Helper classes and functions

../_images/utils.svg

Sphere

struct pcm::utils::Sphere

POD describing a sphere.

Author

Roberto Di Remigio

Date

2011, 2016

Public Functions

void scale(double scaling)

Scale sphere to other units.

Atom

struct pcm::utils::Atom

A POD describing an atom.

Author

Roberto Di Remigio

Date

2011, 2016

Public Members

double charge

Atomic charge

double mass

Atomic mass

double radius

Atomic radius

double radiusScaling

Scaling of the atomic radius

Eigen::Vector3d position

Position of the atom

std::string element

Name of the element

std::string symbol

Atomic symbol

ChargeDistribution

struct pcm::utils::ChargeDistribution

POD representing a classical charge distribution.

Author

Roberto Di Remigio

Date

2016

Public Members

Eigen::VectorXd monopoles

Monopoles

Eigen::Matrix3Xd monopolesSites

Monopoles sites

Eigen::Matrix3Xd dipoles

Dipoles

Eigen::Matrix3Xd dipolesSites

Dipoles sites

Eigen::VectorXd FQChi

FQ electronegativities

Eigen::VectorXd FQEta

FQ hardnesses

Eigen::Matrix3Xd FQSites

FQ sites

Molecule

class pcm::Molecule

Class representing a molecule or general aggregate of atoms.

This class is based on the similar class available in the Mints library of Psi4

Author

Roberto Di Remigio

Date

2014

Unnamed Group

Molecule &operator=(const Molecule &other)

Operators Assignment operator.

Public Functions

Molecule()

Default constructor Initialize a dummy molecule, e.g. as placeholder, see ICavity.cpp loadCavity method.

Molecule(int nat, const Eigen::VectorXd &chg, const Eigen::VectorXd &masses, const Eigen::Matrix3Xd &geo, const std::vector<Atom> &at, const std::vector<Sphere> &sph)

Constructor from full molecular data.

This initializes the molecule in C1 symmetry

Parameters
  • [in] nat: number of atoms

  • [in] chg: vector of atomic charges

  • [in] masses: vector of atomic masses

  • [in] geo: molecular geometry (format nat*3)

  • [in] at: vector of Atom objects

  • [in] sph: vector of Sphere objects

Molecule(int nat, const Eigen::VectorXd &chg, const Eigen::VectorXd &masses, const Eigen::Matrix3Xd &geo, const std::vector<Atom> &at, const std::vector<Sphere> &sph, int nr_gen, std::array<int, 3> gens)

Constructor from full molecular data, plus number of generators and generators.

This initializes the molecule in the symmetry prescribed by nr_gen and gen. See documentation of the

Symmetry object for the conventions.
Parameters
  • [in] nat: number of atoms

  • [in] chg: vector of atomic charges

  • [in] masses: vector of atomic masses

  • [in] geo: molecular geometry (format nat*3)

  • [in] at: vector of Atom objects

  • [in] sph: vector of Sphere objects

  • [in] nr_gen: number of molecular point group generators

  • [in] gen: molecular point group generators

Molecule(int nat, const Eigen::VectorXd &chg, const Eigen::VectorXd &masses, const Eigen::Matrix3Xd &geo, const std::vector<Atom> &at, const std::vector<Sphere> &sph, const Symmetry &pg)

Constructor from full molecular data and point group.

This initializes the molecule in the symmetry prescribed by pg.

Parameters
  • [in] nat: number of atoms

  • [in] chg: vector of atomic charges

  • [in] masses: vector of atomic masses

  • [in] geo: molecular geometry (format nat*3)

  • [in] at: vector of Atom objects

  • [in] sph: vector of Sphere objects

  • [in] pg: the molecular point group (a Symmetry object)

Molecule(const std::vector<Sphere> &sph)

Constructor from list of spheres.

Molecule is treated as an aggregate of spheres. We do not have information on the atomic species involved in the aggregate. Charges are set to 1.0; masses are set based on the radii; geometry is set from the list of spheres. All the atoms are dummy atoms. The point group is C1.

Warning

This constructor is to be used exclusively when initializing the Molecule in EXPLICIT mode, i.e. when the user specifies explicitly spheres centers and radii.

Parameters
  • [in] sph: list of spheres

Molecule(const Molecule &other)

Copy constructor.

void translate(const Eigen::Vector3d &translationVector)

Given a vector, carries out translation of the molecule.

Parameters
  • translationVector: The translation vector.

void moveToCOM()

Performs translation to the Center of Mass Frame.

void rotate(const Eigen::Matrix3d &rotationMatrix)

Given a matrix, carries out rotation of the molecule.

Parameters
  • rotationMatrix: The matrix representing the rotation.

void moveToPAF()

Performs rotation to the Principal Axes Frame.

Private Members

size_t nAtoms_

The number of atoms in the molecule.

Eigen::VectorXd charges_

A vector of dimension (# atoms) containing the charges.

Eigen::VectorXd masses_

A vector of dimension (# atoms) containing the masses.

Eigen::Matrix3Xd geometry_

Molecular geometry, in cartesian coordinates. The dimensions are (# atoms * 3) Units are Bohr.

std::vector<Atom> atoms_

A container for all the atoms composing the molecule.

std::vector<Sphere> spheres_

A container for the spheres composing the molecule.

rotorType rotor_

The molecular rotor type.

Symmetry pointGroup_

The molecular point group.

Solvent

struct pcm::utils::Solvent

POD describing a solvent.

A

Solvent object contains all the solvent-related experimental data needed to set up the Green’s functions and the non-electrostatic terms calculations.
Author

Roberto Di Remigio

Date

2011, 2016

Public Members

std::string name

Solvent name

double epsStatic

Static permittivity, in AU

double epsDynamic

Optical permittivity, in AU

double probeRadius

Radius of the spherical probe mimicking the solvent, in Angstrom

Symmetry

class Symmetry

Contains very basic info about symmetry (only Abelian groups)

Just a wrapper around a vector containing the generators of the group

Author

Roberto Di Remigio

Date

2014

Private Members

int nrGenerators_ = {0}

Number of generators

std::array<int, 3> generators_ = {0}

Generators

int nrIrrep_ = {1}

Number of irreps

Mathematical utilities

namespace pcm

PCMSolver, an API for the Polarizable Continuum Model Copyright (C) 2020 Roberto Di Remigio, Luca Frediani and collaborators.

This file is part of PCMSolver.

PCMSolver is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

PCMSolver is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.

You should have received a copy of the GNU Lesser General Public License along with PCMSolver. If not, see http://www.gnu.org/licenses/.

For information on the complete list of contributors to the PCMSolver API, see: http://pcmsolver.readthedocs.io/

PCMSolver, an API for the Polarizable Continuum Model Copyright (C) 2020 Roberto Di Remigio, Luca Frediani and contributors.

This file is part of PCMSolver.

PCMSolver is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

PCMSolver is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.

You should have received a copy of the GNU Lesser General Public License along with PCMSolver. If not, see http://www.gnu.org/licenses/.

For information on the complete list of contributors to the PCMSolver API, see: http://pcmsolver.readthedocs.io/

namespace utils

Functions

template<size_t nBits>
int parity(std::bitset<nBits> bitrep)

Calculate the parity of the bitset as defined by: bitrep[0] XOR bitrep[1] XOR … XOR bitrep[nBits-1]

Parameters
  • [in] bitrep: a bitset

Template Parameters
  • nBits: lenght of the input bitset

double parity(unsigned int i)

Returns parity of input integer. The parity is defined as the result of using XOR on the bitrep of the given integer. For example: 2 -> 010 -> 0^1^0 = 1 -> -1.0 6 -> 110 -> 1^1^0 = 0 -> 1.0

Parameters
  • [in] i: an integer, usually an index for an irrep or a symmetry operation

It can also be interpreted as the action of a given operation on the Cartesian axes: zyx 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

bool isZero(double value, double threshold)

Returns true if value is less or equal to threshold

Parameters
  • [in] value: the value to be checked

  • [in] threshold: the threshold

bool numericalZero(double value)

Returns true if value is less than 1.0e-14

Parameters
  • [in] value: the value to be checked

template<typename T>
int sign(T val)

This function implements the signum function and returns the sign of the passed value: -1, 0 or 1

Parameters
  • [in] val: value whose sign should be determined

Template Parameters
  • T: of the parameter val

void symmetryBlocking(Eigen::MatrixXd &matrix, PCMSolverIndex cavitySize, PCMSolverIndex ntsirr, int nr_irrep)
void symmetryPacking(std::vector<Eigen::MatrixXd> &blockedMatrix, const Eigen::MatrixXd &fullMatrix, int dimBlock, int nrBlocks)

Parameters
  • [out] blockedMatrix: the result of packing fullMatrix

  • [in] fullMatrix: the matrix to be packed

  • [in] dimBlock: the dimension of the square blocks

  • [in] nrBlocks: the number of square blocks

template<typename Derived>
void hermitivitize(Eigen::MatrixBase<Derived> &obj_)

Given obj_ returns 0.5 * (obj_ + obj_^dagger)

Note

We check if a matrix or vector was given, since in the latter case we only want the complex conjugation operation to happen.

Parameters
  • [out] obj_: the Eigen object to be hermitivitized

Template Parameters
  • Derived: the numeric type of obj_ elements

void eulerRotation(Eigen::Matrix3d &R_, const Eigen::Vector3d &eulerAngles_)

Build rotation matrix between two reference frames given the Euler angles.

We assume the convention

\( R = Z_3 X_2 Z_1 \) for the ordering of the extrinsic elemental rotations (see http://en.wikipedia.org/wiki/Euler_angles) The Euler angles are given in the order \( \phi, \theta, \psi \). If we write \( c_i, s_i \,\, i = 1, 3 \) for their cosines and sines the rotation matrix will be:
\[\begin{split} R = \begin{pmatrix} c_1c_3 - s_1c_2s_3 & -s_1c_3 - c_1c_2s_3 & s_2s_3 \\ c_1s_3 + s_1c_2c_3 & -s_1s_3 + c_1c_2c_3 & -s_2c_3 \\ s_1s_2 & c_1s_2 & c_2 \end{pmatrix} \end{split}\]
Eigen’s geometry module is used to calculate the rotation matrix
Parameters
  • [out] R_: the rotation matrix

  • [in] eulerAngles_: the Euler angles, in degrees, describing the rotation

double linearInterpolation(const double point, const std::vector<double> &grid, const std::vector<double> &function)

Return value of function defined on grid at an arbitrary point.

This function finds the nearest values for the given point and performs a linear interpolation.

Warning

This function assumes that grid has already been sorted!

Parameters
  • [in] point: where the function has to be evaluated

  • [in] grid: holds points on grid where function is known

  • [in] function: holds known function values

double splineInterpolation(const double point, const std::vector<double> &grid, const std::vector<double> &function)

Return value of function defined on grid at an arbitrary point.

This function finds the nearest values for the given point and performs a cubic spline interpolation.

Warning

This function assumes that grid has already been sorted!

Parameters
  • [in] point: where the function has to be evaluated

  • [in] grid: holds points on grid where function is known

  • [in] function: holds known function values

template<typename Derived>
void print_eigen_matrix(const Eigen::MatrixBase<Derived> &matrix, const std::string &fname)

Prints Eigen object (matrix or vector) to file.

Note

This is for debugging only, the format is in fact rather ugly. Row index Column index Matrix entry 0 0 0.0000

Parameters
  • [in] matrix: Eigen object

  • [in] fname: name of the file

Template Parameters
  • Derived: template parameters of the MatrixBase object

Eigen::MatrixXd prune_zero_columns(const Eigen::MatrixXd &incoming, const Eigen::Matrix<bool, 1, Eigen::Dynamic> &filter)

Prune zero columns from matrix.

Outgoing matrix has the same number of rows as the incoming.

Parameters
  • [in] incoming: Matrix to be pruned

  • [in] filter: indexing array for pruning

Eigen::VectorXd prune_vector(const Eigen::VectorXd &incoming, const Eigen::Matrix<bool, 1, Eigen::Dynamic> &filter)

Prune zero elements from Vector.

Parameters
  • [in] incoming: VectorXd to be pruned

  • [in] filter: indexing array for pruning

namespace cnpy
namespace custom

Custom overloads for cnpy load and save functions

Functions

template<typename Scalar, int Rows, int Cols>
void npy_save(const std::string &fname, const Eigen::Matrix<Scalar, Rows, Cols> &obj)

Save Eigen object to NumPy array file.

Parameters
  • fname: name of the NumPy array file

  • obj: Eigen object to be saved, either a matrix or a vector

Template Parameters
  • Scalar: the data type of the matrix to be returned. Default is double

  • Rows: number of rows in the Eigen object. Default is dynamic e

  • Cols: number of columns in the Eigen object. Default is dynamic

template<typename Scalar, int Rows, int Cols>
void npz_save(const std::string &fname, const std::string &name, const Eigen::Matrix<Scalar, Rows, Cols> &obj, bool overwrite = false)

Save Eigen object to a compressed NumPy file.

Parameters
  • fname: name of the compressed NumPy file

  • name: tag for the given object in the compressed NumPy file

  • obj: Eigen object to be saved, either a matrix or a vector

  • overwrite: if file exists, overwrite. Appends by default.

Template Parameters
  • Scalar: the data type of the matrix to be returned. Default is double

  • Rows: number of rows in the Eigen object. Default is dynamic

  • Cols: number of columns in the Eigen object. Default is dynamic

template<typename Scalar>
Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> npy_to_eigen(const NpyArray &npy_array)

Load NpyArray object into Eigen object.

Todo:

Extend to read in also data in row-major (C) storage order

Return

An Eigen object (matrix or vector) with the data

Warning

We check that the rank of the object read is not more than 2 Eigen cannot handle general tensors.

Parameters
  • npy_array: the NpyArray object

Template Parameters
  • Scalar: the data type of the matrix to be returned. Default is double

template<typename Scalar>
Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> npy_load(const std::string &fname)

Load NumPy array file into Eigen object.

Todo:

Extend to read in also data in row-major (C) storage order

Return

An Eigen object (matrix or vector) with the data

Parameters
  • fname: name of the NumPy array file

Template Parameters
  • Scalar: the data type of the matrix to be returned. Default is double