fmm_sph_rotate_oz_adj_work Subroutine

private subroutine fmm_sph_rotate_oz_adj_work(p, vcos, vsin, alpha, src, beta, dst)

Rotate spherical harmonics around OZ axis in an opposite direction

Compute the following matrix-vector product: \f[ \mathrm{dst} = \beta \mathrm{dst} + \alpha R \mathrm{src}, \f] where \f$ \mathrm{dst} \f$ is a vector of coefficients of output spherical harmonics corresponding to a new cartesion system of coordinates, \f$ \mathrm{src} \f$ is a vector of coefficients of input spherical harmonics corresponding to the standard cartesian system of coordinates and \f$ R \f$ is a matrix of rotation of coordinates around OZ axis on angle \f$ \phi \f$, presented by \f$ \cos(m \phi) \f$ and \f$ \sin(m \phi) \f$.

@param[in] p: Maximal order of spherical harmonics @param[in] vcos: Vector \f$ { \cos(m \phi) }{m=0}^p \f$ @param[in] vsin: Vector \f$ { \sin(m \phi) }^p \f$ @param[in] alpha: Scalar multiplier for src @param[in] src: Coefficients of initial spherical harmonics @param[in] beta: Scalar multipler for dst @param[inout] dst: Coefficients of rotated spherical harmonics

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: p
real(kind=rp), intent(in) :: vcos(p+1)
real(kind=rp), intent(in) :: vsin(p+1)
real(kind=rp), intent(in) :: alpha
real(kind=rp), intent(in) :: src((p+1)*(p+1))
real(kind=rp), intent(in) :: beta
real(kind=rp), intent(inout) :: dst((p+1)*(p+1))

Called by

proc~~fmm_sph_rotate_oz_adj_work~~CalledByGraph proc~fmm_sph_rotate_oz_adj_work fmm_sph_rotate_oz_adj_work proc~fmm_m2m_rotation_work fmm_m2m_rotation_work proc~fmm_m2m_rotation_work->proc~fmm_sph_rotate_oz_adj_work proc~fmm_l2l_rotation_work fmm_l2l_rotation_work proc~fmm_l2l_rotation_work->proc~fmm_sph_rotate_oz_adj_work proc~fmm_m2l_rotation_work fmm_m2l_rotation_work proc~fmm_m2l_rotation_work->proc~fmm_sph_rotate_oz_adj_work proc~fmm_m2m fmm_m2m proc~fmm_m2m->proc~fmm_m2m_rotation_work proc~fmm_l2l fmm_l2l proc~fmm_l2l->proc~fmm_l2l_rotation_work proc~fmm_m2l fmm_m2l proc~fmm_m2l->proc~fmm_m2l_rotation_work proc~tree_p2m tree_p2m proc~tree_p2m->proc~fmm_m2m proc~cart_propfar_at_ipart cart_propfar_at_ipart proc~cart_propfar_at_ipart->proc~fmm_l2l proc~tree_m2m tree_m2m proc~tree_m2m->proc~fmm_m2m proc~tree_m2l tree_m2l proc~tree_m2l->proc~fmm_m2l proc~tree_l2l tree_l2l proc~tree_l2l->proc~fmm_l2l proc~cart_propnear_at_ipart cart_propnear_at_ipart proc~cart_propnear_at_ipart->proc~fmm_m2l proc~fmm_solve_for_multipoles fmm_solve_for_multipoles proc~fmm_solve_for_multipoles->proc~tree_p2m proc~fmm_solve_for_multipoles->proc~tree_m2m proc~fmm_solve_for_multipoles->proc~tree_m2l proc~fmm_solve_for_multipoles->proc~tree_l2l proc~elec_prop_d2m elec_prop_D2M proc~elec_prop_d2m->proc~cart_propfar_at_ipart proc~prepare_fmm_ipd prepare_fmm_ipd proc~elec_prop_d2m->proc~prepare_fmm_ipd proc~field_extd2d field_extD2D proc~field_extd2d->proc~cart_propfar_at_ipart proc~prepare_fmm_ext_ipd prepare_fmm_ext_ipd proc~field_extd2d->proc~prepare_fmm_ext_ipd proc~elec_prop_d2d elec_prop_D2D proc~elec_prop_d2d->proc~cart_propfar_at_ipart proc~elec_prop_d2d->proc~prepare_fmm_ipd proc~elec_prop_m2m elec_prop_M2M proc~elec_prop_m2m->proc~cart_propfar_at_ipart proc~preapare_fmm_static preapare_fmm_static proc~elec_prop_m2m->proc~preapare_fmm_static proc~elec_prop_m2d elec_prop_M2D proc~elec_prop_m2d->proc~cart_propfar_at_ipart proc~elec_prop_m2d->proc~preapare_fmm_static proc~cart_prop_at_ipart cart_prop_at_ipart proc~cart_prop_at_ipart->proc~cart_propfar_at_ipart proc~cart_prop_at_ipart->proc~cart_propnear_at_ipart proc~fmm_solve fmm_solve proc~fmm_solve->proc~tree_m2m proc~fmm_solve->proc~tree_m2l proc~fmm_solve->proc~tree_l2l proc~prepare_fmm_ext_ipd->proc~fmm_solve_for_multipoles proc~preapare_fmm_static->proc~fmm_solve_for_multipoles proc~prepare_polelec prepare_polelec proc~prepare_polelec->proc~elec_prop_d2m proc~prepare_polelec->proc~elec_prop_d2d proc~prepare_polelec->proc~elec_prop_m2d proc~prepare_fixedelec prepare_fixedelec proc~prepare_fixedelec->proc~elec_prop_m2m proc~tmatvec_otf TMatVec_otf proc~tmatvec_otf->proc~field_extd2d proc~prepare_fmm_ipd->proc~prepare_fmm_ext_ipd proc~ommp_get_polelec_energy ommp_get_polelec_energy proc~ommp_get_polelec_energy->proc~prepare_polelec proc~energy_mm_pol energy_MM_pol proc~ommp_get_polelec_energy->proc~energy_mm_pol proc~ommp_set_external_field ommp_set_external_field proc~ommp_set_external_field->proc~prepare_polelec proc~polelec_geomgrad polelec_geomgrad proc~polelec_geomgrad->proc~prepare_polelec proc~energy_mm_pol->proc~prepare_polelec proc~fixedelec_geomgrad fixedelec_geomgrad proc~fixedelec_geomgrad->proc~prepare_fixedelec proc~energy_mm_mm energy_MM_MM proc~energy_mm_mm->proc~prepare_fixedelec proc~ommp_get_full_ele_energy ommp_get_full_ele_energy proc~ommp_get_full_ele_energy->proc~ommp_get_polelec_energy proc~ommp_get_fixedelec_energy ommp_get_fixedelec_energy proc~ommp_get_full_ele_energy->proc~ommp_get_fixedelec_energy proc~c_ommp_get_polelec_energy C_ommp_get_polelec_energy proc~c_ommp_get_polelec_energy->proc~ommp_get_polelec_energy proc~ommp_set_external_field_nomm ommp_set_external_field_nomm proc~ommp_set_external_field_nomm->proc~ommp_set_external_field proc~c_ommp_set_external_field C_ommp_set_external_field proc~c_ommp_set_external_field->proc~ommp_set_external_field proc~c_ommp_set_external_field_nomm C_ommp_set_external_field_nomm proc~c_ommp_set_external_field_nomm->proc~ommp_set_external_field proc~ommp_polelec_geomgrad ommp_polelec_geomgrad proc~ommp_polelec_geomgrad->proc~polelec_geomgrad proc~ommp_full_geomgrad ommp_full_geomgrad proc~ommp_full_geomgrad->proc~polelec_geomgrad proc~ommp_full_geomgrad->proc~fixedelec_geomgrad proc~ommp_fixedelec_geomgrad ommp_fixedelec_geomgrad proc~ommp_fixedelec_geomgrad->proc~fixedelec_geomgrad proc~ommp_get_fixedelec_energy->proc~energy_mm_mm proc~ommp_get_full_energy ommp_get_full_energy proc~ommp_get_full_energy->proc~ommp_get_full_ele_energy proc~c_ommp_polelec_geomgrad C_ommp_polelec_geomgrad proc~c_ommp_polelec_geomgrad->proc~ommp_polelec_geomgrad proc~c_ommp_full_geomgrad C_ommp_full_geomgrad proc~c_ommp_full_geomgrad->proc~ommp_full_geomgrad proc~c_ommp_fixedelec_geomgrad C_ommp_fixedelec_geomgrad proc~c_ommp_fixedelec_geomgrad->proc~ommp_fixedelec_geomgrad proc~c_ommp_get_full_ele_energy C_ommp_get_full_ele_energy proc~c_ommp_get_full_ele_energy->proc~ommp_get_full_ele_energy proc~c_ommp_get_fixedelec_energy C_ommp_get_fixedelec_energy proc~c_ommp_get_fixedelec_energy->proc~ommp_get_fixedelec_energy proc~c_ommp_get_full_energy C_ommp_get_full_energy proc~c_ommp_get_full_energy->proc~ommp_get_full_energy

Contents


Source Code

    subroutine fmm_sph_rotate_oz_adj_work(p, vcos, vsin, alpha, src, beta, dst)
        ! Inputs
        integer, intent(in) :: p
        real(rp), intent(in) :: vcos(p+1), vsin(p+1), alpha, src((p+1)*(p+1)), beta
        ! Output
        real(rp), intent(inout) :: dst((p+1)*(p+1))
        ! Local variables
        integer :: l, m, ind
        real(rp) :: v1, v2, v3, v4
        ! In case alpha is 0.0 just scale output
        if (abs(alpha) < eps_rp) then
            ! Set output to 0.0 if beta is also 0.0
            if (abs(beta) < eps_rp) then
                dst = 0.0
            else
                dst = beta * dst
            end if
            ! Exit subroutine
            return
        end if
        ! Now alpha is non-0.0
        ! In case beta is 0.0 output is just overwritten without being read
        if (abs(beta) < eps_rp) then
            ! l = 0
            dst(1) = alpha*src(1)
            ! l > 0
            do l = 1, p
                ind = l*l + l + 1
                ! m = 0
                dst(ind) = alpha*src(ind)
                ! m != 0
                do m = 1, l
                    v1 = src(ind+m)
                    v2 = src(ind-m)
                    v3 = vcos(1+m)
                    v4 = vsin(1+m)
                    ! m > 0
                    dst(ind+m) = alpha * (v1*v3+v2*v4)
                    ! m < 0
                    dst(ind-m) = alpha * (v2*v3-v1*v4)
                end do
            end do
        else
            ! l = 0
            dst(1) = beta*dst(1) + alpha*src(1)
            ! l > 0
            do l = 1, p
                ind = l*l + l + 1
                ! m = 0
                dst(ind) = beta*dst(ind) + alpha*src(ind)
                ! m != 0
                do m = 1, l
                    v1 = src(ind+m)
                    v2 = src(ind-m)
                    v3 = vcos(1+m)
                    v4 = vsin(1+m)
                    ! m > 0
                    dst(ind+m) = beta*dst(ind+m) + alpha*(v1*v3+v2*v4)
                    ! m < 0
                    dst(ind-m) = beta*dst(ind-m) + alpha*(v2*v3-v1*v4)
                end do
            end do
        end if
    end subroutine fmm_sph_rotate_oz_adj_work