ddx 0.6.0
Libary for domain-decomposition methods for polarizable continuum models
ddx.h
Go to the documentation of this file.
1#ifndef DDX_H
2#define DDX_H
3
14#ifdef __cplusplus
15extern "C" {
16#endif
17
20
24void ddx_get_banner(char* banner, int maxlen);
25
32int ddx_supported_lebedev_grids(int maxlen, int* grids);
33
42void ddx_scaled_ylm(const void* ddx, int lmax, const double* x, int sphere, double* ylm);
44
47
50
54int ddx_get_error_flag(const void* error);
55
59void ddx_get_error_message(const void* error, char* message, int maxlen);
60
62
65
67void* ddx_allocate_electrostatics(void* ddx, void* error);
68
70void ddx_multipole_electrostatics(void* ddx, int nsph, int nmultipoles,
71 const double* multipoles,
72 void* electrostatics, void* error);
73
75void ddx_deallocate_electrostatics(void* electrostatics, void* error);
76
78
81
114void* ddx_allocate_model(int model, int enable_force, double solvent_epsilon,
115 double solvent_kappa, double eta, double shift, int lmax,
116 int n_lebedev, int incore, int maxiter, int jacobi_n_diis,
117 int enable_fmm, int fmm_multipole_lmax, int fmm_local_lmax,
118 int n_proc, int n_spheres, const double* sphere_centres,
119 const double* sphere_radii, int length_logfile,
120 const char* logfile, void* error);
121
123void ddx_deallocate_model(void* ddx, void* error);
124
126
129
132void ddx_get_logfile(const void* ddx, char* logfile, int maxlen);
133
134// Generated block, see scripts/generate_cinterface_setup.py.
136int ddx_get_enable_fmm(const void* ddx);
137
139int ddx_get_enable_force(const void* ddx);
140
142int ddx_get_jacobi_n_diis(const void* ddx);
143
145int ddx_get_lmax(const void* ddx);
146
148int ddx_get_incore(const void* ddx);
149
151int ddx_get_maxiter(const void* ddx);
152
154int ddx_get_model(const void* ddx);
155
157int ddx_get_n_lebedev(const void* ddx);
158
160int ddx_get_n_proc(const void* ddx);
161
163int ddx_get_n_spheres(const void* ddx);
164
167
170
172double ddx_get_solvent_epsilon(const void* ddx);
173
175double ddx_get_eta(const void* ddx);
176
178double ddx_get_solvent_kappa(const void* ddx);
179
181double ddx_get_shift(const void* ddx);
182
186void ddx_get_sphere_charges(const void* ddx, int nsph, double* c_charge);
187
191void ddx_get_sphere_centres(const void* ddx, int nsph, double* c_csph);
192
196void ddx_get_sphere_radii(const void* ddx, int nsph, double* c_rsph);
197
199int ddx_get_n_basis(const void* ddx);
200
202int ddx_get_n_cav(const void* ddx);
203
207void ddx_get_cavity(const void* ddx, int ncav, double* c_ccav);
208// end generated block
210
213
214void* ddx_allocate_state(const void* ddx, void* error);
215
217void ddx_deallocate_state(void* state, void* error);
218
229void ddx_get_x(const void* state, int nbasis, int nsph, double* x);
230
233int ddx_get_x_niter(const void* state);
234
245void ddx_get_s(const void* state, int nbasis, int nsph, double* s);
246
249int ddx_get_s_niter(const void* state);
250
257void ddx_get_xi(const void* state, const void* ddx, int ncav, double* xi);
258
268void ddx_get_zeta_dip(const void* state, const void* ddx, int ncav, double* zeta_dip);
270
273
286double ddx_ddsolve(const void* ddx, void* state, const void* electrostatics,
287 int nbasis, int nsph, const double* psi, double tol,
288 double* forces, const int read_guess, void* error);
289
300void ddx_setup(const void* ddx, void* state, const void* electrostatics,
301 int nbasis, int nsph, const double* psi, void* error);
302
306void ddx_fill_guess(const void* ddx, void* state, double tol, void* error);
307
311void ddx_fill_guess_adjoint(const void* ddx, void* state, double tol, void* error);
312
318void ddx_solve(const void* ddx, void* state, double tol, void* error);
319
325void ddx_solve_adjoint(const void* ddx, void* state, double tol, void* error);
326
332double ddx_energy(const void* ddx, void* state, void* error);
333
342void ddx_solvation_force_terms(const void* ddx, void* state, void* electrostatics,
343 int nsph, double* forces, void* error);
344
360double ddx_ddrun(const void* ddx, void* state, void* electrostatics, int nbasis,
361 int nsph, double* psi, const double tol, double* forces,
362 int read_guess, void* error);
363
365
368
378void ddx_cosmo_setup(const void* ddx, void* state, int ncav, int nbasis, int nsph,
379 const double* psi, const double* phi_cav, void* error);
380
384void ddx_cosmo_guess(const void* ddx, void* state, void* error);
385
389void ddx_cosmo_guess_adjoint(const void* ddx, void* state, void* error);
390
396void ddx_cosmo_solve(const void* ddx, void* state, double tol, void* error);
397
403void ddx_cosmo_solve_adjoint(const void* ddx, void* state, double tol, void* error);
404
411double ddx_cosmo_energy(const void* ddx, void* state, void* error);
412
421void ddx_cosmo_solvation_force_terms(const void* ddx, void* state, int nsph,
422 int ncav, const double* e_cav,
423 double* forces, void* error);
424
435void ddx_pcm_setup(const void* ddx, void* state, int ncav, int nbasis, int nsph,
436 const double* psi, const double* phi_cav, void* error);
437
441void ddx_pcm_guess(const void* ddx, void* state, void* error);
442
446void ddx_pcm_guess_adjoint(const void* ddx, void* state, void* error);
447
453void ddx_pcm_solve(const void* ddx, void* state, double tol, void* error);
454
460void ddx_pcm_solve_adjoint(const void* ddx, void* state, double tol, void* error);
461
466double ddx_pcm_energy(const void* ddx, void* state, void* error);
467
476void ddx_pcm_solvation_force_terms(const void* ddx, void* state, int nsph,
477 int ncav, const double* e_cav,
478 double* forces, void* error);
479
492void ddx_lpb_setup(const void* ddx, void* state, int ncav, int nbasis, int nsph,
493 const double* psi, const double* phi_cav, const double* e_cav,
494 void* error);
495
500void ddx_lpb_guess(const void* ddx, void* state, double tol, void* error);
501
506void ddx_lpb_guess_adjoint(const void* ddx, void* state, double tol, void* error);
507
513void ddx_lpb_solve(const void* ddx, void* state, double tol, void* error);
514
520void ddx_lpb_solve_adjoint(const void* ddx, void* state, double tol, void* error);
521
526double ddx_lpb_energy(const void* ddx, void* state, void* error);
527
536void ddx_lpb_solvation_force_terms(const void* ddx, void* state, int nsph,
537 int ncav, const double* g_cav,
538 double* forces, const void* error);
539
541
544
557void ddx_multipole_electrostatics_0(const void* ddx, int nsph, int ncav, int nmultipoles,
558 const double* multipoles, double* phi_cav, void* error);
559
574void ddx_multipole_electrostatics_1(const void* ddx, int nsph, int ncav, int nmultipoles,
575 const double* multipoles, double* phi_cav,
576 double* e_cav, void* error);
577
595void ddx_multipole_electrostatics_2(const void* ddx, int nsph, int ncav, int nmultipoles,
596 const double* multipoles, double* phi_cav,
597 double* e_cav, double* g_cav, void* error);
598
611void ddx_multipole_psi(const void* ddx, int nbasis, int nsph, int nmultipoles,
612 const double* multipoles, double* psi, void* error);
613
626void ddx_multipole_force_terms(const void* ddx, void* state, int nsph,
627 int nmultipoles, const double* multipoles,
628 double* forces, void* error);
629
632
633#ifdef __cplusplus
634}
635#endif
636
637#endif /* DDX_H */
void ddx_get_sphere_centres(const void *ddx, int nsph, double *c_csph)
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_get_logfile(const void *ddx, char *logfile, int maxlen)
void ddx_get_cavity(const void *ddx, int ncav, double *c_ccav)
void ddx_multipole_electrostatics_0(const void *ddx, int nsph, int ncav, int nmultipoles, const double *multipoles, double *phi_cav, void *error)
void ddx_lpb_guess_adjoint(const void *ddx, void *state, double tol, void *error)
void ddx_fill_guess(const void *ddx, void *state, double tol, void *error)
void ddx_get_s(const void *state, int nbasis, int nsph, double *s)
void ddx_lpb_solvation_force_terms(const void *ddx, void *state, int nsph, int ncav, const double *g_cav, double *forces, const void *error)
void * ddx_allocate_state(const void *ddx, void *error)
void ddx_solve_adjoint(const void *ddx, void *state, double tol, void *error)
double ddx_cosmo_energy(const void *ddx, void *state, void *error)
double ddx_get_solvent_kappa(const void *ddx)
void ddx_pcm_solve_adjoint(const void *ddx, void *state, double tol, void *error)
void ddx_cosmo_solve_adjoint(const void *ddx, void *state, double tol, void *error)
int ddx_get_n_cav(const void *ddx)
int ddx_get_error_flag(const void *error)
double ddx_pcm_energy(const void *ddx, void *state, void *error)
void ddx_cosmo_guess_adjoint(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)
void ddx_multipole_electrostatics(void *ddx, int nsph, int nmultipoles, const double *multipoles, void *electrostatics, 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_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_cosmo_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_adjoint(const void *ddx, void *state, void *error)
int ddx_get_n_lebedev(const void *ddx)
double ddx_get_solvent_epsilon(const void *ddx)
void ddx_get_error_message(const void *error, char *message, int maxlen)
int ddx_get_n_proc(const void *ddx)
void ddx_multipole_psi(const void *ddx, int nbasis, int nsph, int nmultipoles, const double *multipoles, double *psi, void *error)
void ddx_cosmo_guess(const void *ddx, void *state, void *error)
int ddx_get_n_spheres(const void *ddx)
void ddx_get_zeta_dip(const void *state, const void *ddx, int ncav, double *zeta_dip)
int ddx_get_s_niter(const void *state)
void ddx_get_banner(char *banner, int maxlen)
int ddx_get_enable_fmm(const void *ddx)
int ddx_get_lmax(const void *ddx)
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_get_x(const void *state, int nbasis, int nsph, double *x)
void ddx_multipole_force_terms(const void *ddx, void *state, int nsph, int nmultipoles, const double *multipoles, double *forces, void *error)
double ddx_lpb_energy(const void *ddx, void *state, void *error)
void ddx_get_xi(const void *state, const void *ddx, int ncav, double *xi)
int ddx_get_incore(const void *ddx)
int ddx_get_x_niter(const void *state)
double ddx_energy(const void *ddx, void *state, void *error)
void ddx_deallocate_electrostatics(void *electrostatics, void *error)
void ddx_deallocate_model(void *ddx, void *error)
void ddx_scaled_ylm(const void *ddx, int lmax, const double *x, int sphere, double *ylm)
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_setup(const void *ddx, void *state, const void *electrostatics, int nbasis, int nsph, const double *psi, void *error)
double ddx_get_eta(const void *ddx)
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_lpb_solve(const void *ddx, void *state, double tol, void *error)
int ddx_get_maxiter(const void *ddx)
int ddx_get_enable_force(const void *ddx)
int ddx_get_model(const void *ddx)
double ddx_get_shift(const void *ddx)
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)
int ddx_supported_lebedev_grids(int maxlen, int *grids)
void ddx_deallocate_state(void *state, void *error)
int ddx_get_n_basis(const void *ddx)
void ddx_lpb_guess(const void *ddx, void *state, double tol, void *error)
void ddx_cosmo_solve(const void *ddx, void *state, double tol, void *error)
void ddx_lpb_solve_adjoint(const void *ddx, void *state, double tol, void *error)
void ddx_get_sphere_radii(const void *ddx, int nsph, double *c_rsph)
void * ddx_allocate_error()
void ddx_fill_guess_adjoint(const void *ddx, void *state, double tol, void *error)
void * ddx_allocate_electrostatics(void *ddx, void *error)
int ddx_get_jacobi_n_diis(const void *ddx)
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_pcm_guess(const void *ddx, void *state, void *error)
void ddx_get_sphere_charges(const void *ddx, int nsph, double *c_charge)
void ddx_solve(const void *ddx, void *state, double tol, void *error)
void ddx_pcm_solve(const void *ddx, void *state, double tol, void *error)
int ddx_get_fmm_multipole_lmax(const void *ddx)
int ddx_get_fmm_local_lmax(const void *ddx)
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)
High-level module of the ddX software.
Definition: ddx.f90:2