Direct M2M translation by 4 rotations and 1 translation
Compute the following matrix-vector product: \f[ \mathrm{dst} = \beta \mathrm{dst} + \alpha M_M \mathrm{src}, \f] where \f$ \mathrm{dst} \f$ is a vector of coefficients of output spherical harmonics, \f$ \mathrm{src} \f$ is a vector of coefficients of input spherical harmonics and \f$ M_M \f$ is a matrix of a multipole-to-multipole translation.
Rotates around OZ and OY axes, translates over OZ and then rotates back around OY and OZ axes.
@param[in] c: Radius-vector from new to old centers of harmonics
@param[in] src_r: Radius of old harmonics
@param[in] dst_r: Radius of new harmonics
@param[in] p: Maximal degree of spherical harmonics
@param[in] vscales: Normalization constants for Y_lm
@param[in] vcnk: Square roots of combinatorial numbers C_n^k
@param[in] alpha: Scalar multiplier for src_m
@param[in] src_m: Expansion in old harmonics
@param[in] beta: Scalar multiplier for dst_m
@param[inout] dst_m: Expansion in new harmonics
@param[out] work: Temporary workspace of a size 6pp+19*p+8
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=rp), | intent(in) | :: | c(3) | |||
integer, | intent(in) | :: | p | |||
real(kind=rp), | intent(in) | :: | vscales((p+1)*(p+1)) | |||
real(kind=rp), | intent(in) | :: | vcnk((2*p+1)*(p+1)) | |||
real(kind=rp), | intent(in) | :: | alpha | |||
real(kind=rp), | intent(in) | :: | src_m((p+1)*(p+1)) | |||
real(kind=rp), | intent(in) | :: | beta | |||
real(kind=rp), | intent(inout) | :: | dst_m((p+1)*(p+1)) | |||
real(kind=rp), | intent(out), | target | :: | work(6*p*p+19*p+8) |
subroutine fmm_m2m_rotation_work(c, p, vscales, vcnk, alpha, &
& src_m, beta, dst_m, work)
! Inputs
integer, intent(in) :: p
real(rp), intent(in) :: c(3), vscales((p+1)*(p+1)), &
& vcnk((2*p+1)*(p+1)), alpha, src_m((p+1)*(p+1)), beta
! Output
real(rp), intent(inout) :: dst_m((p+1)*(p+1))
! Temporary workspace
real(rp), intent(out), target :: work(6*p*p + 19*p + 8)
! Local variables
real(rp) :: rho, ctheta, stheta, cphi, sphi
integer :: m, n
! Pointers for temporary values of harmonics
real(rp), pointer :: tmp_m(:), tmp_m2(:), vcos(:), vsin(:)
! Convert Cartesian coordinates into spherical
call carttosph(c, rho, ctheta, stheta, cphi, sphi)
! If no need for rotations, just do translation along z
if (abs(stheta) < eps_rp) then
! Workspace here is 2*(p+1)
call fmm_m2m_ztranslate_work(c(3), p, vscales, vcnk, &
& alpha, src_m, beta, dst_m, work)
return
end if
! Prepare pointers
m = (p+1)**2
n = 4*m + 5*p ! 4*p*p + 13*p + 4
tmp_m(1:m) => work(n+1:n+m) ! 5*p*p + 15*p + 5
n = n + m
tmp_m2(1:m) => work(n+1:n+m) ! 6*p*p + 17*p + 6
n = n + m
m = p + 1
vcos => work(n+1:n+m) ! 6*p*p + 18*p + 7
n = n + m
vsin => work(n+1:n+m) ! 6*p*p + 19*p + 8
! Compute arrays of cos and sin that are needed for rotations of harmonics
call trgev(cphi, sphi, p, vcos, vsin)
! Rotate around OZ axis (work array might appear in the future)
call fmm_sph_rotate_oz_adj_work(p, vcos, vsin, alpha, src_m, 0.0_rp, tmp_m)
! Perform rotation in the OXZ plane, work size is 4*p*p+13*p+4
call fmm_sph_rotate_oxz_work(p, ctheta, -stheta, 1.0_rp, tmp_m, 0.0_rp, &
& tmp_m2, work)
! OZ translation, workspace here is 2*(p+1)
call fmm_m2m_ztranslate_work(rho, p, vscales, vcnk, 1.0_rp, &
& tmp_m2, 0.0_rp, tmp_m, work)
! Backward rotation in the OXZ plane, work size is 4*p*p+13*p+4
call fmm_sph_rotate_oxz_work(p, ctheta, stheta, 1.0_rp, tmp_m, 0.0_rp, tmp_m2, &
& work)
! Backward rotation around OZ axis (work array might appear in the future)
call fmm_sph_rotate_oz_work(p, vcos, vsin, 1.0_rp, tmp_m2, beta, dst_m)
end subroutine fmm_m2m_rotation_work