ddx 0.6.0
Libary for domain-decomposition methods for polarizable continuum models
Fortran interface: core routines

Functions/Subroutines

subroutine ddx::ddfromfile (fname, ddx_data, tol, charges, ddx_error)
 Read the configuration from ddX input file and return a ddx_data structure. More...
 
subroutine ddx::ddinit (model, nsph, coords, radii, eps, ddx_data, ddx_error, force, kappa, eta, shift, lmax, ngrid, incore, maxiter, jacobi_ndiis, enable_fmm, pm, pl, nproc, logfile, adjoint, eps_int)
 Wrapper to the initialization routine which supports optional arguments. This will make future maintenance easier. More...
 
subroutine ddx::ddrun (ddx_data, state, electrostatics, psi, tol, esolv, ddx_error, force, read_guess)
 Main solver routine. More...
 
subroutine ddx::setup (params, constants, workspace, state, electrostatics, psi, ddx_error)
 Setup the state for the different models. More...
 
subroutine ddx::fill_guess (params, constants, workspace, state, tol, ddx_error)
 Do a guess for the primal linear system for the different models. More...
 
subroutine ddx::fill_guess_adjoint (params, constants, workspace, state, tol, ddx_error)
 Do a guess for the adjoint linear system for the different models. More...
 
subroutine ddx::solve (params, constants, workspace, state, tol, ddx_error)
 Solve the primal linear system for the different models. More...
 
subroutine ddx::solve_adjoint (params, constants, workspace, state, tol, ddx_error)
 Solve the adjoint linear system for the different models. More...
 
subroutine ddx::energy (params, constants, workspace, state, solvation_energy, ddx_error)
 Compute the energy for the different models. More...
 
subroutine ddx::solvation_force_terms (params, constants, workspace, state, electrostatics, force, ddx_error)
 Compute the solvation terms of the forces (solute aspecific) for the different models. This must be summed to the solute specific term to get the full forces. More...
 
subroutine ddx_core::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 ddx_core::deallocate_electrostatics (electrostatics, ddx_error)
 Deallocate the electrostatic properties. More...
 
subroutine ddx_core::allocate_state (params, constants, state, ddx_error)
 Initialize the ddx_state object. More...
 
subroutine ddx_core::deallocate_state (state, ddx_error)
 Deallocate the ddx_state object. More...
 
subroutine ddx_core::get_banner (string)
 Print the ddX logo. More...
 

Detailed Description

Function/Subroutine Documentation

◆ ddfromfile()

subroutine ddx::ddfromfile ( character(len=*), intent(in)  fname,
type(ddx_type), intent(out)  ddx_data,
real(dp), intent(out)  tol,
real(dp), dimension(:), intent(out), allocatable  charges,
type(ddx_error_type), intent(inout)  ddx_error 
)

Read the configuration from ddX input file and return a ddx_data structure.

Parameters
[in]fnameFilename containing all the required info
[out]ddx_dataObject containing all inputs
[out]toltolerance for iterative solvers
[out]chargescharge array, size(nsph)
[in,out]ddx_errorddX error

Definition at line 25 of file ddx.f90.

◆ ddinit()

subroutine ddx::ddinit ( integer, intent(in)  model,
integer, intent(in)  nsph,
real(dp), dimension(3, nsph), intent(in)  coords,
real(dp), dimension(nsph), intent(in)  radii,
real(dp), intent(in)  eps,
type(ddx_type), intent(out), target  ddx_data,
type(ddx_error_type), intent(inout)  ddx_error,
integer, intent(in), optional  force,
real(dp), intent(in), optional  kappa,
real(dp), intent(in), optional  eta,
real(dp), intent(in), optional  shift,
integer, intent(in), optional  lmax,
integer, intent(in), optional  ngrid,
integer, intent(in), optional  incore,
integer, intent(in), optional  maxiter,
integer, intent(in), optional  jacobi_ndiis,
integer, intent(in), optional  enable_fmm,
integer, intent(in), optional  pm,
integer, intent(in), optional  pl,
integer, intent(in), optional  nproc,
character(len=255), intent(in), optional  logfile,
integer, intent(in), optional  adjoint,
real(dp), intent(in), optional  eps_int 
)

Wrapper to the initialization routine which supports optional arguments. This will make future maintenance easier.

Parameters
[in]model1 for COSMO, 2 for PCM and 3 for LPB
[in]nsphnumber of spheres n > 0
[in]coordscoordinates of the spheres, size (3, nsph)
[in]radiiradii of the spheres, size (nsph)
[in]epsrelative dielectric permittivity eps > 1
[out]ddx_dataContainer for ddX data structures
[in,out]ddx_errorddX error
[in,optional]force: 1 if forces are required and 0 otherwise
[in,optional]kappa: Debye-H"{u}ckel parameter
[in,optional]eta: regularization parameter 0 < eta <= 1
[in,optional]shift: shift of characteristic function -1 for interior, 0 for centered and 1 for outer regularization
[in,optional]lmax: Maximal degree of modeling spherical harmonics, lmax >= 0
[in,optional]ngrid: Number of Lebedev grid points ngrid >= 0
[in,optional]incore: handling of the sparse matrices, 1 for precomputing them and keeping them in memory, 0 for assembling the matrix-vector products on-the-fly
[in,optional]maxiter: Maximum number of iterations for the iterative solver, maxiter > 0
[in,optional]jacobi_ndiis: Number of extrapolation points for the Jacobi/DIIS solver, ndiis >= 0
[in,optional]enable_fmm: 1 to use FMM acceleration and 0 otherwise
[in,optional]pm: Maximal degree of multipole spherical harmonics. Ignored in the case fmm=0. Value -1 means no far-field FMM interactions are computed, pm >= -1
[in,optional]pl: Maximal degree of local spherical harmonics. Ignored in the case fmm=0. Value -1 means no far-field FMM interactions are computed, pl >= -1
[in,optional]nproc: Number of OpenMP threads, nproc >= 0.
[in,optional]logfile: file name for log information.

Definition at line 243 of file ddx.f90.

◆ ddrun()

subroutine ddx::ddrun ( type(ddx_type), intent(inout)  ddx_data,
type(ddx_state_type), intent(inout)  state,
type(ddx_electrostatics_type), intent(in)  electrostatics,
real(dp), dimension(ddx_data % constants % nbasis, ddx_data % params % nsph), intent(in)  psi,
real(dp), intent(in)  tol,
real(dp), intent(out)  esolv,
type(ddx_error_type), intent(inout)  ddx_error,
real(dp), dimension(3, ddx_data % params % nsph), intent(out), optional  force,
logical, intent(in), optional  read_guess 
)

Main solver routine.

Solves the solvation problem, computes the energy, and if required computes the forces.

Parameters
[in]ddx_dataddX object with all input information
[in,out]stateddx state (contains RHSs and solutions)
[in]electrostaticselectrostatic property container
[in]psiRHS of the adjoint problem
[in]toltolerance for the linear system solvers
[out]esolvsolvation energy
[in,out]ddx_errorddX error
[out]forceAnalytical forces (optional argument, only if required)
[in]read_guessoptional argument, if true read the guess from the state object

Definition at line 362 of file ddx.f90.

◆ setup()

subroutine ddx::setup ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(inout)  constants,
type(ddx_workspace_type), intent(inout)  workspace,
type(ddx_state_type), intent(inout)  state,
type(ddx_electrostatics_type), intent(in)  electrostatics,
real(dp), dimension(constants % nbasis, params % nsph), intent(in)  psi,
type(ddx_error_type), intent(inout)  ddx_error 
)

Setup the state for the different models.

Parameters
[in]paramsUser specified parameters
[in]constantsPrecomputed constants
[in,out]workspacePreallocated workspaces
[in,out]stateddx state (contains solutions and RHSs)
[in]electrostaticselectrostatic property container
[in]psiRepresentation of the solute potential in spherical harmonics, size (nbasis, nsph)
[in,out]ddx_errorddX error

Definition at line 471 of file ddx.f90.

◆ fill_guess()

subroutine ddx::fill_guess ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(inout)  constants,
type(ddx_workspace_type), intent(inout)  workspace,
type(ddx_state_type), intent(inout)  state,
real(dp), intent(in)  tol,
type(ddx_error_type), intent(inout)  ddx_error 
)

Do a guess for the primal linear system for the different models.

Parameters
[in]paramsUser specified parameters
[in,out]constantsPrecomputed constants
[in,out]workspacePreallocated workspaces
[in,out]stateddx state (contains solutions and RHSs)
[in]toltolerance
[in,out]ddx_errorddX error

Definition at line 509 of file ddx.f90.

◆ fill_guess_adjoint()

subroutine ddx::fill_guess_adjoint ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(inout)  constants,
type(ddx_workspace_type), intent(inout)  workspace,
type(ddx_state_type), intent(inout)  state,
real(dp), intent(in)  tol,
type(ddx_error_type), intent(inout)  ddx_error 
)

Do a guess for the adjoint linear system for the different models.

Parameters
[in]paramsUser specified parameters
[in,out]constantsPrecomputed constants
[in,out]workspacePreallocated workspaces
[in,out]stateddx state (contains solutions and RHSs)
[in]toltolerance
[in,out]ddx_errorddX error

Definition at line 541 of file ddx.f90.

◆ solve()

subroutine ddx::solve ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(inout)  constants,
type(ddx_workspace_type), intent(inout)  workspace,
type(ddx_state_type), intent(inout)  state,
real(dp), intent(in)  tol,
type(ddx_error_type), intent(inout)  ddx_error 
)

Solve the primal linear system for the different models.

Parameters
[in]params: General options
[in]constants: Precomputed constants
[in,out]workspace: Preallocated workspaces
[in,out]state: Solutions, guesses and relevant quantities
[in]tol: Tolerance for the iterative solvers
[in,out]ddx_errorddX error

Definition at line 573 of file ddx.f90.

◆ solve_adjoint()

subroutine ddx::solve_adjoint ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(inout)  constants,
type(ddx_workspace_type), intent(inout)  workspace,
type(ddx_state_type), intent(inout)  state,
real(dp), intent(in)  tol,
type(ddx_error_type), intent(inout)  ddx_error 
)

Solve the adjoint linear system for the different models.

Parameters
[in]params: General options
[in]constants: Precomputed constants
[in,out]workspace: Preallocated workspaces
[in,out]state: Solutions, guesses and relevant quantities
[in]tol: Tolerance for the iterative solvers
[in,out]ddx_errorddX error

Definition at line 605 of file ddx.f90.

◆ energy()

subroutine ddx::energy ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
type(ddx_workspace_type), intent(in)  workspace,
type(ddx_state_type), intent(in)  state,
real(dp), intent(out)  solvation_energy,
type(ddx_error_type), intent(inout)  ddx_error 
)

Compute the energy for the different models.

Parameters
[in]paramsGeneral options
[in]constantsPrecomputed constants
[in,out]workspacePreallocated workspaces
[in]stateddx state (contains solutions and RHSs)
[out]solvation_energyresulting energy
[in,out]ddx_errorddX error

Definition at line 637 of file ddx.f90.

◆ solvation_force_terms()

subroutine ddx::solvation_force_terms ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
type(ddx_workspace_type), intent(inout)  workspace,
type(ddx_state_type), intent(inout)  state,
type(ddx_electrostatics_type), intent(in)  electrostatics,
real(dp), dimension(3, params % nsph), intent(out)  force,
type(ddx_error_type), intent(inout)  ddx_error 
)

Compute the solvation terms of the forces (solute aspecific) for the different models. This must be summed to the solute specific term to get the full forces.

Parameters
[in]paramsGeneral options
[in]constantsPrecomputed constants
[in,out]workspacePreallocated workspaces
[in,out]stateSolutions and relevant quantities
[in]electrostaticsElectrostatic properties container.
[out]forceGeometrical contribution to the forces
[in,out]ddx_errorddX error

Definition at line 675 of file ddx.f90.

◆ allocate_electrostatics()

subroutine ddx_core::allocate_electrostatics ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
type(ddx_electrostatics_type), intent(out)  electrostatics,
type(ddx_error_type), intent(inout)  ddx_error 
)

Given the chosen model, find the required electrostatic properties, and allocate the arrays for them in the container.

Parameters
[in]paramsddx parameters
[in]constantsddx constants
[out]electrostatics,ddXelectrostatic properties container
[in,out]ddx_errorddX error

Definition at line 340 of file ddx_core.f90.

◆ deallocate_electrostatics()

subroutine ddx_core::deallocate_electrostatics ( type(ddx_electrostatics_type), intent(inout)  electrostatics,
type(ddx_error_type), intent(inout)  ddx_error 
)

Deallocate the electrostatic properties.

Parameters
[out]electrostatics,ddXelectrostatic properties container
[in,out]ddx_errorddX error

Definition at line 406 of file ddx_core.f90.

◆ allocate_state()

subroutine ddx_core::allocate_state ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
type(ddx_state_type), intent(out)  state,
type(ddx_error_type), intent(inout)  ddx_error 
)

Initialize the ddx_state object.

Parameters
[in]paramsUser specified parameters
[in]constantsPrecomputed constants
[in,out]stateddx state (contains solutions and RHSs)
[in,out]ddx_errorddX error

Definition at line 446 of file ddx_core.f90.

◆ deallocate_state()

subroutine ddx_core::deallocate_state ( type(ddx_state_type), intent(inout)  state,
type(ddx_error_type), intent(inout)  ddx_error 
)

Deallocate the ddx_state object.

Parameters
[in,out]stateddx state (contains solutions and RHSs)
[in,out]ddx_errorddX error

Definition at line 774 of file ddx_core.f90.

◆ get_banner()

subroutine ddx_core::get_banner ( character (len=2047), intent(out)  string)

Print the ddX logo.

Parameters
[out]stringcontainer for the logo

Definition at line 2788 of file ddx_core.f90.