Helper classes and functions

../_images/utils.png

Input

Warning

doxygenclass: Cannot find class “Input” in doxygen xml output for project “PCMSolver” from directory: xml

Sphere

struct Sphere

POD describing a sphere.

Author
Roberto Di Remigio
Date
2011, 2016

Atom

struct Atom

A POD describing an atom.

Author
Roberto Di Remigio
Date
2011, 2016

Molecule

Warning

doxygenclass: Cannot find class “Molecule” in doxygen xml output for project “PCMSolver” from directory: xml

Solvent

struct 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

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

Public Functions

Symmetry()

Default constructor sets up C1 point group.

Private Members

int nrGenerators_

Number of generators

int generators_[3]

Generators

int nrIrrep_

Number of irreps

Mathematical utilities

namespace pcm

PCMSolver, an API for the Polarizable Continuum Model Copyright (C) 2017 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/

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
  • 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
  • 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
  • value: the value to be checked
  • threshold: the threshold

bool numericalZero(double value)

Returns true if value is less than 1.0e-14

Parameters
  • value: the value to be checked

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
  • blockedMatrix: the result of packing fullMatrix
  • fullMatrix: the matrix to be packed
  • dimBlock: the dimension of the square blocks
  • 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
  • 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
  • R_: the rotation matrix
  • 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
  • point: where the function has to be evaluated
  • grid: holds points on grid where function is known
  • 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
  • point: where the function has to be evaluated
  • grid: holds points on grid where function is known
  • 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
  • matrix: Eigen object
  • fname: name of the file
Template Parameters
  • Derived: template parameters of the MatrixBase object

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.

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.

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