make_vcnk Subroutine

private subroutine make_vcnk(dmax)

Compute FMM-related constants

@param[in] dmax: Maximal degree of spherical harmonics to be evaluated. dmax >= 0 @param[in] pm: Maximal degree of the multipole expansion. pm >= 0. @param[in] pl: Maximal degree of the local expansion. pl >= 0. @param[out] vcnk: Array of squre roots of combinatorial factors C_n^k. Dimension is (2*dmax+1)*(dmax+1). @param[out] m2l_ztranslate_coef: Constants for M2L translation over OZ axis. Dimension is (pm+1, pl+1, pl+1). @param[out] m2l_ztranslate_coef: Constants for adjoint M2L translation over OZ axis. Dimension is (pl+1, pl+1, pm+1).

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: dmax

Called by

proc~~make_vcnk~~CalledByGraph proc~make_vcnk make_vcnk proc~prepare_fmmm_constants prepare_fmmm_constants proc~prepare_fmmm_constants->proc~make_vcnk proc~make_m2l_ztranslate_coef make_m2l_ztranslate_coef proc~prepare_fmmm_constants->proc~make_m2l_ztranslate_coef proc~make_m2l_ztranslate_coef->proc~make_vcnk proc~fmm_init fmm_init proc~fmm_init->proc~prepare_fmmm_constants proc~field_extd2d field_extD2D proc~field_extd2d->proc~fmm_init proc~fmm_coordinates_update fmm_coordinates_update proc~fmm_coordinates_update->proc~fmm_init proc~tmatvec_otf TMatVec_otf proc~tmatvec_otf->proc~field_extd2d proc~c_ommp_set_fmm_distance C_ommp_set_fmm_distance proc~c_ommp_set_fmm_distance->proc~fmm_coordinates_update proc~update_coordinates update_coordinates proc~update_coordinates->proc~fmm_coordinates_update proc~c_ommp_set_fmm_min_cell_size C_ommp_set_fmm_min_cell_size proc~c_ommp_set_fmm_min_cell_size->proc~fmm_coordinates_update proc~mmpol_prepare mmpol_prepare proc~mmpol_prepare->proc~fmm_coordinates_update proc~c_ommp_update_coordinates C_ommp_update_coordinates proc~c_ommp_update_coordinates->proc~update_coordinates proc~mmpol_init_from_mmp mmpol_init_from_mmp proc~mmpol_init_from_mmp->proc~mmpol_prepare proc~ommp_system_from_qm_helper ommp_system_from_qm_helper proc~ommp_system_from_qm_helper->proc~mmpol_prepare proc~mmpol_init_from_xyz mmpol_init_from_xyz proc~mmpol_init_from_xyz->proc~mmpol_prepare proc~ommp_init_mmp ommp_init_mmp proc~ommp_init_mmp->proc~mmpol_init_from_mmp proc~c_ommp_system_from_qm_helper C_ommp_system_from_qm_helper proc~c_ommp_system_from_qm_helper->proc~ommp_system_from_qm_helper proc~ommp_init_xyz ommp_init_xyz proc~ommp_init_xyz->proc~mmpol_init_from_xyz proc~c_ommp_init_mmp C_ommp_init_mmp proc~c_ommp_init_mmp->proc~ommp_init_mmp proc~c_ommp_init_xyz C_ommp_init_xyz proc~c_ommp_init_xyz->proc~ommp_init_xyz

Contents

Source Code


Source Code

subroutine make_vcnk(dmax)
    ! Inputs
    integer, intent(in) :: dmax
    ! Local variables

    integer(ip) :: i, j, indi

    if(dmax <= vcnk_dmax) then
        ! Job is already done
        return
    end if 

    if(allocated(vcnk)) then
        ! vscales is allocated but it should be expanded.
        ! Just remove everything and restart from scratch
        deallocate(vcnk)
    end if

    vcnk_dmax = dmax
    allocate(vcnk((2*dmax+1)*(dmax+1)))

    vcnk(1) = 1.0_rp
    do i = 2, 2*dmax+1
        ! Offset to the C_{i-2}^{i-2}, next item to be stored is C_{i-1}^0
        indi = (i-1) * i / 2
        ! C_{i-1}^0 = 1
        vcnk(indi+1) = 1.0_rp
        ! C_{i-1}^{i-1} = 1
        vcnk(indi+i) = 1.0_rp
        ! C_{i-1}^{j-1} = C_{i-2}^{j-1} + C_{i-2}^{j-2}
        ! Offset to C_{i-3}^{i-3} is indi-i+1
        do j = 2, i-1
            vcnk(indi+j) = vcnk(indi-i+j+1) + vcnk(indi-i+j)
        end do
    end do
    ! Get square roots of C_n^k. sqrt(1.0_rp) is 1.0_rp, so no need to update C_n^0
    ! and C_n^n
    do i = 3, 2*dmax+1
        indi = (i-1) * i / 2
        do j = 2, i-1
            vcnk(indi+j) = sqrt(vcnk(indi+j))
        end do
    end do
end subroutine make_vcnk