ddx 0.6.0
Libary for domain-decomposition methods for polarizable continuum models
ddx_core Module Reference

Core routines and parameters of the ddX software. More...

Data Types

type  ddx_electrostatics_type
 Container for the electrostatic properties. Since different methods require different electrostatic properties, defining this kind of general purpose container simplify the interfaces as we don't need different interfaces for different methods. More...
 
type  ddx_state_type
 This defined type contains the primal and adjoint RHSs, the solution of the primal and adjoint linear systems, useful intermediates for the computation of the forces, and the information about the convergence of the linear system solver (time, number of iterations, relative difference at each iteration). More...
 
type  ddx_type
 Main ddX type that stores all the required information. Container for the params, contants and workspace derived types. More...
 

Functions/Subroutines

subroutine allocate_model (nsph, x, y, z, rvdw, model, lmax, ngrid, force, fmm, pm, pl, se, eta, eps, kappa, matvecmem, maxiter, jacobi_ndiis, nproc, output_filename, ddx_data, ddx_error)
 Initialize ddX input with a full set of parameters. More...
 
subroutine deallocate_model (ddx_data, ddx_error)
 Deallocate object with corresponding data. More...
 
subroutine allocate_electrostatics (params, constants, electrostatics, ddx_error)
 Given the chosen model, find the required electrostatic properties, and allocate the arrays for them in the container. More...
 
subroutine deallocate_electrostatics (electrostatics, ddx_error)
 Deallocate the electrostatic properties. More...
 
subroutine allocate_state (params, constants, state, ddx_error)
 Initialize the ddx_state object. More...
 
subroutine deallocate_state (state, ddx_error)
 Deallocate the ddx_state object. More...
 
subroutine print_spherical (iunit, label, nbasis, lmax, ncol, icol, x)
 Print array of spherical harmonics. More...
 
subroutine print_nodes (iunit, label, ngrid, ncol, icol, x)
 Print array of quadrature points. More...
 
subroutine print_ddvector (ddx_data, label, vector)
 Print dd Solution vector. More...
 
subroutine ddintegrate (nsph, nbasis, ngrid, vwgrid, ldvwgrid, x_grid, x_lm)
 Integrate against spherical harmonics. More...
 
subroutine dbasis (params, constants, x, basloc, dbsloc, vplm, vcos, vsin)
 Compute first derivatives of spherical harmonics. More...
 
real(dp) function intmlp (params, constants, t, sigma, basloc)
 TODO. More...
 
subroutine wghpot (ncav, phi_cav, nsph, ngrid, ui, phi_grid, g)
 Weigh potential at cavity points by characteristic function TODO use cavity points in CSR format. More...
 
subroutine hsnorm (lmax, nbasis, u, unorm)
 Compute the local Sobolev H^(-1/2)-norm on one sphere of u. More...
 
real(dp) function hnorm (lmax, nbasis, nsph, x)
 Compute the global Sobolev H^(-1/2)-norm of x. More...
 
real(dp) function rmsnorm (lmax, nbasis, nsph, x)
 
subroutine rmsvec (n, v, vrms, vmax)
 compute root-mean-square and max norm More...
 
subroutine adjrhs (params, constants, isph, xi, vlm, work)
 
subroutine calcv (params, constants, isph, pot, sigma, work)
 
subroutine ddeval_grid (params, constants, alpha, x_sph, beta, x_grid)
 Evaluate values of spherical harmonics at Lebedev grid points. More...
 
subroutine ddeval_grid_work (nbasis, ngrid, nsph, vgrid, ldvgrid, alpha, x_sph, beta, x_grid)
 Evaluate values of spherical harmonics at Lebedev grid points. More...
 
subroutine ddintegrate_sph (params, constants, alpha, x_grid, beta, x_sph)
 Integrate values at grid points into spherical harmonics. More...
 
subroutine ddintegrate_sph_work (nbasis, ngrid, nsph, vwgrid, ldvwgrid, alpha, x_grid, beta, x_sph)
 Integrate values at grid points into spherical harmonics. More...
 
subroutine ddcav_to_grid (params, constants, x_cav, x_grid)
 Unwrap values at cavity points into values at all grid points. More...
 
subroutine ddcav_to_grid_work (ngrid, nsph, ncav, icav_ia, icav_ja, x_cav, x_grid)
 Unwrap values at cavity points into values at all grid points. More...
 
subroutine ddproject_cav (params, constants, s, xi)
 Integrate by a characteristic function at Lebedev grid points \xi(n,i) = sum w_n U_n^i Y_l^m(s_n) [S_i]_l^m l,m. More...
 
subroutine tree_m2m_rotation (params, constants, node_m)
 Transfer multipole coefficients over a tree. More...
 
subroutine tree_m2m_rotation_work (params, constants, node_m, work)
 Transfer multipole coefficients over a tree. More...
 
subroutine tree_m2m_bessel_rotation (params, constants, node_m)
 Transfer multipole coefficients over a tree. More...
 
subroutine tree_m2m_bessel_rotation_work (params, constants, node_m)
 Transfer multipole coefficients over a tree. More...
 
subroutine tree_m2m_rotation_adj (params, constants, node_m)
 Adjoint transfer multipole coefficients over a tree. More...
 
subroutine tree_m2m_rotation_adj_work (params, constants, node_m, work)
 Adjoint transfer multipole coefficients over a tree. More...
 
subroutine tree_m2m_bessel_rotation_adj (params, constants, node_m)
 Adjoint transfer multipole coefficients over a tree. More...
 
subroutine tree_m2m_bessel_rotation_adj_work (params, constants, node_m)
 Adjoint transfer multipole coefficients over a tree. More...
 
subroutine tree_l2l_rotation (params, constants, node_l)
 Transfer local coefficients over a tree. More...
 
subroutine tree_l2l_rotation_work (params, constants, node_l, work)
 Transfer local coefficients over a tree. More...
 
subroutine tree_l2l_bessel_rotation (params, constants, node_l)
 Transfer local coefficients over a tree. More...
 
subroutine tree_l2l_bessel_rotation_work (params, constants, node_l)
 Transfer local coefficients over a tree. More...
 
subroutine tree_l2l_rotation_adj (params, constants, node_l)
 Adjoint transfer local coefficients over a tree. More...
 
subroutine tree_l2l_rotation_adj_work (params, constants, node_l, work)
 Adjoint transfer local coefficients over a tree. More...
 
subroutine tree_l2l_bessel_rotation_adj (params, constants, node_l)
 Adjoint transfer local coefficients over a tree. More...
 
subroutine tree_m2l_rotation (params, constants, node_m, node_l)
 Transfer multipole local coefficients into local over a tree. More...
 
subroutine tree_m2l_bessel_rotation (params, constants, node_m, node_l)
 Transfer multipole local coefficients into local over a tree. More...
 
subroutine tree_m2l_bessel_rotation_adj (params, constants, node_l, node_m)
 Adjoint transfer multipole local coefficients into local over a tree. More...
 
subroutine tree_m2l_bessel_rotation_adj_work (params, constants, node_l, node_m)
 Adjoint transfer multipole local coefficients into local over a tree. More...
 
subroutine tree_m2l_rotation_adj (params, constants, node_l, node_m)
 Adjoint transfer multipole local coefficients into local over a tree. More...
 
subroutine tree_l2p (params, constants, alpha, node_l, beta, grid_v, sph_l)
 TODO. More...
 
subroutine tree_l2p_bessel (params, constants, alpha, node_l, beta, grid_v)
 TODO. More...
 
subroutine tree_l2p_adj (params, constants, alpha, grid_v, beta, node_l, sph_l)
 TODO. More...
 
subroutine tree_l2p_bessel_adj (params, constants, alpha, grid_v, beta, node_l)
 TODO. More...
 
subroutine tree_m2p (params, constants, p, alpha, sph_m, beta, grid_v)
 TODO. More...
 
subroutine tree_m2p_bessel (params, constants, p, alpha, sph_p, sph_m, beta, grid_v)
 TODO. More...
 
subroutine tree_m2p_adj (params, constants, p, alpha, grid_v, beta, sph_m)
 TODO. More...
 
subroutine tree_m2p_bessel_adj (params, constants, p, alpha, grid_v, beta, sph_p, sph_m)
 TODO. More...
 
subroutine tree_m2p_bessel_nodiag_adj (params, constants, p, alpha, grid_v, beta, sph_p, sph_m)
 TODO. More...
 
subroutine tree_grad_m2m (params, constants, sph_m, sph_m_grad, work)
 TODO. More...
 
subroutine tree_grad_l2l (params, constants, node_l, sph_l_grad, work)
 TODO. More...
 
subroutine get_banner (string)
 Print the ddX logo. More...
 
subroutine cav_to_spherical (params, constants, workspace, property_cav, property_sph)
 Transform a function defined at the exposed cavity points (cav) to a spherical harmonics expansion. Note that the function is also multiplied by the characteristic function U. More...
 

Detailed Description

Core routines and parameters of the ddX software.

Function/Subroutine Documentation

◆ allocate_model()

subroutine ddx_core::allocate_model ( integer, intent(in)  nsph,
real(dp), dimension(nsph), intent(in)  x,
real(dp), dimension(nsph), intent(in)  y,
real(dp), dimension(nsph), intent(in)  z,
real(dp), dimension(nsph), intent(in)  rvdw,
integer, intent(in)  model,
integer, intent(in)  lmax,
integer, intent(in)  ngrid,
integer, intent(in)  force,
integer, intent(in)  fmm,
integer, intent(in)  pm,
integer, intent(in)  pl,
real(dp), intent(in)  se,
real(dp), intent(in)  eta,
real(dp), intent(in)  eps,
real(dp), intent(in)  kappa,
integer, intent(in)  matvecmem,
integer, intent(in)  maxiter,
integer, intent(in)  jacobi_ndiis,
integer, intent(in)  nproc,
character(len=255), intent(in)  output_filename,
type(ddx_type), intent(out), target  ddx_data,
type(ddx_error_type), intent(inout)  ddx_error 
)

Initialize ddX input with a full set of parameters.

Parameters
[in]nsphNumber of atoms. n > 0.
[in]x\( x \) coordinates of atoms. Dimension is (n)
[in]y\( y \) coordinates of atoms. Dimension is (n)
[in]z\( z \) coordinates of atoms. Dimension is (n)
[in]rvdwVan-der-Waals radii of atoms. Dimension is (n)
[in]modelChoose model: 1 for COSMO, 2 for PCM and 3 for LPB
[in]lmaxMaximal degree of modeling spherical harmonics. lmax >= 0
[in]ngridNumber of Lebedev grid points ngrid >= 0
[in]force1 if forces are required and 0 otherwise
[in]fmm1 to use FMM acceleration and 0 otherwise
[in]pmMaximal degree of multipole spherical harmonics. Ignored in the case fmm=0. Value -1 means no far-field FMM interactions are computed. pm >= -1.
[in]plMaximal degree of local spherical harmonics. Ignored in the case fmm=0. Value -1 means no far-field FMM interactions are computed. pl >= -1.
[in]seShift of characteristic function. -1 for interior, 0 for centered and 1 for outer regularization
[in]etaRegularization parameter. 0 < eta <= 1.
[in]epsRelative dielectric permittivity. eps > 1.
[in]kappaDebye-H"{u}ckel parameter
[in]matvecmemhandling of the sparse matrices. 1 for precomputin them and keeping them in memory, 0 for assembling the matrix-vector product on-the-fly.
[in]maxiterMaximum number of iterations for an iterative solver. maxiter > 0.
[in]ndiisNumber of extrapolation points for Jacobi/DIIS solver. ndiis >= 0.
[in,out]nprocNumber of OpenMP threads to be used where applicable.
[out]ddx_dataObject containing all inputs
[in,out]ddx_errorddX error

Definition at line 258 of file ddx_core.f90.

◆ deallocate_model()

subroutine ddx_core::deallocate_model ( type(ddx_type), intent(inout)  ddx_data,
type(ddx_error_type), intent(inout)  ddx_error 
)

Deallocate object with corresponding data.

Parameters
[in,out]ddx_dataobject to deallocate
[in,out]ddx_errorddX error

Definition at line 320 of file ddx_core.f90.

◆ print_spherical()

subroutine ddx_core::print_spherical ( integer, intent(in)  iunit,
character (len=*), intent(in)  label,
integer, intent(in)  nbasis,
integer, intent(in)  lmax,
integer, intent(in)  ncol,
integer, intent(in)  icol,
real(dp), dimension(nbasis, ncol), intent(in)  x 
)

Print array of spherical harmonics.

Prints (nbasis, ncol) array

Parameters
[in]labelLabel to print
[in]nbasisNumber of rows of input x. nbasis >= 0
[in]lmaxMaximal degree of corresponding spherical harmonics. (lmax+1)**2 = nbasis
[in]ncolNumber of columns to print
[in]icolThis number is only for printing purposes
[in]xActual data to print

Definition at line 996 of file ddx_core.f90.

◆ print_nodes()

subroutine ddx_core::print_nodes ( integer, intent(in)  iunit,
character (len=*), intent(in)  label,
integer, intent(in)  ngrid,
integer, intent(in)  ncol,
integer, intent(in)  icol,
real(dp), dimension(ngrid, ncol), intent(in)  x 
)

Print array of quadrature points.

Prints (ngrid, ncol) array

Parameters
[in]labelLabel to print
[in]ngridNumber of rows of input x. ngrid >= 0
[in]ncolNumber of columns to print
[in]icolThis number is only for printing purposes
[in]xActual data to print

Definition at line 1053 of file ddx_core.f90.

◆ print_ddvector()

subroutine ddx_core::print_ddvector ( type(ddx_type), intent(in)  ddx_data,
character(len=*)  label,
real(dp), dimension(ddx_data % constants % nbasis, ddx_data % params % nsph)  vector 
)

Print dd Solution vector.

Parameters
[in]label: Label to print
[in]vectorVector to print

Definition at line 1099 of file ddx_core.f90.

◆ ddintegrate()

subroutine ddx_core::ddintegrate ( integer, intent(in)  nsph,
integer, intent(in)  nbasis,
integer, intent(in)  ngrid,
real(dp), dimension(ldvwgrid, ngrid), intent(in)  vwgrid,
integer, intent(in)  ldvwgrid,
real(dp), dimension(ngrid, nsph), intent(in)  x_grid,
real(dp), dimension(nbasis, nsph), intent(out)  x_lm 
)

Integrate against spherical harmonics.

Integrate by Lebedev spherical quadrature. This function can be simply substituted by a matrix-vector product.

Parameters
[in]nsphNumber of spheres. nsph >= 0.
[in]nbasisNumber of spherical harmonics. nbasis >= 0.
[in]ngridNumber of Lebedev grid points. ngrid >= 0.
[in]vwgridValues of spherical harmonics at Lebedev grid points, multiplied by weights of grid points. Dimension is (ldvwgrid, ngrid).
[in]ldvwgridLeading dimension of vwgrid.
[in]x_gridInput values at grid points of the sphere. Dimension is (ngrid, nsph).
[out]x_lmOutput spherical harmonics. Dimension is (nbasis, nsph).

Definition at line 1132 of file ddx_core.f90.

◆ dbasis()

subroutine ddx_core::dbasis ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension(3), intent(in)  x,
real(dp), dimension(constants % nbasis), intent(inout)  basloc,
real(dp), dimension(3, constants % nbasis), intent(inout)  dbsloc,
real(dp), dimension(constants % nbasis), intent(inout)  vplm,
real(dp), dimension(params % lmax+1), intent(inout)  vcos,
real(dp), dimension(params % lmax+1), intent(inout)  vsin 
)

Compute first derivatives of spherical harmonics.

Parameters
[in]x
[out]basloc
[out]dbsloc
[out]vplm
[out]vcos
[out]vsin

TODO: rewrite code and fill description. Computing sqrt(one-cthe*cthe) reduces effective range of input double precision values. cthe*cthe for cthe=1d+155 is NaN.

Definition at line 1160 of file ddx_core.f90.

◆ intmlp()

real(dp) function ddx_core::intmlp ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), intent(in)  t,
real(dp), dimension(constants % nbasis), intent(in)  sigma,
real(dp), dimension(constants % nbasis), intent(in)  basloc 
)

TODO.

Definition at line 1291 of file ddx_core.f90.

◆ wghpot()

subroutine ddx_core::wghpot ( integer, intent(in)  ncav,
real(dp), dimension(ncav), intent(in)  phi_cav,
integer, intent(in)  nsph,
integer, intent(in)  ngrid,
real(dp), dimension(ngrid, nsph), intent(in)  ui,
real(dp), dimension(ngrid, nsph), intent(out)  phi_grid,
real(dp), dimension(ngrid, nsph), intent(out)  g 
)

Weigh potential at cavity points by characteristic function TODO use cavity points in CSR format.

Definition at line 1337 of file ddx_core.f90.

◆ hsnorm()

subroutine ddx_core::hsnorm ( integer, intent(in)  lmax,
integer, intent(in)  nbasis,
real(dp), dimension(nbasis), intent(in)  u,
real(dp), intent(inout)  unorm 
)

Compute the local Sobolev H^(-1/2)-norm on one sphere of u.

Definition at line 1368 of file ddx_core.f90.

◆ hnorm()

real(dp) function ddx_core::hnorm ( integer, intent(in)  lmax,
integer, intent(in)  nbasis,
integer, intent(in)  nsph,
real(dp), dimension(nbasis, nsph), intent(in)  x 
)

Compute the global Sobolev H^(-1/2)-norm of x.

Definition at line 1389 of file ddx_core.f90.

◆ rmsnorm()

real(dp) function ddx_core::rmsnorm ( integer, intent(in)  lmax,
integer, intent(in)  nbasis,
integer, intent(in)  nsph,
real(dp), dimension(nbasis, nsph), intent(in)  x 
)

Definition at line 1406 of file ddx_core.f90.

◆ rmsvec()

subroutine ddx_core::rmsvec ( integer, intent(in)  n,
real(dp), dimension(n), intent(in)  v,
real(dp), intent(inout)  vrms,
real(dp), intent(inout)  vmax 
)

compute root-mean-square and max norm

Definition at line 1421 of file ddx_core.f90.

◆ adjrhs()

subroutine ddx_core::adjrhs ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
integer, intent(in)  isph,
real(dp), dimension(params % ngrid, params % nsph), intent(in)  xi,
real(dp), dimension(constants % nbasis), intent(inout)  vlm,
real(dp), dimension(params % lmax+1), intent(inout)  work 
)

Definition at line 1438 of file ddx_core.f90.

◆ calcv()

subroutine ddx_core::calcv ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
integer, intent(in)  isph,
real(dp), dimension(params % ngrid), intent(inout)  pot,
real(dp), dimension(constants % nbasis, params % nsph), intent(in)  sigma,
real(dp), dimension(params % lmax+1), intent(inout)  work 
)

Definition at line 1472 of file ddx_core.f90.

◆ ddeval_grid()

subroutine ddx_core::ddeval_grid ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), intent(in)  alpha,
real(dp), dimension(constants % nbasis, params % nsph), intent(in)  x_sph,
real(dp), intent(in)  beta,
real(dp), dimension(params % ngrid, params % nsph), intent(inout)  x_grid 
)

Evaluate values of spherical harmonics at Lebedev grid points.

Definition at line 1519 of file ddx_core.f90.

◆ ddeval_grid_work()

subroutine ddx_core::ddeval_grid_work ( integer, intent(in)  nbasis,
integer, intent(in)  ngrid,
integer, intent(in)  nsph,
real(dp), dimension(ldvgrid, ngrid), intent(in)  vgrid,
integer, intent(in)  ldvgrid,
real(dp), intent(in)  alpha,
real(dp), dimension(nbasis, nsph), intent(in)  x_sph,
real(dp), intent(in)  beta,
real(dp), dimension(ngrid, nsph), intent(inout)  x_grid 
)

Evaluate values of spherical harmonics at Lebedev grid points.

Definition at line 1537 of file ddx_core.f90.

◆ ddintegrate_sph()

subroutine ddx_core::ddintegrate_sph ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), intent(in)  alpha,
real(dp), dimension(params % ngrid, params % nsph), intent(in)  x_grid,
real(dp), intent(in)  beta,
real(dp), dimension(constants % nbasis, params % nsph), intent(inout)  x_sph 
)

Integrate values at grid points into spherical harmonics.

Definition at line 1554 of file ddx_core.f90.

◆ ddintegrate_sph_work()

subroutine ddx_core::ddintegrate_sph_work ( integer, intent(in)  nbasis,
integer, intent(in)  ngrid,
integer, intent(in)  nsph,
real(dp), dimension(ldvwgrid, ngrid), intent(in)  vwgrid,
integer, intent(in)  ldvwgrid,
real(dp), intent(in)  alpha,
real(dp), dimension(ngrid, nsph), intent(in)  x_grid,
real(dp), intent(in)  beta,
real(dp), dimension(nbasis, nsph), intent(inout)  x_sph 
)

Integrate values at grid points into spherical harmonics.

Definition at line 1570 of file ddx_core.f90.

◆ ddcav_to_grid()

subroutine ddx_core::ddcav_to_grid ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension(constants % ncav), intent(in)  x_cav,
real(dp), dimension(params % ngrid, params % nsph), intent(inout)  x_grid 
)

Unwrap values at cavity points into values at all grid points.

Definition at line 1588 of file ddx_core.f90.

◆ ddcav_to_grid_work()

subroutine ddx_core::ddcav_to_grid_work ( integer, intent(in)  ngrid,
integer, intent(in)  nsph,
integer, intent(in)  ncav,
integer, dimension(nsph+1), intent(in)  icav_ia,
integer, dimension(ncav), intent(in)  icav_ja,
real(dp), dimension(ncav), intent(in)  x_cav,
real(dp), dimension(ngrid, nsph), intent(out)  x_grid 
)

Unwrap values at cavity points into values at all grid points.

Definition at line 1604 of file ddx_core.f90.

◆ ddproject_cav()

subroutine ddx_core::ddproject_cav ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension(constants%nbasis, params%nsph), intent(in)  s,
real(dp), dimension(constants%ncav), intent(out)  xi 
)

Integrate by a characteristic function at Lebedev grid points \xi(n,i) = sum w_n U_n^i Y_l^m(s_n) [S_i]_l^m l,m.

Definition at line 1633 of file ddx_core.f90.

◆ tree_m2m_rotation()

subroutine ddx_core::tree_m2m_rotation ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pm+1)**2, constants % nclusters), intent(inout)  node_m 
)

Transfer multipole coefficients over a tree.

Definition at line 1655 of file ddx_core.f90.

◆ tree_m2m_rotation_work()

subroutine ddx_core::tree_m2m_rotation_work ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pm+1)**2, constants % nclusters), intent(inout)  node_m,
real(dp), dimension(6*params % pm**2 + 19*params % pm + 8), intent(out)  work 
)

Transfer multipole coefficients over a tree.

Definition at line 1671 of file ddx_core.f90.

◆ tree_m2m_bessel_rotation()

subroutine ddx_core::tree_m2m_bessel_rotation ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pm+1)**2, constants % nclusters), intent(inout)  node_m 
)

Transfer multipole coefficients over a tree.

Definition at line 1716 of file ddx_core.f90.

◆ tree_m2m_bessel_rotation_work()

subroutine ddx_core::tree_m2m_bessel_rotation_work ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pm+1)**2, constants % nclusters), intent(inout)  node_m 
)

Transfer multipole coefficients over a tree.

Definition at line 1732 of file ddx_core.f90.

◆ tree_m2m_rotation_adj()

subroutine ddx_core::tree_m2m_rotation_adj ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pm+1)**2, constants % nclusters), intent(inout)  node_m 
)

Adjoint transfer multipole coefficients over a tree.

Definition at line 1777 of file ddx_core.f90.

◆ tree_m2m_rotation_adj_work()

subroutine ddx_core::tree_m2m_rotation_adj_work ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pm+1)**2, constants % nclusters), intent(inout)  node_m,
real(dp), dimension(6*params % pm**2 + 19*params % pm + 8), intent(out)  work 
)

Adjoint transfer multipole coefficients over a tree.

Definition at line 1793 of file ddx_core.f90.

◆ tree_m2m_bessel_rotation_adj()

subroutine ddx_core::tree_m2m_bessel_rotation_adj ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pm+1)**2, constants % nclusters), intent(inout)  node_m 
)

Adjoint transfer multipole coefficients over a tree.

Definition at line 1821 of file ddx_core.f90.

◆ tree_m2m_bessel_rotation_adj_work()

subroutine ddx_core::tree_m2m_bessel_rotation_adj_work ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pm+1)**2, constants % nclusters), intent(inout)  node_m 
)

Adjoint transfer multipole coefficients over a tree.

Definition at line 1836 of file ddx_core.f90.

◆ tree_l2l_rotation()

subroutine ddx_core::tree_l2l_rotation ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(inout)  node_l 
)

Transfer local coefficients over a tree.

Definition at line 1865 of file ddx_core.f90.

◆ tree_l2l_rotation_work()

subroutine ddx_core::tree_l2l_rotation_work ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(inout)  node_l,
real(dp), dimension(6*params % pl**2 + 19*params % pl + 8), intent(out)  work 
)

Transfer local coefficients over a tree.

Definition at line 1881 of file ddx_core.f90.

◆ tree_l2l_bessel_rotation()

subroutine ddx_core::tree_l2l_bessel_rotation ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(inout)  node_l 
)

Transfer local coefficients over a tree.

Definition at line 1912 of file ddx_core.f90.

◆ tree_l2l_bessel_rotation_work()

subroutine ddx_core::tree_l2l_bessel_rotation_work ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(inout)  node_l 
)

Transfer local coefficients over a tree.

Definition at line 1928 of file ddx_core.f90.

◆ tree_l2l_rotation_adj()

subroutine ddx_core::tree_l2l_rotation_adj ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(inout)  node_l 
)

Adjoint transfer local coefficients over a tree.

Definition at line 1958 of file ddx_core.f90.

◆ tree_l2l_rotation_adj_work()

subroutine ddx_core::tree_l2l_rotation_adj_work ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(inout)  node_l,
real(dp), dimension(6*params % pl**2 + 19*params % pl + 8), intent(out)  work 
)

Adjoint transfer local coefficients over a tree.

Definition at line 1974 of file ddx_core.f90.

◆ tree_l2l_bessel_rotation_adj()

subroutine ddx_core::tree_l2l_bessel_rotation_adj ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(inout)  node_l 
)

Adjoint transfer local coefficients over a tree.

Definition at line 2015 of file ddx_core.f90.

◆ tree_m2l_rotation()

subroutine ddx_core::tree_m2l_rotation ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pm+1)**2, constants % nclusters), intent(in)  node_m,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(out)  node_l 
)

Transfer multipole local coefficients into local over a tree.

Definition at line 2045 of file ddx_core.f90.

◆ tree_m2l_bessel_rotation()

subroutine ddx_core::tree_m2l_bessel_rotation ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pm+1)**2, constants % nclusters), intent(in)  node_m,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(out)  node_l 
)

Transfer multipole local coefficients into local over a tree.

Definition at line 2094 of file ddx_core.f90.

◆ tree_m2l_bessel_rotation_adj()

subroutine ddx_core::tree_m2l_bessel_rotation_adj ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(in)  node_l,
real(dp), dimension((params % pm+1)**2, constants % nclusters), intent(inout)  node_m 
)

Adjoint transfer multipole local coefficients into local over a tree.

Definition at line 2145 of file ddx_core.f90.

◆ tree_m2l_bessel_rotation_adj_work()

subroutine ddx_core::tree_m2l_bessel_rotation_adj_work ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(in)  node_l,
real(dp), dimension((params % pm+1)**2, constants % nclusters), intent(out)  node_m 
)

Adjoint transfer multipole local coefficients into local over a tree.

Definition at line 2161 of file ddx_core.f90.

◆ tree_m2l_rotation_adj()

subroutine ddx_core::tree_m2l_rotation_adj ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(in)  node_l,
real(dp), dimension((params % pm+1)**2, constants % nclusters), intent(out)  node_m 
)

Adjoint transfer multipole local coefficients into local over a tree.

Definition at line 2207 of file ddx_core.f90.

◆ tree_l2p()

subroutine ddx_core::tree_l2p ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), intent(in)  alpha,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(in)  node_l,
real(dp), intent(in)  beta,
real(dp), dimension(params % ngrid, params % nsph), intent(inout)  grid_v,
real(dp), dimension((params % pl+1)**2, params % nsph), intent(out)  sph_l 
)

TODO.

Definition at line 2251 of file ddx_core.f90.

◆ tree_l2p_bessel()

subroutine ddx_core::tree_l2p_bessel ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), intent(in)  alpha,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(in)  node_l,
real(dp), intent(in)  beta,
real(dp), dimension(params % ngrid, params % nsph), intent(inout)  grid_v 
)

TODO.

Definition at line 2290 of file ddx_core.f90.

◆ tree_l2p_adj()

subroutine ddx_core::tree_l2p_adj ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), intent(in)  alpha,
real(dp), dimension(params % ngrid, params % nsph), intent(in)  grid_v,
real(dp), intent(in)  beta,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(inout)  node_l,
real(dp), dimension((params % pl+1)**2, params % nsph), intent(out)  sph_l 
)

TODO.

Definition at line 2325 of file ddx_core.f90.

◆ tree_l2p_bessel_adj()

subroutine ddx_core::tree_l2p_bessel_adj ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), intent(in)  alpha,
real(dp), dimension(params % ngrid, params % nsph), intent(in)  grid_v,
real(dp), intent(in)  beta,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(inout)  node_l 
)

TODO.

Definition at line 2363 of file ddx_core.f90.

◆ tree_m2p()

subroutine ddx_core::tree_m2p ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
integer, intent(in)  p,
real(dp), intent(in)  alpha,
real(dp), dimension((p+1)**2, params % nsph), intent(in)  sph_m,
real(dp), intent(in)  beta,
real(dp), dimension(params % ngrid, params % nsph), intent(inout)  grid_v 
)

TODO.

Definition at line 2400 of file ddx_core.f90.

◆ tree_m2p_bessel()

subroutine ddx_core::tree_m2p_bessel ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
integer, intent(in)  p,
real(dp), intent(in)  alpha,
integer, intent(in)  sph_p,
real(dp), dimension((sph_p+1)**2, params % nsph), intent(in)  sph_m,
real(dp), intent(in)  beta,
real(dp), dimension(params % ngrid, params % nsph), intent(inout)  grid_v 
)

TODO.

Definition at line 2450 of file ddx_core.f90.

◆ tree_m2p_adj()

subroutine ddx_core::tree_m2p_adj ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
integer, intent(in)  p,
real(dp), intent(in)  alpha,
real(dp), dimension(params % ngrid, params % nsph), intent(in)  grid_v,
real(dp), intent(in)  beta,
real(dp), dimension((p+1)**2, params % nsph), intent(inout)  sph_m 
)

TODO.

Definition at line 2502 of file ddx_core.f90.

◆ tree_m2p_bessel_adj()

subroutine ddx_core::tree_m2p_bessel_adj ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
integer, intent(in)  p,
real(dp), intent(in)  alpha,
real(dp), dimension(params % ngrid, params % nsph), intent(in)  grid_v,
real(dp), intent(in)  beta,
integer, intent(in)  sph_p,
real(dp), dimension((sph_p+1)**2, params % nsph), intent(inout)  sph_m 
)

TODO.

Definition at line 2553 of file ddx_core.f90.

◆ tree_m2p_bessel_nodiag_adj()

subroutine ddx_core::tree_m2p_bessel_nodiag_adj ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
integer, intent(in)  p,
real(dp), intent(in)  alpha,
real(dp), dimension(params % ngrid, params % nsph), intent(in)  grid_v,
real(dp), intent(in)  beta,
integer, intent(in)  sph_p,
real(dp), dimension((sph_p+1)**2, params % nsph), intent(inout)  sph_m 
)

TODO.

Definition at line 2602 of file ddx_core.f90.

◆ tree_grad_m2m()

subroutine ddx_core::tree_grad_m2m ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension(constants % nbasis, params % nsph), intent(in)  sph_m,
real(dp), dimension((params % lmax+2)**2, 3, params % nsph), intent(inout)  sph_m_grad,
real(dp), dimension((params % lmax+2)**2, params % nsph), intent(inout)  work 
)

TODO.

Definition at line 2649 of file ddx_core.f90.

◆ tree_grad_l2l()

subroutine ddx_core::tree_grad_l2l ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), dimension((params % pl+1)**2, constants % nclusters), intent(in)  node_l,
real(dp), dimension((params % pl+1)**2, 3, params % nsph), intent(out)  sph_l_grad,
real(dp), dimension((params % pl+1)**2, params % nsph), intent(out)  work 
)

TODO.

Definition at line 2713 of file ddx_core.f90.

◆ cav_to_spherical()

subroutine ddx_core::cav_to_spherical ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
type(ddx_workspace_type), intent(inout)  workspace,
real(dp), dimension(constants % ncav), intent(in)  property_cav,
real(dp), dimension(constants % nbasis, params % nsph), intent(out)  property_sph 
)

Transform a function defined at the exposed cavity points (cav) to a spherical harmonics expansion. Note that the function is also multiplied by the characteristic function U.

Parameters
[in]paramsddx parameters
[in]constantsddx constants
[in,out]workspaceddx workspace
[in]property_cavproperty defined at the exposed cavity points, size (ncav)
[out]property_sphproperty as a spherical harmonics expansion, size (nbasis, nsph)

Definition at line 2836 of file ddx_core.f90.