ddx 0.6.0
Libary for domain-decomposition methods for polarizable continuum models
C Interface

Generic quantities

void ddx_get_banner (char *banner, int maxlen)
 
int ddx_supported_lebedev_grids (int maxlen, int *grids)
 
void ddx_scaled_ylm (const void *ddx, int lmax, const double *x, int sphere, double *ylm)
 

Allocate and manage the ddX error object

void * ddx_allocate_error ()
 
int ddx_get_error_flag (const void *error)
 
void ddx_get_error_message (const void *error, char *message, int maxlen)
 

Allocate and manage the ddX electrostatics object

void * ddx_allocate_electrostatics (void *ddx, void *error)
 
void ddx_multipole_electrostatics (void *ddx, int nsph, int nmultipoles, const double *multipoles, void *electrostatics, void *error)
 
void ddx_deallocate_electrostatics (void *electrostatics, void *error)
 

Allocate and manage the ddx model object

void * ddx_allocate_model (int model, int enable_force, double solvent_epsilon, double solvent_kappa, double eta, double shift, int lmax, int n_lebedev, int incore, int maxiter, int jacobi_n_diis, int enable_fmm, int fmm_multipole_lmax, int fmm_local_lmax, int n_proc, int n_spheres, const double *sphere_centres, const double *sphere_radii, int length_logfile, const char *logfile, void *error)
 
void ddx_deallocate_model (void *ddx, void *error)
 

Getters for properties of the state object

void ddx_get_logfile (const void *ddx, char *logfile, int maxlen)
 
int ddx_get_enable_fmm (const void *ddx)
 
int ddx_get_enable_force (const void *ddx)
 
int ddx_get_jacobi_n_diis (const void *ddx)
 
int ddx_get_lmax (const void *ddx)
 
int ddx_get_incore (const void *ddx)
 
int ddx_get_maxiter (const void *ddx)
 
int ddx_get_model (const void *ddx)
 
int ddx_get_n_lebedev (const void *ddx)
 
int ddx_get_n_proc (const void *ddx)
 
int ddx_get_n_spheres (const void *ddx)
 
int ddx_get_fmm_local_lmax (const void *ddx)
 
int ddx_get_fmm_multipole_lmax (const void *ddx)
 
double ddx_get_solvent_epsilon (const void *ddx)
 
double ddx_get_eta (const void *ddx)
 
double ddx_get_solvent_kappa (const void *ddx)
 
double ddx_get_shift (const void *ddx)
 
void ddx_get_sphere_charges (const void *ddx, int nsph, double *c_charge)
 
void ddx_get_sphere_centres (const void *ddx, int nsph, double *c_csph)
 
void ddx_get_sphere_radii (const void *ddx, int nsph, double *c_rsph)
 
int ddx_get_n_basis (const void *ddx)
 
int ddx_get_n_cav (const void *ddx)
 
void ddx_get_cavity (const void *ddx, int ncav, double *c_ccav)
 

Allocate and manage the state object

void * ddx_allocate_state (const void *ddx, void *error)
 
void ddx_deallocate_state (void *state, void *error)
 
void ddx_get_x (const void *state, int nbasis, int nsph, double *x)
 
int ddx_get_x_niter (const void *state)
 
void ddx_get_s (const void *state, int nbasis, int nsph, double *s)
 
int ddx_get_s_niter (const void *state)
 
void ddx_get_xi (const void *state, const void *ddx, int ncav, double *xi)
 
void ddx_get_zeta_dip (const void *state, const void *ddx, int ncav, double *zeta_dip)
 

Model nonspecific setup and solution routines

double ddx_ddsolve (const void *ddx, void *state, const void *electrostatics, int nbasis, int nsph, const double *psi, double tol, double *forces, const int read_guess, void *error)
 
void ddx_setup (const void *ddx, void *state, const void *electrostatics, int nbasis, int nsph, const double *psi, void *error)
 
void ddx_fill_guess (const void *ddx, void *state, double tol, void *error)
 
void ddx_fill_guess_adjoint (const void *ddx, void *state, double tol, void *error)
 
void ddx_solve (const void *ddx, void *state, double tol, void *error)
 
void ddx_solve_adjoint (const void *ddx, void *state, double tol, void *error)
 
double ddx_energy (const void *ddx, void *state, void *error)
 
void ddx_solvation_force_terms (const void *ddx, void *state, void *electrostatics, int nsph, double *forces, void *error)
 
double ddx_ddrun (const void *ddx, void *state, void *electrostatics, int nbasis, int nsph, double *psi, const double tol, double *forces, int read_guess, void *error)
 

Problem setup and solution routines

void ddx_cosmo_setup (const void *ddx, void *state, int ncav, int nbasis, int nsph, const double *psi, const double *phi_cav, void *error)
 
void ddx_cosmo_guess (const void *ddx, void *state, void *error)
 
void ddx_cosmo_guess_adjoint (const void *ddx, void *state, void *error)
 
void ddx_cosmo_solve (const void *ddx, void *state, double tol, void *error)
 
void ddx_cosmo_solve_adjoint (const void *ddx, void *state, double tol, void *error)
 
double ddx_cosmo_energy (const void *ddx, void *state, void *error)
 
void ddx_cosmo_solvation_force_terms (const void *ddx, void *state, int nsph, int ncav, const double *e_cav, double *forces, void *error)
 
void ddx_pcm_setup (const void *ddx, void *state, int ncav, int nbasis, int nsph, const double *psi, const double *phi_cav, void *error)
 
void ddx_pcm_guess (const void *ddx, void *state, void *error)
 
void ddx_pcm_guess_adjoint (const void *ddx, void *state, void *error)
 
void ddx_pcm_solve (const void *ddx, void *state, double tol, void *error)
 
void ddx_pcm_solve_adjoint (const void *ddx, void *state, double tol, void *error)
 
double ddx_pcm_energy (const void *ddx, void *state, void *error)
 
void ddx_pcm_solvation_force_terms (const void *ddx, void *state, int nsph, int ncav, const double *e_cav, double *forces, void *error)
 
void ddx_lpb_setup (const void *ddx, void *state, int ncav, int nbasis, int nsph, const double *psi, const double *phi_cav, const double *e_cav, void *error)
 
void ddx_lpb_guess (const void *ddx, void *state, double tol, void *error)
 
void ddx_lpb_guess_adjoint (const void *ddx, void *state, double tol, void *error)
 
void ddx_lpb_solve (const void *ddx, void *state, double tol, void *error)
 
void ddx_lpb_solve_adjoint (const void *ddx, void *state, double tol, void *error)
 
double ddx_lpb_energy (const void *ddx, void *state, void *error)
 
void ddx_lpb_solvation_force_terms (const void *ddx, void *state, int nsph, int ncav, const double *g_cav, double *forces, const void *error)
 

Multipolar solutes

void ddx_multipole_electrostatics_0 (const void *ddx, int nsph, int ncav, int nmultipoles, const double *multipoles, double *phi_cav, void *error)
 
void ddx_multipole_electrostatics_1 (const void *ddx, int nsph, int ncav, int nmultipoles, const double *multipoles, double *phi_cav, double *e_cav, void *error)
 
void ddx_multipole_electrostatics_2 (const void *ddx, int nsph, int ncav, int nmultipoles, const double *multipoles, double *phi_cav, double *e_cav, double *g_cav, void *error)
 
void ddx_multipole_psi (const void *ddx, int nbasis, int nsph, int nmultipoles, const double *multipoles, double *psi, void *error)
 
void ddx_multipole_force_terms (const void *ddx, void *state, int nsph, int nmultipoles, const double *multipoles, double *forces, void *error)
 

Detailed Description

Header file defining the C interface of ddx.

Function Documentation

◆ ddx_get_banner()

void ddx_get_banner ( char *  banner,
int  maxlen 
)

Fill the string banner with an informative header about ddx. At most maxlen characters will be copied, so use this to specify the number of characters the char array can hold.

◆ ddx_supported_lebedev_grids()

int ddx_supported_lebedev_grids ( int  maxlen,
int *  grids 
)

Fill the grids array with the number of grid points of the supported lebedev grids. At most maxlen elements will be touched, so pass the size of the allocated integer array as maxlen. Returns the number of entries in grids actually altered. E.g. if the returned value is 5 the valid entries in the grid array after the call are grid[0], ..., grid[4].

◆ ddx_scaled_ylm()

void ddx_scaled_ylm ( const void *  ddx,
int  lmax,
const double *  x,
int  sphere,
double *  ylm 
)

With reference to a atomic sphere sphere of radius r centred at a compute:

4π |x - a|^l
------ ----------- Y_l^m(|x - a|)
2l + 1 r^l

where lmax should be identical to the value stored inside ddx or less.

◆ ddx_allocate_error()

void * ddx_allocate_error ( )

Allocate the ddX error.

◆ ddx_get_error_flag()

int ddx_get_error_flag ( const void *  error)

Get the error flag stored inside the model (0 is no error). It is the responsibility of the user of DDX to regularly check for a non-zero error flag and take appropriate actions if an error is present.

◆ ddx_get_error_message()

void ddx_get_error_message ( const void *  error,
char *  message,
int  maxlen 
)

Store the error message corresponding to the most recent error inside the passed char array. at most maxlen characters are copied (so set this to the number of allocated characters in the passed array to avoid a buffer overflow.

◆ ddx_allocate_electrostatics()

void * ddx_allocate_electrostatics ( void *  ddx,
void *  error 
)

Allocate the ddX electrostatics object.

◆ ddx_multipole_electrostatics()

void ddx_multipole_electrostatics ( void *  ddx,
int  nsph,
int  nmultipoles,
const double *  multipoles,
void *  electrostatics,
void *  error 
)

Compute the electrostatic properties for a multipolar distribution

◆ ddx_deallocate_electrostatics()

void ddx_deallocate_electrostatics ( void *  electrostatics,
void *  error 
)

Deallocate the ddX electrostatics object.

◆ ddx_allocate_model()

void * ddx_allocate_model ( int  model,
int  enable_force,
double  solvent_epsilon,
double  solvent_kappa,
double  eta,
double  shift,
int  lmax,
int  n_lebedev,
int  incore,
int  maxiter,
int  jacobi_n_diis,
int  enable_fmm,
int  fmm_multipole_lmax,
int  fmm_local_lmax,
int  n_proc,
int  n_spheres,
const double *  sphere_centres,
const double *  sphere_radii,
int  length_logfile,
const char *  logfile,
void *  error 
)

Allocate a ddx model object.

Parameters
modelInteger describing the model to use.
enable_forceIf 1 enable force calculations
solvent_epsilonRelative dielectric permittivity
solvent_kappaDebye-Hückel parameter (inverse screening length)
etaRegularization parameter for the regularised characteristic function chi_eta. Range [0,1]
shiftShift for the regularized characteristic function chi_eta
lmaxMaximal degree of modelling spherical harmonics
n_lebedevNumber of Lebedev grid points to use
incoreIf 1 store more large objects in memory
maxiterMaximal number of iterations
jacobi_n_diisNumber of iterates stored in the DIIS space for acceleration
enable_fmmIf 1 enable the fast-multipole method for computations
fmm_multipole_lmaxMaximal degree of multipole spherical harmonics, ignored in case !enable_fmm. Value -1 means no far-field FFM interactions are computed.
fmm_local_lmaxMaximal degree of local spherical harmonics, ignored in case use_fmm=false. Value -1 means no local FFM interactions are computed.
n_procNumber of processors to use
n_spheresNumber of cavity spheres
sphere_centresThe centres of the cavity spheres as a column-major (3, n_spheres) array.
sphere_radiiThe radii of the spheres.
length_logfileLength of the logfile name (0 if no log). Note that logs can be quite verbose for larger runs and should only be used for debugging.
logfileThe logfile name (empty if no log)
errorThe ddX error object

◆ ddx_deallocate_model()

void ddx_deallocate_model ( void *  ddx,
void *  error 
)

Deallocate the model object

◆ ddx_get_logfile()

void ddx_get_logfile ( const void *  ddx,
char *  logfile,
int  maxlen 
)

Copy the logfile used by the model to the passed logfile array. At most maxlen characters are copied.

◆ ddx_get_enable_fmm()

int ddx_get_enable_fmm ( const void *  ddx)

Return the current value of enable_fmm stored in the model.

◆ ddx_get_enable_force()

int ddx_get_enable_force ( const void *  ddx)

Return the current value of enable_force stored in the model.

◆ ddx_get_jacobi_n_diis()

int ddx_get_jacobi_n_diis ( const void *  ddx)

Return the current value of jacobi_n_diis stored in the model.

◆ ddx_get_lmax()

int ddx_get_lmax ( const void *  ddx)

Return the current value of lmax stored in the model.

◆ ddx_get_incore()

int ddx_get_incore ( const void *  ddx)

Return the current value of incore stored in the model.

◆ ddx_get_maxiter()

int ddx_get_maxiter ( const void *  ddx)

Return the current value of maxiter stored in the model.

◆ ddx_get_model()

int ddx_get_model ( const void *  ddx)

Return the current value of model stored in the model.

◆ ddx_get_n_lebedev()

int ddx_get_n_lebedev ( const void *  ddx)

Return the current value of n_lebedev stored in the model.

◆ ddx_get_n_proc()

int ddx_get_n_proc ( const void *  ddx)

Return the current value of n_proc stored in the model.

◆ ddx_get_n_spheres()

int ddx_get_n_spheres ( const void *  ddx)

Return the current value of n_spheres stored in the model.

◆ ddx_get_fmm_local_lmax()

int ddx_get_fmm_local_lmax ( const void *  ddx)

Return the current value of fmm_local_lmax stored in the model.

◆ ddx_get_fmm_multipole_lmax()

int ddx_get_fmm_multipole_lmax ( const void *  ddx)

Return the current value of fmm_multipole_lmax stored in the model.

◆ ddx_get_solvent_epsilon()

double ddx_get_solvent_epsilon ( const void *  ddx)

Return the current value of solvent_epsilon stored in the model.

◆ ddx_get_eta()

double ddx_get_eta ( const void *  ddx)

Return the current value of eta stored in the model.

◆ ddx_get_solvent_kappa()

double ddx_get_solvent_kappa ( const void *  ddx)

Return the current value of solvent_kappa stored in the model.

◆ ddx_get_shift()

double ddx_get_shift ( const void *  ddx)

Return the current value of shift stored in the model.

◆ ddx_get_sphere_charges()

void ddx_get_sphere_charges ( const void *  ddx,
int  nsph,
double *  c_charge 
)

Store the current values of sphere_charges into the passed pointer location. The pointer should refer to the memory location of an array of size (nsph). Data will be stored in column-major format.

◆ ddx_get_sphere_centres()

void ddx_get_sphere_centres ( const void *  ddx,
int  nsph,
double *  c_csph 
)

Store the current values of sphere_centres into the passed pointer location. The pointer should refer to the memory location of an array of size (nsph). Data will be stored in column-major format.

◆ ddx_get_sphere_radii()

void ddx_get_sphere_radii ( const void *  ddx,
int  nsph,
double *  c_rsph 
)

Store the current values of sphere_radii into the passed pointer location. The pointer should refer to the memory location of an array of size (nsph). Data will be stored in column-major format.

◆ ddx_get_n_basis()

int ddx_get_n_basis ( const void *  ddx)

Return the current value of n_basis stored in the model.

◆ ddx_get_n_cav()

int ddx_get_n_cav ( const void *  ddx)

Return the current value of n_cav stored in the model.

◆ ddx_get_cavity()

void ddx_get_cavity ( const void *  ddx,
int  ncav,
double *  c_ccav 
)

Store the current values of cavity into the passed pointer location. The pointer should refer to the memory location of an array of size (ncav). Data will be stored in column-major format.

◆ ddx_allocate_state()

void * ddx_allocate_state ( const void *  ddx,
void *  error 
)

Allocate an empty state from a model

◆ ddx_deallocate_state()

void ddx_deallocate_state ( void *  state,
void *  error 
)

Deallocate a state object

◆ ddx_get_x()

void ddx_get_x ( const void *  state,
int  nbasis,
int  nsph,
double *  x 
)

Get the solution of the forward problem stored inside the state as a (nbasis, nsph) array in column-major ordering.

Parameters
stateDDX state
nbasisNumber of basis functions used by DDX
nsphNumber of cavity spheres
xOutput pointer
Note
This function just returns the data stored in the state, which only becomes valid once the appropriate ddx_cosmo_solve or ddx_pcm_solve have been called for the state.

◆ ddx_get_x_niter()

int ddx_get_x_niter ( const void *  state)

Get the number of iterations which where needed to obtain the solution returned by ddx_get_x

◆ ddx_get_s()

void ddx_get_s ( const void *  state,
int  nbasis,
int  nsph,
double *  s 
)

Get the solution of the adjoint problem stored inside the state as a (nbasis, nsph) array in column-major ordering.

Parameters
stateDDX state
nbasisNumber of basis functions used by DDX
nsphNumber of cavity spheres
sOutput pointer
Note
This function just returns the data stored in the state, which only becomes valid once the appropriate ddx_cosmo_solve_adjoint or ddx_pcm_solve_adjoint have been called for the state.

◆ ddx_get_s_niter()

int ddx_get_s_niter ( const void *  state)

Get the number of iterations which where needed to obtain the adjoint solution returned by ddx_get_s

◆ ddx_get_xi()

void ddx_get_xi ( const void *  state,
const void *  ddx,
int  ncav,
double *  xi 
)

Return the xi, the adjoint solution s projected onto the cavity points.

Parameters
stateDDX state
ddxDDX model
ncavNumber of cavity points
xiArray to be filled with ncav

◆ ddx_get_zeta_dip()

void ddx_get_zeta_dip ( const void *  state,
const void *  ddx,
int  ncav,
double *  zeta_dip 
)

Return the zeta_dip, the dipolar contribution of the LPB adjoint solution, at the cavity points. The pointer should refer to the memory location of an array of size (3, ncav). Data will be stored in column-major format.

Parameters
stateDDX state
ddxDDX model
ncavNumber of cavity points
zeta_dipArray to be filled with ncav

◆ ddx_ddsolve()

double ddx_ddsolve ( const void *  ddx,
void *  state,
const void *  electrostatics,
int  nbasis,
int  nsph,
const double *  psi,
double  tol,
double *  forces,
const int  read_guess,
void *  error 
)

Solve everything and return the energy.

Parameters
ddxDDX model
stateDDX state
electrostaticsDDX electrostatic properties container
nbasisNumber of basis functions used by DDX
nsphNumber of cavity spheres
psiPsi array (nbasis, nsph)-shaped array (in column-major ordering)
tolTolerance for the linear system solver
forcesForce array
read_guessFlag for guess, if different from zero the guess is read from the state
errorDDX error

◆ ddx_setup()

void ddx_setup ( const void *  ddx,
void *  state,
const void *  electrostatics,
int  nbasis,
int  nsph,
const double *  psi,
void *  error 
)

In-place adjust the guess inside the state. Setup a problem in the passed state.

Parameters
ddxDDX model
stateDDX state
electrostaticsDDX electrostatic properties container
nbasisNumber of basis functions used by DDX
nsphNumber of cavity spheres
psiPsi array (nbasis, nsph)-shaped array (in column-major ordering)
errorDDX error

◆ ddx_fill_guess()

void ddx_fill_guess ( const void *  ddx,
void *  state,
double  tol,
void *  error 
)

In-place adjust the guess inside the state. Avoid calling this step if you want to use the currently stored solution as an initial guess

◆ ddx_fill_guess_adjoint()

void ddx_fill_guess_adjoint ( const void *  ddx,
void *  state,
double  tol,
void *  error 
)

In-place adjust the adjoint guess inside the state. problem. Avoid calling this step if you want to use the currently stored solution as an initial guess

◆ ddx_solve()

void ddx_solve ( const void *  ddx,
void *  state,
double  tol,
void *  error 
)

Solve the forward problem.

Parameters
stateDDX state
ddxDDX model
tolTolerance up to which the problem is solved
errorDDX error

◆ ddx_solve_adjoint()

void ddx_solve_adjoint ( const void *  ddx,
void *  state,
double  tol,
void *  error 
)

Solve the adjoint COSMO problem.

Parameters
stateDDX state
ddxDDX model
tolTolerance up to which the problem is solved
errorDDX error

◆ ddx_energy()

double ddx_energy ( const void *  ddx,
void *  state,
void *  error 
)

Compute the solvation energy.

Parameters
stateDDX state
ddxDDX model
errorDDX ERROR

◆ ddx_solvation_force_terms()

void ddx_solvation_force_terms ( const void *  ddx,
void *  state,
void *  electrostatics,
int  nsph,
double *  forces,
void *  error 
)

Compute the solvation force terms.

Parameters
ddxDDX model
stateDDX state
electrostaticsDDX electrostatic properties container
nsphNumber of cavity spheres
forcesOutput force array (3, nsph) in column-major order
errorDDX error

◆ ddx_ddrun()

double ddx_ddrun ( const void *  ddx,
void *  state,
void *  electrostatics,
int  nbasis,
int  nsph,
double *  psi,
const double  tol,
double *  forces,
int  read_guess,
void *  error 
)

Run all the steps of a ddX calculation and return energy and, if requested, forces.

Parameters
ddxDDX model
stateDDX state
electrostaticsDDX electrostatic properties container
nbasisNumber of basis functions used by DDX
nsphNumber of cavity spheres
psiPsi array (nbasis, nsph)-shaped array (in column-major ordering)
tolTolerance up to which the problem is solved
forcesOutput force array (3, nsph) in column-major order
read_guessGuess control flag: 0 for default guess, 1 for using the content of DDX state for the guess.
errorDDX error

◆ ddx_cosmo_setup()

void ddx_cosmo_setup ( const void *  ddx,
void *  state,
int  ncav,
int  nbasis,
int  nsph,
const double *  psi,
const double *  phi_cav,
void *  error 
)

Setup a COSMO problem in the passed state.

Parameters
ddxDDX model
stateDDX state
ncavNumber of cavity points
nbasisNumber of basis functions used by DDX
nsphNumber of cavity spheres
psiPsi array (nbasis, nsph)-shaped array (in column-major ordering)
phi_cavPhi array adjoint (ncav, )-shaped array

◆ ddx_cosmo_guess()

void ddx_cosmo_guess ( const void *  ddx,
void *  state,
void *  error 
)

In-place adjust the guess inside the state, getting ready to solve a COSMO problem. Avoid calling this step if you want to use the currently stored solution as an initial guess

◆ ddx_cosmo_guess_adjoint()

void ddx_cosmo_guess_adjoint ( const void *  ddx,
void *  state,
void *  error 
)

In-place adjust the guess inside the state, getting ready to solve an adjoint COSMO problem. Avoid calling this step if you want to use the currently stored solution as an initial guess

◆ ddx_cosmo_solve()

void ddx_cosmo_solve ( const void *  ddx,
void *  state,
double  tol,
void *  error 
)

Solve the COSMO problem.

Parameters
stateDDX state
ddxDDX model
tolTolerance up to which the problem is solved
errorDDX error

◆ ddx_cosmo_solve_adjoint()

void ddx_cosmo_solve_adjoint ( const void *  ddx,
void *  state,
double  tol,
void *  error 
)

Solve the adjoint COSMO problem.

Parameters
stateDDX state
ddxDDX model
tolTolerance up to which the problem is solved
errorDDX error

◆ ddx_cosmo_energy()

double ddx_cosmo_energy ( const void *  ddx,
void *  state,
void *  error 
)

Compute COSMO energy (without any scaling by (epsilon - 1) / epsilon or similar).

Parameters
stateDDX state
ddxDDX model
errorDDX ERROR

◆ ddx_cosmo_solvation_force_terms()

void ddx_cosmo_solvation_force_terms ( const void *  ddx,
void *  state,
int  nsph,
int  ncav,
const double *  e_cav,
double *  forces,
void *  error 
)

Compute the COSMO force terms.

Parameters
ddxDDX model
stateDDX state
nsphNumber of cavity spheres
e_cavElectric field
forcesOutput force array (3, nsph) in column-major order
errorDDX error

◆ ddx_pcm_setup()

void ddx_pcm_setup ( const void *  ddx,
void *  state,
int  ncav,
int  nbasis,
int  nsph,
const double *  psi,
const double *  phi_cav,
void *  error 
)

Setup a PCM problem in the passed state.

Parameters
ddxDDX model
stateDDX state
ncavNumber of cavity points
nbasisNumber of basis functions used by DDX
nsphNumber of cavity spheres
psiPsi array (nbasis, nsph)-shaped array (in column-major ordering)
phi_cavPhi array (ncav, )-shaped array
errorDDX error

◆ ddx_pcm_guess()

void ddx_pcm_guess ( const void *  ddx,
void *  state,
void *  error 
)

In-place adjust the guess inside the state, getting ready to solve a PCM problem. Avoid calling this step if you want to use the currently stored solution as an initial guess

◆ ddx_pcm_guess_adjoint()

void ddx_pcm_guess_adjoint ( const void *  ddx,
void *  state,
void *  error 
)

In-place adjust the guess inside the state, getting ready to solve an adjoint PCM problem. Avoid calling this step if you want to use the currently stored solution as an initial guess

◆ ddx_pcm_solve()

void ddx_pcm_solve ( const void *  ddx,
void *  state,
double  tol,
void *  error 
)

Solve the forward PCM problem.

Parameters
stateDDX state
ddxDDX model
errorDDX error
tolTolerance up to which the problem is solved

◆ ddx_pcm_solve_adjoint()

void ddx_pcm_solve_adjoint ( const void *  ddx,
void *  state,
double  tol,
void *  error 
)

Solve the adjoint PCM problem.

Parameters
stateDDX state
ddxDDX model
errorDDX error
tolTolerance up to which the problem is solved

◆ ddx_pcm_energy()

double ddx_pcm_energy ( const void *  ddx,
void *  state,
void *  error 
)

Compute PCM energy

Parameters
stateDDX state
errorDDX error
ddxDDX model

◆ ddx_pcm_solvation_force_terms()

void ddx_pcm_solvation_force_terms ( const void *  ddx,
void *  state,
int  nsph,
int  ncav,
const double *  e_cav,
double *  forces,
void *  error 
)

Compute the PCM force terms.

Parameters
ddxDDX model
stateDDX state
nsphNumber of cavity spheres
e_cavElectric field
forcesOutput force array (3, nsph) in column-major order
errorDDX error

◆ ddx_lpb_setup()

void ddx_lpb_setup ( const void *  ddx,
void *  state,
int  ncav,
int  nbasis,
int  nsph,
const double *  psi,
const double *  phi_cav,
const double *  e_cav,
void *  error 
)

Setup a LPB problem in the passed state.

Parameters
ddxDDX model
stateDDX state
ncavNumber of cavity points
nbasisNumber of basis functions used by DDX
nsphNumber of cavity spheres
psiPsi array (nbasis, nsph)-shaped array (in column-major ordering)
phi_cavPhi array (ncav, )-shaped array
e_cavElectric field generated by multipoles: (3, ncav)-shaped array (column major)
errorDDX error

◆ ddx_lpb_guess()

void ddx_lpb_guess ( const void *  ddx,
void *  state,
double  tol,
void *  error 
)

In-place adjust the guess inside the state, getting ready to solve a LPB problem. Avoid calling this step if you want to use the currently stored solution as an initial guess. tol is the tolerance for solving a simplified initial-guess problem. The same tolerance as for the ddx_lpb_solve call should be chosen.

◆ ddx_lpb_guess_adjoint()

void ddx_lpb_guess_adjoint ( const void *  ddx,
void *  state,
double  tol,
void *  error 
)

In-place adjust the guess inside the state, getting ready to solve an adjoint LPB problem. Avoid calling this step if you want to use the currently stored solution as an initial guess. tol is the tolerance for solving a simplified initial-guess problem. The same tolerance as for the ddx_lpb_solve call should be chosen.

◆ ddx_lpb_solve()

void ddx_lpb_solve ( const void *  ddx,
void *  state,
double  tol,
void *  error 
)

Solve the forward LPB problem.

Parameters
stateDDX state
ddxDDX model
errorDDX error
tolTolerance up to which the problem is solved

◆ ddx_lpb_solve_adjoint()

void ddx_lpb_solve_adjoint ( const void *  ddx,
void *  state,
double  tol,
void *  error 
)

Solve the adjoint LPB problem.

Parameters
stateDDX state
ddxDDX model
errorDDX error
tolTolerance up to which the problem is solved

◆ ddx_lpb_energy()

double ddx_lpb_energy ( const void *  ddx,
void *  state,
void *  error 
)

Compute LPB energy

Parameters
stateDDX state
errorDDX error
ddxDDX model

◆ ddx_lpb_solvation_force_terms()

void ddx_lpb_solvation_force_terms ( const void *  ddx,
void *  state,
int  nsph,
int  ncav,
const double *  g_cav,
double *  forces,
const void *  error 
)

Compute the LPB force terms.

Parameters
ddxDDX model
stateDDX state
nsphNumber of cavity spheres
g_cavElectric field gradient
forcesOutput force array (3, nsph) in column-major order
errorDDX error

◆ ddx_multipole_electrostatics_0()

void ddx_multipole_electrostatics_0 ( const void *  ddx,
int  nsph,
int  ncav,
int  nmultipoles,
const double *  multipoles,
double *  phi_cav,
void *  error 
)

Build potential generated by a multipolar charge distribution.

Parameters
ddxDDX model
nsphNumber of cavity spheres
ncavNumber of cavity points
nmultipolesTotal number of multipoles. If multipoles up to mmax are used, this is (mmax+1)^2. E.g. for charges, dipoles and quadrupoles (mmax=2), this is 9.
multipolesMultipoles (nmultipoles, nsph)-shaped array in column-major order.
phi_cavPotential (phi) generated by multipoles: (ncav, )-shaped array
errorDDX error

◆ ddx_multipole_electrostatics_1()

void ddx_multipole_electrostatics_1 ( const void *  ddx,
int  nsph,
int  ncav,
int  nmultipoles,
const double *  multipoles,
double *  phi_cav,
double *  e_cav,
void *  error 
)

Build potential and electric field generated by a multipolar charge distribution.

Parameters
ddxDDX model
nsphNumber of cavity spheres
ncavNumber of cavity points
nmultipolesTotal number of multipoles. If multipoles up to mmax are used, this is (mmax+1)^2. E.g. for charges, dipoles and quadrupoles (mmax=2), this is 9.
multipolesMultipoles (nmultipoles, nsph)-shaped array in column-major order.
phi_cavPotential (phi) generated by multipoles: (ncav, )-shaped array
e_cavElectric field generated by multipoles: (3, ncav)-shaped array (column major)
errorDDX error

◆ ddx_multipole_electrostatics_2()

void ddx_multipole_electrostatics_2 ( const void *  ddx,
int  nsph,
int  ncav,
int  nmultipoles,
const double *  multipoles,
double *  phi_cav,
double *  e_cav,
double *  g_cav,
void *  error 
)

Build potential, electric field and field gradient generated by a multipolar charge distribution.

Parameters
ddxDDX model
nsphNumber of cavity spheres
ncavNumber of cavity points
nmultipolesTotal number of multipoles. If multipoles up to mmax are used, this is (mmax+1)^2. E.g. for charges, dipoles and quadrupoles (mmax=2), this is 9.
multipolesMultipoles (nmultipoles, nsph)-shaped array in column-major order.
phi_cavPotential (phi) generated by multipoles: (ncav, )-shaped array
e_cavElectric field generated by multipoles: (3, ncav)-shaped array (column major)
g_cavElectric field gradient gen. by multipoles: (3, 3, ncav)-shaped array (column major)
errorDDX error

◆ ddx_multipole_psi()

void ddx_multipole_psi ( const void *  ddx,
int  nbasis,
int  nsph,
int  nmultipoles,
const double *  multipoles,
double *  psi,
void *  error 
)

Build the Psi generated by a multipolar charge distribution

Parameters
ddxDDX model
nbasisNumber of basis functions used by DDX
nsphNumber of cavity spheres
nmultipolesTotal number of multipoles. If multipoles up to mmax are used, this is (mmax+1)^2. E.g. for charges, dipoles and quadrupoles (mmax=2), this is 9.
multipolesMultipoles (nmultipoles, nsph)-shaped array in column-major order.
psiPsi array (nbasis, nsph)-shaped array (column major)
errorDDX error

◆ ddx_multipole_force_terms()

void ddx_multipole_force_terms ( const void *  ddx,
void *  state,
int  nsph,
int  nmultipoles,
const double *  multipoles,
double *  forces,
void *  error 
)

Compute the force terms generated by a multipolar charge distribution

Parameters
ddxDDX model
stateDDX state
nsphNumber of cavity spheres
nmultipolesTotal number of multipoles. If multipoles up to mmax are used, this is (mmax+1)^2. E.g. for charges, dipoles and quadrupoles (mmax=2), this is 9.
multipolesMultipoles (nmultipoles, nsph)-shaped array in column-major order.
forcesOutput force array (3, nsph) in column-major order
errorDDX error