make_vscales Subroutine

private subroutine make_vscales(p)

Compute scaling factors of real normalized spherical harmonics

Output values of scaling factors of \f$ Y_\ell^m \f$ harmonics are filled only for non-negative \f$ m \f$ since scaling factor of \f$ Y_\ell^{-m} \f$ is the same as scaling factor of \f$ Y_\ell^m \f$.

@param[in] p: Maximal degree of spherical harmonics. p >= 0 @param[out] vscales: Array of scaling factors. Dimension is (p+1)**2 @param[out] vscales: Array of values 4pi/(2l+1). Dimension is p+1 @param[out] vscales_rel: Array of relative scaling factors. Dimension is (p+1)**, vscales2.

Arguments

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

Called by

proc~~make_vscales~~CalledByGraph proc~make_vscales make_vscales proc~prepare_fmmm_constants prepare_fmmm_constants proc~prepare_fmmm_constants->proc~make_vscales 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_vscales(p)
        ! Input
        integer, intent(in) :: p
        real(rp) :: tmp, twolp1, tmp2
        integer :: l, ind, m

        if(p <= vscales_p) then
            ! Job is already done
            return
        end if 

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

        vscales_p = p
        allocate(vscales((vscales_p+1)**2))
        allocate(vscales_rel((vscales_p+1)**2))

        twolp1 = 1.0_rp
        do l = 0, p
            ! m = 0
            ind = l*l + l + 1
            tmp = 4.0*pi / twolp1
            tmp2 = tmp
            tmp = sqrt(tmp)
            vscales_rel(ind) = tmp
            vscales(ind) = 1.0_rp / tmp
            twolp1 = twolp1 + 2.0_rp
            tmp = vscales(ind) * sqrt(2.0)
            ! m != 0
            do m = 1, l
                tmp = -tmp / sqrt(dble((l-m+1)*(l+m)))
                vscales(ind+m) = tmp
                vscales(ind-m) = tmp
                vscales_rel(ind+m) = tmp * tmp2
                vscales_rel(ind-m) = vscales_rel(ind+m)
            end do
        end do
    end subroutine make_vscales