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

Run-time constants. More...

Data Types

type  ddx_constants_type
 Container for precomputed constants. More...
 

Functions/Subroutines

subroutine constants_init (params, constants, ddx_error)
 Compute all necessary constants. More...
 
subroutine build_itrnl (constants, params, ddx_error)
 Build the transposed neighbor list. More...
 
subroutine build_l (constants, params, ddx_error)
 Allocate and build the ddCOSMO sparse matrix, only if incore is set. More...
 
subroutine build_b (constants, params, ddx_error)
 Allocate and build the HSP sparse matrix, only if incore is set. More...
 
subroutine mkpmat (params, constants, isph, pmat)
 Computation of P_chi. More...
 
subroutine constants_geometry_init (params, constants, ddx_error)
 Initialize geometry-related constants like list of neighbouring spheres. More...
 
subroutine neighbor_list_init (params, constants, ddx_error)
 Build the neighbor list using a N^2 code. More...
 
subroutine neighbor_list_init_fmm (params, constants, ddx_error)
 Build the neighbor list using a linear scaling code (only if FMMs are enabled) More...
 
real(dp) function fsw (t, se, eta)
 Switching function. More...
 
real(dp) function dfsw (t, se, eta)
 Derivative of a switching function. More...
 
subroutine mkprec (lmax, nbasis, nsph, ngrid, eps, ui, wgrid, vgrid, vgrid_nbasis, rx_prc, ddx_error)
 Compute preconditioner assemble the diagonal blocks of the reps matrix then invert them to build the preconditioner. More...
 
subroutine tree_rib_build (nsph, csph, rsph, order, cluster, children, parent, cnode, rnode, snode, ddx_error)
 Build a recursive inertial binary tree. More...
 
subroutine tree_rib_node_bisect (nsph, csph, n, order, div, ddx_error)
 Divide given cluster of spheres into two subclusters by inertial bisection. More...
 
subroutine tree_get_farnear_work (n, children, cnode, rnode, lwork, iwork, jwork, work, nnfar, nfar, nnnear, nnear)
 Find near and far admissible pairs of tree nodes and store it in work array. More...
 
subroutine tree_get_farnear (jwork, lwork, work, n, nnfar, nfar, sfar, far, nnnear, nnear, snear, near)
 Get near and far admissible pairs from work array of tree_get_farnear_work Works only for binary tree. More...
 
subroutine constants_free (constants, ddx_error)
 Deallocate the constants. More...
 

Detailed Description

Run-time constants.

Function/Subroutine Documentation

◆ constants_init()

subroutine ddx_constants::constants_init ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(out)  constants,
type(ddx_error_type), intent(inout)  ddx_error 
)

Compute all necessary constants.

Parameters
[in]paramsObject containing all inputs.
[out]constantsObject containing all constants.
[in,out]ddx_errorddX error

Definition at line 236 of file ddx_constants.f90.

◆ build_itrnl()

subroutine ddx_constants::build_itrnl ( type(ddx_constants_type), intent(inout)  constants,
type(ddx_params_type), intent(in)  params,
type(ddx_error_type), intent(inout)  ddx_error 
)

Build the transposed neighbor list.

Definition at line 555 of file ddx_constants.f90.

◆ build_l()

subroutine ddx_constants::build_l ( type(ddx_constants_type), intent(inout)  constants,
type(ddx_params_type), intent(in)  params,
type(ddx_error_type), intent(inout)  ddx_error 
)

Allocate and build the ddCOSMO sparse matrix, only if incore is set.

Definition at line 580 of file ddx_constants.f90.

◆ build_b()

subroutine ddx_constants::build_b ( type(ddx_constants_type), intent(inout)  constants,
type(ddx_params_type), intent(in)  params,
type(ddx_error_type), intent(inout)  ddx_error 
)

Allocate and build the HSP sparse matrix, only if incore is set.

Definition at line 652 of file ddx_constants.f90.

◆ mkpmat()

subroutine ddx_constants::mkpmat ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
integer, intent(in)  isph,
real(dp), dimension(constants % nbasis, constants % nbasis0), intent(inout)  pmat 
)

Computation of P_chi.

Definition at line 731 of file ddx_constants.f90.

◆ constants_geometry_init()

subroutine ddx_constants::constants_geometry_init ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(inout)  constants,
type(ddx_error_type), intent(inout)  ddx_error 
)

Initialize geometry-related constants like list of neighbouring spheres.

This routine does not check ddx_error state of params or constants as it is intended for a low-level use.

Parameters
[in]paramsObject containing all inputs.
[in,out]constantsObject containing all constants.
[in,out]ddx_errorddX error

Definition at line 767 of file ddx_constants.f90.

◆ neighbor_list_init()

subroutine ddx_constants::neighbor_list_init ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(inout)  constants,
type(ddx_error_type), intent(inout)  ddx_error 
)

Build the neighbor list using a N^2 code.

Definition at line 1074 of file ddx_constants.f90.

◆ neighbor_list_init_fmm()

subroutine ddx_constants::neighbor_list_init_fmm ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(inout)  constants,
type(ddx_error_type), intent(inout)  ddx_error 
)

Build the neighbor list using a linear scaling code (only if FMMs are enabled)

Definition at line 1159 of file ddx_constants.f90.

◆ fsw()

real(dp) function ddx_constants::fsw ( real(dp), intent(in)  t,
real(dp), intent(in)  se,
real(dp), intent(in)  eta 
)

Switching function.

This is an implementation of \( \chi(t) \) with a shift \( se \):

\[ \chi_\eta(t) = \left\{ \begin{array}{ll} 0, & \text{if} \quad x \geq 1 \\ p_\eta(x), & \text{if} \quad 1-\eta < x < 1 \\ 1, & \text{if} \quad x \leq 1-\eta \end{array} \right. \]

where \( x = t - \frac{1+se}{2} \eta \) is a shifted coordinate and

\[ p_\eta(x) = \frac{1}{\eta^5} (1-t)^3 (6t^2 + (15\eta-12)t + (10\eta^2 -15\eta+6)) \]

is a smoothing polynomial of the 5th degree. In the case shift \( se=1 \) the switch function \( \chi_\eta(t) \) is 1 for \( t \in [0,1] \), is 0 for \( t \in [1+\eta, \infty) \) and varies in \( [1, 1+\eta] \) which is an external shift. In the case shift \( se=-1 \) the switch function \( \chi_\eta(t) \) is 1 for \( t \in [0,1-\eta] \), is 0 for \( t \in [1, \infty) \) and varies in \( [1-\eta, 1] \) which is an internal shift. In the case shift \( se=0 \) the switch function \( \chi_\eta(t) \) is 1 for \( t \in [0,1-\eta/2] \), is 0 for \( t \in [1+\eta/2, \infty) \) and varies in \( [1-\eta/2, 1+\eta/2] \) which is a centered shift.

Parameters
[in]treal non-negative input value
[in]seshift
[in]etaregularization parameter \( \eta \)

Definition at line 1276 of file ddx_constants.f90.

◆ dfsw()

real(dp) function ddx_constants::dfsw ( real(dp), intent(in)  t,
real(dp), intent(in)  se,
real(dp), intent(in)  eta 
)

Derivative of a switching function.

This is an implementation of \( \chi'(t) \) with a shift \( se \):

\[ \chi_\eta'(t) = \left\{ \begin{array}{ll} 0, & \text{if} \quad x \geq 1 \\ p'_\eta(x), & \text{if} \quad 1-\eta < x < 1 \\ 0, & \text{if} \quad x \leq 1-\eta \end{array} \right. \]

where \( x = t - \frac{1+se}{2} \eta \) is a shifted coordinate and

\[ p'_\eta(x) = -\frac{30}{\eta^5} (1-t)^2 (t-1+\eta)^2 \]

is a derivative of the smoothing polynomial. In the case shift \( se=1 \) the derivative \( \chi_\eta'(t) \) is 0 for \( t \in [0,1] \cup [1+\eta, \infty) \) and varies in \( [1, 1+\eta] \) which is an external shift. In the case shift \( se=-1 \) the derivative \( \chi_\eta'(t) \) is 0 for \( t \in [0,1-\eta] \cup [1, \infty) \) and varies in \( [1-\eta, 1] \) which is an internal shift. In the case shift \( se=0 \) the derivative \( \chi_\eta'(t) \) is 0 for \( t \in [0,1-\eta/2] \cup [1+\eta/2, \infty) \) and varies in \( [1-\eta/2, 1+\eta/2] \) which is a centered shift.

Parameters
[in]treal non-negative input value
[in]seshift
[in]etaregularization parameter \( \eta \)

Definition at line 1331 of file ddx_constants.f90.

◆ mkprec()

subroutine ddx_constants::mkprec ( integer, intent(in)  lmax,
integer, intent(in)  nbasis,
integer, intent(in)  nsph,
integer, intent(in)  ngrid,
real(dp), intent(in)  eps,
real(dp), dimension(ngrid, nsph), intent(in)  ui,
real(dp), dimension(ngrid), intent(in)  wgrid,
real(dp), dimension(vgrid_nbasis, ngrid), intent(in)  vgrid,
integer, intent(in)  vgrid_nbasis,
real(dp), dimension(nbasis, nbasis, nsph), intent(out)  rx_prc,
type(ddx_error_type), intent(inout)  ddx_error 
)

Compute preconditioner assemble the diagonal blocks of the reps matrix then invert them to build the preconditioner.

Definition at line 1357 of file ddx_constants.f90.

◆ tree_rib_build()

subroutine ddx_constants::tree_rib_build ( integer, intent(in)  nsph,
real(dp), dimension(3, nsph), intent(in)  csph,
real(dp), dimension(nsph), intent(in)  rsph,
integer, dimension(nsph), intent(out)  order,
integer, dimension(2, 2*nsph-1), intent(out)  cluster,
integer, dimension(2, 2*nsph-1), intent(out)  children,
integer, dimension(2*nsph-1), intent(out)  parent,
real(dp), dimension(3, 2*nsph-1), intent(out)  cnode,
real(dp), dimension(2*nsph-1), intent(out)  rnode,
integer, dimension(nsph), intent(out)  snode,
type(ddx_error_type), intent(inout)  ddx_error 
)

Build a recursive inertial binary tree.

Uses inertial bisection in a recursive manner until each leaf node has only one input sphere inside. Number of tree nodes is always 2*nsph-1.

Parameters
[in]nsphNumber of input spheres
[in]csphCenters of input spheres
[in]rsphRadii of input spheres
[out]orderOrdering of input spheres to make spheres of one cluster contiguous.
[out]clusterIndexes in order array of the first and the last spheres of each cluster
[out]childrenIndexes of the first and the last children nodes. If both indexes are equal this means there is only 1 child node and if value of the first child node is 0, there are no children nodes.
[out]parentParent of each cluster. Value 0 means there is no parent node which corresponds to the root node (with index 1).
[out]cnodeCenter of a bounding sphere of each node
[out]rnodeRadius of a bounding sphere of each node
[out]snodeArray of leaf nodes containing input spheres
[in,out]ddx_errorddX error

Definition at line 1448 of file ddx_constants.f90.

◆ tree_rib_node_bisect()

subroutine ddx_constants::tree_rib_node_bisect ( integer, intent(in)  nsph,
real(dp), dimension(3, nsph), intent(in)  csph,
integer, intent(in)  n,
integer, dimension(n), intent(inout)  order,
integer, intent(out)  div,
type(ddx_error_type), intent(inout)  ddx_error 
)

Divide given cluster of spheres into two subclusters by inertial bisection.

Parameters
[in]nsphNumber of all input spheres
[in]csphCenters of all input spheres
[in]nNumber of spheres in a given cluster
[in,out]orderIndexes of spheres in a given cluster. On exit, indexes order(1:div) correspond to the first subcluster and indexes order(div+1:n) correspond to the second subcluster.
[out]divBreak point of order array between two clusters.
[in,out]ddx_errorddX error

Definition at line 1575 of file ddx_constants.f90.

◆ tree_get_farnear_work()

subroutine ddx_constants::tree_get_farnear_work ( integer, intent(in)  n,
integer, dimension(2, n), intent(in)  children,
real(dp), dimension(3, n), intent(in)  cnode,
real(dp), dimension(n), intent(in)  rnode,
integer, intent(in)  lwork,
integer, intent(inout)  iwork,
integer, intent(inout)  jwork,
integer, dimension(3, lwork), intent(inout)  work,
integer, intent(out)  nnfar,
integer, dimension(n), intent(out)  nfar,
integer, intent(out)  nnnear,
integer, dimension(n), intent(out)  nnear 
)

Find near and far admissible pairs of tree nodes and store it in work array.

Parameters
[in]nnumber of nodes
[in]childrenfirst and last children of each cluster. Value 0 means no children (leaf node).
[in]cnodecenter of bounding sphere of each cluster (node) of tree
[in]rnoderadius of bounding sphere of each cluster (node) of tree
[in]lworksize of work array in dimension 2
[in,out]iworkindex of current pair of nodes that needs to be checked for admissibility. Must be 0 for the first call of this subroutine. If on exit iwork is less or equal to jwork, that means lwork was too small, please reallocate work array and copy all the values into new array and then run procedure again.
[in,out]jworkamount of stored possible admissible pairs of nodes. Please read iwork comments.
[in,out]workall the far and near pairs will be stored here
[out]nnfartotal amount of far admissible pairs. valid only if iwork is greater than jwork on exit.
[out]nfaramount of far admissible pairs for each node. valid only if iwork is greater than jwork on exit.
[out]nnneartotal amount of near admissible pairs. valid only if iwork is greater than jwork on exit
[out]nnearamount of near admissible pairs for each node. valid only if iwork is greater than jwork on exit

Definition at line 1683 of file ddx_constants.f90.

◆ tree_get_farnear()

subroutine ddx_constants::tree_get_farnear ( integer, intent(in)  jwork,
integer, intent(in)  lwork,
integer, dimension(3, lwork), intent(in)  work,
integer, intent(in)  n,
integer, intent(in)  nnfar,
integer, dimension(n), intent(in)  nfar,
integer, dimension(n+1), intent(out)  sfar,
integer, dimension(nnfar), intent(out)  far,
integer, intent(in)  nnnear,
integer, dimension(n), intent(in)  nnear,
integer, dimension(n+1), intent(out)  snear,
integer, dimension(nnnear), intent(out)  near 
)

Get near and far admissible pairs from work array of tree_get_farnear_work Works only for binary tree.

Definition at line 1776 of file ddx_constants.f90.

◆ constants_free()

subroutine ddx_constants::constants_free ( type(ddx_constants_type), intent(out)  constants,
type(ddx_error_type), intent(inout)  ddx_error 
)

Deallocate the constants.

Parameters
[out]constantsPrecomputed constants
[in,out]ddx_errorddX error

Definition at line 1835 of file ddx_constants.f90.