Interfacing a QM program and PCMSolver

For the impatients: tl;dr

In these examples, we want to show how every function in the API works. If your program is written in Fortran, head over to Interfacing with a Fortran host If your program is written in C/C++, head over to Interfacing with a C host

How PCMSolver handles potentials and charges: surface functions

Electrostatic potential vectors and the corresponding apparent surface charge vectors are handled internally as surface functions. The actual values are stored into Eigen vectors and saved into a map. The mapping is between the name of the surface function, given by the programmer writing the interface to the library, and the vector holding the values.

What you should care about: API functions

These are the contents of the pcmsolver.h file defining the public API of the PCMSolver library. The Fortran bindings for the API are in the pcmsolver.f90 file. The indexing of symmetry operations and their mapping to a bitstring is explained in the following Table. This is important when passing symmetry information to the pcmsolver_new() function.

Symmetry operations indexing within the module
Index zyx Generator Parity
0 000 E 1.0
1 001 Oyz -1.0
2 010 Oxz -1.0
3 011 C2z 1.0
4 100 Oxy -1.0
5 101 C2y 1.0
6 110 C2x 1.0
7 111 i -1.0

Defines

PCMSolver_EXPORT
pcmsolver_bool_t_DEFINED

Typedefs

typedef bool pcmsolver_bool_t
typedef pcmsolver_context_t

Workaround to have pcmsolver_context_t available to C

typedef HostWriter

Flushes module output to host program

Parameters
  • message: contents of the module output

Enums

enum pcmsolver_reader_t

Input processing strategies.

Values:

PCMSOLVER_READER_OWN

Module reads input on its own

PCMSOLVER_READER_HOST

Module receives input from host

Functions

pcmsolver_context_t *pcmsolver_new(pcmsolver_reader_t input_reading, int nr_nuclei, double charges[], double coordinates[], int symmetry_info[], struct PCMInput *host_input, HostWriter writer)

Creates a new PCM context object.

The molecular point group information is passed as an array of 4 integers: number of generators, first, second and third generator respectively. Generators map to integers as in table :ref:symmetry-ops

Parameters
  • input_reading: input processing strategy
  • nr_nuclei: number of atoms in the molecule
  • charges: atomic charges
  • coordinates: atomic coordinates
  • symmetry_info: molecular point group information
  • host_input: input to the module, as read by the host
  • writer: flush-to-host function

pcmsolver_context_t *pcmsolver_new_v1112(pcmsolver_reader_t input_reading, int nr_nuclei, double charges[], double coordinates[], int symmetry_info[], const char *parsed_fname, struct PCMInput *host_input, HostWriter writer)

Creates a new PCM context object, updated in v1.1.12.

The molecular point group information is passed as an array of 4 integers: number of generators, first, second and third generator respectively. Generators map to integers as in table :ref:symmetry-ops

Parameters
  • input_reading: input processing strategy
  • nr_nuclei: number of atoms in the molecule
  • charges: atomic charges
  • coordinates: atomic coordinates
  • symmetry_info: molecular point group information
  • parsed_fname: name of the input file parsed by pcmsolver.py
  • host_input: input to the module, as read by the host
  • writer: flush-to-host function

void pcmsolver_delete(pcmsolver_context_t *context)

Deletes a PCM context object.

Parameters
  • context: the PCM context object to be deleted

pcmsolver_bool_t pcmsolver_is_compatible_library(void)

Whether the library is compatible with the header file Checks that the compiled library and header file version match. Host should abort when that is not the case.

Warning
This function should be called before instantiating any PCM context objects.

void pcmsolver_print(pcmsolver_context_t *context)

Prints citation and set up information.

Parameters
  • context: the PCM context object

int pcmsolver_get_cavity_size(pcmsolver_context_t *context)

Getter for the number of finite elements composing the molecular cavity.

Return
the size of the cavity
Parameters
  • context: the PCM context object

int pcmsolver_get_irreducible_cavity_size(pcmsolver_context_t *context)

Getter for the number of irreducible finite elements composing the molecular cavity.

Return
the number of irreducible finite elements
Parameters
  • context: the PCM context object

void pcmsolver_get_centers(pcmsolver_context_t *context, double centers[])

Getter for the centers of the finite elements composing the molecular cavity.

Parameters
  • context: the PCM context object
  • centers: array holding the coordinates of the finite elements centers

void pcmsolver_get_center(pcmsolver_context_t *context, int its, double center[])

Getter for the center of the i-th finite element.

Parameters
  • context: the PCM context object
  • its: index of the finite element
  • center: array holding the coordinates of the finite element center

void pcmsolver_get_areas(pcmsolver_context_t *context, double areas[])

Getter for the areas/weights of the finite elements.

Parameters
  • context: the PCM context object
  • areas: array holding the weights/areas of the finite elements

void pcmsolver_compute_asc(pcmsolver_context_t *context, const char *mep_name, const char *asc_name, int irrep)

Computes ASC given a MEP and the desired irreducible representation.

Parameters
  • context: the PCM context object
  • mep_name: label of the MEP surface function
  • asc_name: label of the ASC surface function
  • irrep: index of the desired irreducible representation The module uses the surface function concept to handle potentials and charges. Given labels for each, this function retrieves the MEP and computes the corresponding ASC.

void pcmsolver_compute_response_asc(pcmsolver_context_t *context, const char *mep_name, const char *asc_name, int irrep)

Computes response ASC given a MEP and the desired irreducible representation.

Parameters
  • context: the PCM context object
  • mep_name: label of the MEP surface function
  • asc_name: label of the ASC surface function
  • irrep: index of the desired irreducible representation If Nonequilibrium = True in the input, calculates a response ASC using the dynamic permittivity. Falls back to the solver with static permittivity otherwise.

double pcmsolver_compute_polarization_energy(pcmsolver_context_t *context, const char *mep_name, const char *asc_name)

Computes the polarization energy.

Return
the polarization energy This function calculates the dot product of the given MEP and ASC vectors.
Parameters
  • context: the PCM context object
  • mep_name: label of the MEP surface function
  • asc_name: label of the ASC surface function

double pcmsolver_get_asc_dipole(pcmsolver_context_t *context, const char *asc_name, double dipole[])

Getter for the ASC dipole.

Return
the ASC dipole, i.e. { ^2}
Parameters
  • context: the PCM context object
  • asc_name: label of the ASC surface function
  • dipole: the Cartesian components of the ASC dipole moment

void pcmsolver_get_surface_function(pcmsolver_context_t *context, int size, double values[], const char *name)

Retrieves data wrapped in a given surface function.

Parameters
  • context: the PCM context object
  • size: the size of the surface function
  • values: the values wrapped in the surface function
  • name: label of the surface function

void pcmsolver_set_surface_function(pcmsolver_context_t *context, int size, double values[], const char *name)

Sets a surface function given data and label.

Parameters
  • context: the PCM context object
  • size: the size of the surface function
  • values: the values to be wrapped in the surface function
  • name: label of the surface function

void pcmsolver_print_surface_function(pcmsolver_context_t *context, const char *name)

Prints surface function contents to host output.

Parameters
  • context: the PCM context object
  • name: label of the surface function

void pcmsolver_save_surface_functions(pcmsolver_context_t *context)

Dumps all currently saved surface functions to NumPy arrays in .npy files.

Parameters
  • context: the PCM context object

void pcmsolver_save_surface_function(pcmsolver_context_t *context, const char *name)

Dumps a surface function to NumPy array in .npy file.

Note
The name parameter is the name of the NumPy array file without .npy extension
Parameters
  • context: the PCM context object
  • name: label of the surface function

void pcmsolver_load_surface_function(pcmsolver_context_t *context, const char *name)

Loads a surface function from a .npy file.

Note
The name parameter is the name of the NumPy array file without .npy extension
Parameters
  • context: the PCM context object
  • name: label of the surface function

void pcmsolver_write_timings(pcmsolver_context_t *context)

Writes timing results for the API.

Parameters
  • context: the PCM context object

Host input forwarding

struct PCMInput

Data structure for host-API input communication.

Forward-declare PCMInput input wrapping struct

Public Members

char cavity_type[8]

Type of cavity requested.

int patch_level

Wavelet cavity mesh patch level.

double coarsity

Wavelet cavity mesh coarsity.

double area

Average tesserae area.

char radii_set[8]

The built-in radii set to be used.

double min_distance

Minimal distance between sampling points.

int der_order

Derivative order for the switching function.

pcmsolver_bool_t scaling

Whether to scale or not the atomic radii.

char restart_name[20]

Name of the .npz file for GePol cavity restart.

double min_radius

Minimal radius for the added spheres.

char solver_type[7]

Type of solver requested.

double correction

Correction in the CPCM apparent surface charge scaling factor.

char solvent[16]

Name of the solvent.

double probe_radius

Radius of the spherical probe mimicking the solvent.

char equation_type[11]

Type of the integral equation to be used.

char inside_type[7]

Type of Green’s function requested inside the cavity.

double outside_epsilon

Value of the static permittivity outside the cavity.

char outside_type[22]

Type of Green’s function requested outside the cavity.

Internal details of the API

class pcm::Meddle

Contains functions exposing an interface to the module internals.

Author
Roberto Di Remigio
Date
2015-2017

Public Functions

Meddle(const Input &input, const HostWriter &writer)

CTOR from Input object.

Warning
This CTOR is meant to be used with the standalone executable only
Parameters
  • input: an Input object
  • writer: the global HostWriter object

Meddle(const std::string &inputFileName, const HostWriter &writer)

CTOR from own input reader.

Warning
This CTOR is meant to be used with the standalone executable only
Parameters
  • inputFileName: name of the parsed, machine-readable input file
  • writer: the global HostWriter object

Meddle(int nr_nuclei, double charges[], double coordinates[], int symmetry_info[], const HostWriter &writer, const std::string &inputFileName = "@pcmsolver.inp")

CTOR from parsed input file name.

Parameters
  • inputFileName: name of the parsed, machine-readable input file
  • nr_nuclei: number of atoms in the molecule
  • charges: atomic charges
  • coordinates: atomic coordinates
  • symmetry_info: molecular point group information
  • writer: the global HostWriter object

Meddle(int nr_nuclei, double charges[], double coordinates[], int symmetry_info[], const PCMInput &host_input, const HostWriter &writer)

Constructor.

The molecular point group information is passed as an array of 4 integers: number of generators, first, second and third generator respectively. Generators map to integers as in table :ref:symmetry-ops

Parameters
  • nr_nuclei: number of atoms in the molecule
  • charges: atomic charges
  • coordinates: atomic coordinates
  • symmetry_info: molecular point group information
  • host_input: input to the module, as read by the host
  • writer: the global HostWriter object

Molecule molecule() const

Getter for the molecule object.

PCMSolverIndex getCavitySize() const

Getter for the number of finite elements composing the molecular cavity.

Return
the size of the cavity

PCMSolverIndex getIrreducibleCavitySize() const

Getter for the number of irreducible finite elements composing the molecular cavity.

Return
the number of irreducible finite elements

void getCenters(double centers[]) const

Getter for the centers of the finite elements composing the molecular cavity.

Parameters
  • centers: array holding the coordinates of the finite elements centers

void getCenter(int its, double center[]) const

Getter for the center of the i-th finite element.

Parameters
  • its: index of the finite element
  • center: array holding the coordinates of the finite element center

Eigen::Matrix3Xd getCenters() const

Getter for the centers of the finite elements composing the molecular cavity.

Return
a matrix with the finite elements centers (dimensions 3 x number of finite elements)

void getAreas(double areas[]) const

Getter for the areas/weights of the finite elements.

Parameters
  • areas: array holding the weights/areas of the finite elements

void computeASC(const char *mep_name, const char *asc_name, int irrep) const

Computes ASC given a MEP and the desired irreducible representation.

Parameters
  • mep_name: label of the MEP surface function
  • asc_name: label of the ASC surface function
  • irrep: index of the desired irreducible representation The module uses the surface function concept to handle potentials and charges. Given labels for each, this function retrieves the MEP and computes the corresponding ASC.

void computeResponseASC(const char *mep_name, const char *asc_name, int irrep) const

Computes response ASC given a MEP and the desired irreducible representation.

Parameters
  • mep_name: label of the MEP surface function
  • asc_name: label of the ASC surface function
  • irrep: index of the desired irreducible representation If Nonequilibrium = True in the input, calculates a response ASC using the dynamic permittivity. Falls back to the solver with static permittivity otherwise.

double computePolarizationEnergy(const char *mep_name, const char *asc_name) const

Computes the polarization energy.

Return
the polarization energy This function calculates the dot product of the given MEP and ASC vectors.
Parameters
  • mep_name: label of the MEP surface function
  • asc_name: label of the ASC surface function

double getASCDipole(const char *asc_name, double dipole[]) const

Getter for the ASC dipole.

Return
the ASC dipole, i.e. { ^2}
Parameters
  • asc_name: label of the ASC surface function
  • dipole: the Cartesian components of the ASC dipole moment

void getSurfaceFunction(PCMSolverIndex size, double values[], const char *name) const

Retrieves data wrapped in a given surface function.

Parameters
  • size: the size of the surface function
  • values: the values wrapped in the surface function
  • name: label of the surface function

void setSurfaceFunction(PCMSolverIndex size, double values[], const char *name) const

Sets a surface function given data and label.

Parameters
  • size: the size of the surface function
  • values: the values to be wrapped in the surface function
  • name: label of the surface function

void printSurfaceFunction(const char *name) const

Prints surface function contents to host output.

Parameters
  • name: label of the surface function

void saveSurfaceFunctions() const

Dumps all currently saved surface functions to NumPy arrays in .npy files.

void saveSurfaceFunction(const char *name) const

Dumps a surface function to NumPy array in .npy file.

Note
The name parameter is the name of the NumPy array file without .npy extension
Parameters
  • name: label of the surface function

void loadSurfaceFunction(const char *name) const

Loads a surface function from a .npy file.

Note
The name parameter is the name of the NumPy array file without .npy extension
Parameters
  • name: label of the surface function

void printInfo() const

Prints citation and set up information.

void writeTimings() const

Writes timing results for the API.

Private Functions

void CTORBody()

Common implemenation for the CTOR-s

void initInput(int nr_nuclei, double charges[], double coordinates[], int symmetry_info[])

Initialize input_.

Parameters
  • nr_nuclei: number of atoms in the molecule
  • charges: atomic charges
  • coordinates: atomic coordinates
  • symmetry_info: molecular point group information

void initCavity()

Initialize cavity_

void initStaticSolver()

Initialize static solver K_0_

void initDynamicSolver()

Initialize dynamic solver K_d_

void mediumInfo(IGreensFunction *gf_i, IGreensFunction *gf_o) const

Collect info on medium

Private Members

Printer hostWriter_

Output redirect-or to host program output

Input input_

Input object

ICavity *cavity_

Cavity

ISolver *K_0_

Solver with static permittivity

ISolver *K_d_

Solver with dynamic permittivity

std::ostringstream infoStream_

PCMSolver set up information

bool hasDynamic_

Whether K_d_ was initialized

SurfaceFunctionMap functions_

SurfaceFunction map

class pcm::Input

A wrapper class for the Getkw Library C++ bindings.

An Input object is to be used as the unique point of access to user-provided input: input > parsed input (Python script) > Input object (contains all the input data) Definition of input parameters is to be done in the Python script and in this class. They must be specified as private data members with public accessor methods (get-ters). Most of the data members are anyway accessed through the input wrapping struct-s In general, no mutator methods (set-ters) should be needed, exceptions to this rule should be carefully considered.

Author
Roberto Di Remigio
Date
2013

Public Functions

Input()

Default constructor.

Input(const std::string &filename)

Constructor from human-readable input file name.

Input(const PCMInput &host_input)

Constructor from host input structs.

std::string units() const

Accessor methods.

Top-level section input

bool scaling() const

Cavity section input.

void molecule(const Molecule &m)

This method sets the molecule and the list of spheres.

Solvent solvent() const

Medium section input.

std::string providedBy() const

Keeps track of who did the parsing: the API or the host program.

CavityData cavityParams() const

Get-ters for input wrapping structs.

Private Functions

void reader(const PCMInput &host_input)

Read host data structures (host-side syntactic input parsing) into Input object. It provides access to a limited number of options only, basically the ones that can be filled into the cavityInput, solverInput and greenInput data structures. Lengths and areas are expected to be in Angstrom/Angstrom^2 and will hence be converted to au/au^2.

Note
Specification of the solvent by name overrides any input given through the greenInput data structure!
Warning
The cavity can only be built in the “Implicit” mode, i.e. by grabbing the coordinates for the sphere centers from the host program. Atomic coordinates are expected to be in au! The “Atoms” and “Explicit” methods are only available using the explicit parsing by our Python script of a separate input file.

void semanticCheck() const

Perform semantic input parsing aka sanity check

Private Members

std::string units_

Units of measure.

int CODATAyear_

Year of the CODATA set to be used.

std::string cavityType_

The type of cavity.

std::string cavFilename_

Filename for the .npz cavity restart file.

double area_

GePol cavity average element area.

bool scaling_

Whether the radii should be scaled by 1.2.

std::string radiiSet_

The set of radii to be used.

std::string radiiSetName_

Collects info on atomic radii set.

double minimalRadius_

Minimal radius of an added sphere.

std::string mode_

How the API should get the coordinates of the sphere centers.

std::vector<int> atoms_

List of selected atoms with custom radii.

std::vector<double> radii_

List of radii attached to the selected atoms.

std::vector<Sphere> spheres_

List of spheres for fully custom cavity generation.

Molecule molecule_

Molecule or atomic aggregate.

Solvent solvent_

The solvent for a vacuum/uniform dielectric run.

bool hasSolvent_

Whether the medium was initialized from a solvent object.

std::string solverType_

The solver type.

double correction_

Correction factor (C-PCM)

bool hermitivitize_

Whether the PCM matrix should be hermitivitized (collocation solvers)

bool isDynamic_

Whether the dynamic PCM matrix should be used.

double probeRadius_

Solvent probe radius.

std::string integratorType_

Type of integrator for the diagonal of the boundary integral operators.

double integratorScaling_

Scaling factor for the diagonal of the approximate collocation boundary integral operators

std::string greenInsideType_

The Green’s function type inside the cavity. It encodes the Green’s function type, derivative calculation strategy and dielectric profile: TYPE_DERIVATIVE_PROFILE

std::string greenOutsideType_

The Green’s function type outside the cavity It encodes the Green’s function type, derivative calculation strategy and dielectric profile: TYPE_DERIVATIVE_PROFILE

double epsilonInside_

Permittivity inside the cavity.

double epsilonStaticOutside_

Static permittivity outside the cavity.

double epsilonDynamicOutside_

Dynamic permittivity outside the cavity.

double epsilonReal_

Real part of the metal NP permittivity.

double epsilonImaginary_

Imaginary part of the metal NP permittivity.

std::vector<double> spherePosition_

Center of the spherical metal NP.

double sphereRadius_

Radius of the spherical metal NP.

double epsilonStatic1_

Diffuse interface: static permittivity inside the interface.

double epsilonDynamic1_

Diffuse interface: dynamic permittivity inside the interface.

double epsilonStatic2_

Diffuse interface: static permittivity outside the interface.

double epsilonDynamic2_

Diffuse interface: dynamic permittivity outside the interface.

double center_

Center of the diffuse interface.

double width_

Width of the diffuse interface.

int maxL_

Maximum angular momentum.

std::vector<double> origin_

Center of the dielectric sphere.

std::vector<double> geometry_

Molecular geometry.

bool MEPfromMolecule_

Whether to calculate the MEP from the molecular geometry.

ChargeDistribution multipoles_

Classical charge distribution of point multipoles.

std::string providedBy_

Who performed the syntactic input parsing.

Friends

std::ostream &operator<<(std::ostream &os, const Input &input)

Operators operator<<