Namespace Scine

namespace Scine

This header file contains functions that allow for common notation for common things that can be done at a different degree of derivatives.

This header contains alias definitions defining which classes to use for the different degrees of derivatives.

namespace Utils

Typedefs

using Bohr = StrongType<double, struct BohrType>
using Angstrom = StrongType<double, struct AngstromType>
using Hartree = StrongType<double, struct HartreeType>
using kJPerMol = StrongType<double, struct kJPerMolType>
using kCalPerMol = StrongType<double, struct kCalPerMolType>
using eV = StrongType<double, struct eVType>
using SecondDerivative = AutomaticDifferentiation::Second3D
using Displacement = Eigen::RowVector3d
using DisplacementCollection = Eigen::Matrix<double, Eigen::Dynamic, 3, Eigen::RowMajor>
using ElementTypeCollection = std::vector<ElementType>
using Gradient = Eigen::RowVector3d
using GradientCollection = Eigen::Matrix<double, Eigen::Dynamic, 3, Eigen::RowMajor>
using HessianMatrix = Eigen::MatrixXd
using Dipole = Eigen::RowVector3d
using DipoleGradient = Eigen::Matrix<double, Eigen::Dynamic, 3, Eigen::RowMajor>
using Position = Eigen::RowVector3d
using PositionCollection = Eigen::Matrix<double, Eigen::Dynamic, 3, Eigen::RowMajor>

Enums

enum Property

The properties contained are assigned a bit. This can be switched on or off to flag the presence or absence of the property.

Values:

Energy = (1 << 0)
Gradients = (1 << 1)
Hessian = (1 << 2)
Dipole = (1 << 3)
DipoleGradient = (1 << 4)
DipoleMatrixAO = (1 << 5)
DipoleMatrixMO = (1 << 6)
OneElectronMatrix = (1 << 7)
TwoElectronMatrix = (1 << 8)
BondOrderMatrix = (1 << 9)
enum StateSize

enum class defining the desired size of the State to save.

Values:

minimal
regular
extensive
enum ElementType

Enum class defining all elements.

Values:

none = 0
H
He
Li
Be
B
C
N
O
F
Ne
Na
Mg
Al
Si
P
S
Cl
Ar
K
Ca
Sc
Ti
V
Cr
Mn
Fe
Co
Ni
Cu
Zn
Ga
Ge
As
Se
Br
Kr
Rb
Sr
Y
Zr
Nb
Mo
Tc
Ru
Rh
Pd
Ag
Cd
In
Sn
Sb
Te
I
Xe
Cs
Ba
La
Ce
Pr
Nd
Pm
Sm
Eu
Gd
Tb
Dy
Ho
Er
Tm
Yb
Lu
Hf
Ta
W
Re
Os
Ir
Pt
Au
Hg
Tl
Pb
Bi
Po
At
Rn
Fr
Ra
Ac
Th
Pa
U
Np
Pu
Am
Cm
Bk
Cf
Es
Fm
Md
No
Lr
Rf
Db
Sg
Bh
Hs
Mt
Ds
Rg
Cn
enum derivOrder

Values:

zero
one
two
enum derivativeType

Values:

none
first
second_atomic
second_full
enum scf_mixer_t

Values:

none
fock_diis
ediis
ediis_diis
fock_simple
charge_simple

Functions

constexpr Property operator|(Property v1, Property v2)

Returns a Property object that is the superset of the two properties given as argument.

Allow combination of properties.

constexpr bool operator&(Property v1, Property v2)

Returns a Property object that is the subset of the two properties given as argument.

Allow to check if there is a flag overlap.

constexpr double toAngstrom(Bohr b)

Conversion from Bohr to Angstrom.

Return

constexpr double Returns the length in Angstrom.

Parameters
  • b: Length in/as Bohr.

constexpr double toBohr(Angstrom a)

Conversion from Angstrom to Bohr.

Return

constexpr double Returns the length in Bohr.

Parameters
  • a: Length in/as Angstrom.

constexpr double toKJPerMol(Hartree h)

Conversion from Hartree to kJ/mol.

Return

constexpr double Returns the energy in kJ/mol.

Parameters
  • h: Energy in/as Hartree.

constexpr double toHartree(kJPerMol kjpm)

Conversion from kJ/mol to Hartree.

Return

constexpr double Returns the energy in Hartree.

Parameters
  • kjpm: Energy in/as kJ/mol.

constexpr double toKCalPerMol(Hartree h)

Conversion from Hartree to kcal/mol.

Return

constexpr double Returns the energy in kcal/mol.

Parameters
  • h: Energy in/as Hartree.

constexpr double toHartree(kCalPerMol kcpm)

Conversion from kcal/mol to Hartree.

Return

constexpr double Returns the energy in Hartree.

Parameters
  • kcpm: Energy in/as kcal/mol.

std::string generateChemicalFormula(const ElementTypeCollection &elements, const std::string &numberPrefix = "", const std::string &numberPostfix = "")

This function returns the elemental composition of a compound based on its ElementTypeCollection.

The convention is the Hill Order System. It starts with C, then H, then all the other elements in alphabetical order. Returns “(empty)” if elements is empty.

Return

std::string The final string.

Parameters
  • elements: The collection of elements to be transcribed.

  • numberPrefix: Characters to put before each number for rich text formats (Latex etc).

  • numberPostfix: Characters to put after each number for rich text formats (Latex etc).

std::string singleElementPartOfFormula(std::string symbol, int numberOccurrences, const std::string &numberPrefix, const std::string &numberPostfix)

This function transcribes a single element into a string.

Return

std::string The final string.

Parameters
  • symbol: The element symbol

  • numberOccurrences: The number of times the element is present in a given structure.

  • numberPrefix: Characters to put before each number for rich text formats (Latex etc).

  • numberPostfix: Characters to put after each number for rich text formats (Latex etc).

template<class MatrixType>
std::pair<int, int> getMatrixDimensions(const MatrixType &matrix)
template<>
std::pair<int, int> getMatrixDimensions(const Eigen::MatrixXd &matrix)
template<>
std::pair<int, int> getMatrixDimensions(const Eigen::Matrix<double, Eigen::Dynamic, 3, Eigen::RowMajor> &matrix)
template<class MatrixType>
double getElement(int row, int col, const MatrixType &matrix)
template<>
double getElement(int row, int col, const Eigen::MatrixXd &matrix)
template<>
double getElement(int row, int col, const Eigen::Matrix<double, Eigen::Dynamic, 3, Eigen::RowMajor> &matrix)
template<>
double getElement(int row, int col, const NormalModesContainer &matrix)
template<class MatrixType>
void matrixPrettyPrint(std::ostream &out, const MatrixType &matrix)
void printElement(std::ostream &out, const Eigen::SparseMatrix<double>::InnerIterator &it)
void matrixPrettyPrint(std::ostream &out, const Eigen::SparseMatrix<double> matrix, double threshold)
void matrixPrettyPrint(std::ostream &out, const NormalModesContainer &container, const ElementTypeCollection &elementTypes)
static void nodeToSettings(Utils::Settings &settings, YAML::Node node, bool allowSuperfluous = false)

Parses data from a YAML::Node (yaml-cpp) into the value collection of a settings object.

Parameters
  • settings: The settings.

  • node: The read yaml node.

  • allowSuperfluous: It true, the function will ignore kys in the node that are not present in the settings object. If false, the function will throw an error if such a case is encountered.

template<typename T>
bool checkMethodType(SinglePointMethod *m)
template<>
MatrixWithDerivatives::Matrix0 &get<Utils::derivOrder::zero>()
template<>
MatrixWithDerivatives::Matrix1 &get<Utils::derivOrder::one>()
template<>
MatrixWithDerivatives::Matrix2 &get<Utils::derivOrder::two>()
template<>
const MatrixWithDerivatives::Matrix0 &get<Utils::derivOrder::zero>() const
template<>
const MatrixWithDerivatives::Matrix1 &get<Utils::derivOrder::one>() const
template<>
const MatrixWithDerivatives::Matrix2 &get<Utils::derivOrder::two>() const
class BondDetector
#include <BondDetector.h>

Detecting bonds from 3D structures using covalent radii.

The radii are stored in Scine::Utils::BondDetectorRadii .

References:

  • DOI: 10.1002/jcc.540120716

  • DOI: 10.1186/1758-2946-4-26

  • DOI: 10.1002/jcc.24309

class BondDetectorRadii
#include <BondDetectorRadii.h>

Holds the covalent radii required in BondDetector.

class BondOrderCollection
#include <BondOrderCollection.h>

Class defining the bond orders between atoms of some molecular system.

class VanDerWaalsBondDetector
#include <VanDerWaalsBondDetector.h>

Implements the function to find bonds from a set of elements and positions, using the van der Waals radii.

class PropertyList
#include <PropertyList.h>

This class defines a list of properties that can be calculated in a single-point calculation.

class PropertyNotPresentException : public exception
#include <Results.h>

Exception thrown if a property is requested and it is not stored.

class Results
#include <Results.h>

Class for the properties obtained in a single-point calculation.

To obtain the properties:

  • ”get” methods return references to the results (also allows simple copy)

  • ”take” methods move the results (which will not be present in Results afterwards).

Implementation note:

  • New voices can be added to the Result. Minimally a set, a get and a has method must be present.

  • Use pImpl idiom to hide boost::optional dependency

class State
#include <State.h>

Base class for the implementation of a generic State.

A state should be viewed as a checkpoint. Implementation note: If a class having an interface needs to save/load a state, then it has to have a polymorphic pointer to StatesHandler in the interface and then populate it with a derived class, specific to the class. The same is true with the polymorphic State class.

Subclassed by Scine::Utils::ExternalQC::OrcaState

class EmptyStatesHandlerContainer : public exception
#include <StatesHandler.h>

This exception is thrown if an empty states handler is popped.

class StatesHandler
#include <StatesHandler.h>

Base class for the implementation of a generic saving and loading capability.

Implementation note: If a class having an interface needs to save/load a state, then it has to populate the polymorphic pointers to the State class with a derived class, specific to the class.

Subclassed by Scine::Utils::ExternalQC::OrcaStatesHandler

class OrcaCalculatorSettings : public Scine::Utils::Settings
#include <OrcaCalculatorSettings.h>

Settings for ORCA calculations.

class HessianUtilities
#include <HessianUtilities.h>

A utility class for Hessians allowing easier access to eigenvalues and eigenvectors of transformed versions.

class NormalMode
#include <NormalMode.h>

Container for a single normal mode.

class NormalModeAnalyzer
#include <NormalModeAnalyzer.h>

Class to calculate the normal modes of a molecule from its hessian matrix.

class NormalModesContainer
#include <NormalModesContainer.h>

Class that holds the and manages the normal modes.

class NumericalHessianCalculator
#include <NumericalHessianCalculator.h>

This class calculates the Hessian and, optionally, the dipole gradient (semi-)numerically.

Useful if analytic second derivatives are not available.

calculateFromEnergyDifferences uses the following, where D = delta / 2: d^2E/dx^2 = 1 / delta^2 * (E(x-2D) - 2 * E(x) + E(x+2D)) d^2E/dxdy = 1 / delta^2 * (E(x+D,y+D) - E(x+D,y-D) - E(x-D,y+D) + E(x-D,y-D))

calculateFromGradientDifferences uses the following, where D = delta / 2: d^2E/dxdy = 1 / (2*delta) * (g_j(i+D,j) - g_j(i-D,j) + g_i(i,j+D) - g_i(i,j-D))

The second formulation is more stable numerically and is used as default.

In order to calculate the dipole gradient, from each displacement the dipole is calculated. dmu/dx_i = 1 / (2*delta) * (mu_x(x_i + delta) - mu_x(x_i - delta)) Only 3N calculations are needed to fill the 3Nx3 dipole gradient matrix, as the dipole comes as a x,y,z vector at each single point calculation.

class Atom
#include <Atom.h>

A basic atom in the chemical sense.

class AtomCollection
#include <AtomCollection.h>

A Collection of Atoms.

Has the same functionality as a std::vector<Atom>, but is implemented as a composition of a ElementTypeCollection and a PositionCollection.

class ElementSymbolNotFound : public runtime_error
#include <ElementInfo.h>

A runtime error specific to an element not being found.

class ElementInfo
#include <ElementInfo.h>

Provides information about elements, such as mass, van-der-Waals radius, etc.

This class only wraps around the actual data and their handling. For the underlying data see ElementInfoData.h and ElementInfoData.cpp.

class StructuralCompletion
#include <StructuralCompletion.h>

Contains function to generate positions from other positions.

F.i. to be applied to add hydrogen atoms to a carbon chain.

template<class OptimizerType>
class AFIROptimizer : public Scine::Utils::AFIROptimizerBase
#include <AFIROptimizer.h>

A version of the GeometryOptimizer that optimizes the underlying structure while applying an additional artificial force for the induction of reactions.

For more details about AFIR see the base class Utils::AFIROptimizerBase .

Template Parameters
  • OptimizerType: Expects any of the Optimizer classes. Note that some special optimizers may not yet be supported or may need additional specialization.

class AFIROptimizerBase
#include <AFIROptimizerBase.h>

The base class for all AFIR optimizers.

The main purpose of this base class is to hide the template parameter(s) of the derived class(es).

This optimizer is a version of the GeometryOptimizer that optimizes the underlying structure while applying an additional artificial force for the induction of reactions.

The artificial force induced reaction (AFIR) optimizer is based on the works by Maeda et. al. The implementation was done following this reference: J Comput Chem., 2018, 39, 233 [DOI: 10.1002/jcc.25106] https://www.ncbi.nlm.nih.gov/pmc/articles/PMC5765425/ (accessed 05.04.2019)

The original references are: J. Chem. Phys., 2010, 132, 241102 J. Chem. Theory Comput., 2011, 7, 2335

Subclassed by Scine::Utils::AFIROptimizer< OptimizerType >

template<class OptimizerType, class ConvergenceCheckType>
class AFIROptimizerSettings : public Scine::Utils::Settings
#include <AFIROptimizerSettings.h>

Settings for an AFIROptimizer.

Uses template arguments in order to automatically include the settings of underlying objects into the given settings.

Template Parameters

class GeometryOptimizerBase
#include <GeometryOptimizer.h>

The base class for all Geometry optimizers.

The purpose of the geometry optimizers is to wrap technical details needed for the actual optimization into a class that is easily used and has a small set of data needed in its constructor, is mainly configured via its settings, and then exposes the optimization through a simple function accepting a geometry (AtomCollection)

The main purpose of this base class is to hide the template parameter(s) of the derived class(es).

Subclassed by Scine::Utils::GeometryOptimizer< OptimizerType >

template<class OptimizerType, class ConvergenceCheckType>
class GeometryOptimizerSettings : public Scine::Utils::Settings
#include <GeometryOptimizer.h>

Settings for a GeometryOptimizer.

Uses template arguments in order to automatically include the settings of underlying objects into the given settings.

Template Parameters

template<class OptimizerType>
class GeometryOptimizer : public Scine::Utils::GeometryOptimizerBase
#include <GeometryOptimizer.h>

Basically just a templated version of the base class GeometryOptimizerBase, where the template defines the actual optimizer used in the geometry optimization.

Template Parameters
  • OptimizerType: Expects any of the Optimizer classes. Note that some special optimizers may not yet be supported or may need additional specialization.

class IRCOptimizerBase
#include <IRCOptimizer.h>

The base class for all IRC optimizers.

The main purpose of this base class is to hide the template parameter(s) of the derived class(es).

Subclassed by Scine::Utils::IRCOptimizer< OptimizerType >

template<class OptimizerType, class ConvergenceCheckType>
class IRCOptimizerSettings : public Scine::Utils::Settings
#include <IRCOptimizer.h>

Settings for an IRCOptimizer.

Uses template arguments in order to automatically include the settings of underlying objects into the given settings.

Template Parameters

template<class OptimizerType>
class IRCOptimizer : public Scine::Utils::IRCOptimizerBase
#include <IRCOptimizer.h>

A version of the GeometryOptimizer that optimizes along an internal reaction coordinate (IRC).

This optimizer mass-weights the actual gradient in order to optimize in the mass-weighted coordinate system.

Template Parameters
  • OptimizerType: Expects any of the Optimizer classes. Note that some special optimizers may not yet be supported or may need additional specialization.

class ChemicalFileHandler
#include <ChemicalFileHandler.h>

Handles the reading and writing of chemical file formats.

Matches suffixes to chemical file formats and then dispatches the calls to appropriate streaming I/O handlers. Catches I/O errors.

class FormattedStreamHandler
#include <FormattedStreamHandler.h>

The interface for all classes handling formatted streaming IO.

Subclassed by Scine::Utils::MOLStreamHandler, Scine::Utils::OpenBabelStreamHandler, Scine::Utils::XYZStreamHandler

class MOLStreamHandler : public Scine::Utils::FormattedStreamHandler
#include <MOLStreamHandler.h>

Handler for MOL stream IO.

This class implements the FormattedStreamHandler interface for use with Core’s module management. The core stream IO functions are available as static functions, too.

class OpenBabelStreamHandler : public Scine::Utils::FormattedStreamHandler
#include <OpenBabelStreamHandler.h>

If obabel is in the PATH, supports IO to range of OpenBabel’s file formats.

Note

Concerning the GPLv2 OpenBabel license, this StreamHandler very much communicates ‘at arms length’ with OpenBabel (https://www.gnu.org/licenses/old-licenses/gpl-2.0-faq.en.html#TOCGPLInProprietarySystem) by merely calling it via system process-level interaction, if present.

class XYZStreamHandler : public Scine::Utils::FormattedStreamHandler
#include <XYZStreamHandler.h>

Handler for XYZ stream IO.

This class implements the FormattedStreamHandler interface for use with Core’s module management. The core stream IO functions are available as static functions, too.

Note

Since the XYZ file format does not include bond information in any part, the read/write pair involving BondOrderCollections just throws in all cases.

class FilesystemHelpers
#include <FilesystemHelpers.h>

This class contains utility functions for dealing with files and directories.

Hides the use of Boost::filesystem.

class IOException : public runtime_error
#include <IOExceptions.h>

Base class for exceptions related to input-output.

Subclassed by Scine::Utils::InvalidAtomCountSpecification, Scine::Utils::InvalidAtomSpecification, Scine::Utils::InvalidFile

class InvalidFile : public Scine::Utils::IOException
#include <IOExceptions.h>

Exception when file cannot be opened.

class InvalidAtomCountSpecification : public Scine::Utils::IOException
#include <IOExceptions.h>

Exception thrown when the first line of a xyz input does not contain the number of atoms.

class InvalidAtomSpecification : public Scine::Utils::IOException
#include <IOExceptions.h>

Exception thrown when a line for an atom of a xyz input is invalid.

class Log
#include <Logger.h>

Static class for logging.

By default, logging is disabled. Logging can be started with startConsoleLogging() or startFileLogging(). It is possible to have several files for logging.

class MolecularTrajectoryIO
#include <MolecularTrajectoryIO.h>

Class for input and output of MolecularTrajectory classes.

class NativeFilenames
#include <NativeFilenames.h>

This class contains utility functions for storing file paths in strings, and cares for cross-platform issues. Hides the use of Boost::filesystem.

class YAMLParsingException : public runtime_error
#include <Yaml.h>

A custom exception for YAML parsing errors.

class AtomicSecondDerivativeCollection
#include <AtomicSecondDerivativeCollection.h>

Container class for second derivatives.

Mainly an override of std::vector<SecondDerivative>. Generally used for storing the gradients of a molecular structure.

class FullSecondDerivativeCollection
#include <FullSecondDerivativeCollection.h>

Container class for full second derivatives (Hessian matrix + first derivatives).

class QuaternionFit
#include <QuaternionFit.h>

Compares two sets of positions and calculates translation and rotation to go from one to the other.

References:

class AvailableMethods
#include <AvailableMethods.h>

This class defines the methods that are available for single-point calculations.

class MethodIdentifier
#include <MethodIdentifier.h>

Identifier for some method type. Fulfills the same role as an enum that can be changed at runtime.

class MethodInitializer
#include <MethodInitializer.h>

MethodInitializer is an abstract base class for the method initializers.

Date

20.11.2015

class MethodSpecifier
#include <MethodSpecifier.h>

This class defines the properties of a method and implements a factory to create instances of it.

class ConvergenceChecker
#include <ConvergenceChecker.h>

Class for a density matrix based convergence checker. It looks at the RMS change in density matrix elements.

class DensityMatrixGuessCalculator
#include <DensityMatrixGuessCalculator.h>

Interface for the calculation of the density matrix guess in SCF calculations.

class ElectronicContributionCalculator
#include <ElectronicContributionCalculator.h>

Interface for the computation of the Fock matrix. TODO: Separate for non-SCF and SCF specializations?

class ElectronicEnergyCalculator
#include <ElectronicEnergyCalculator.h>

Interface for the calculation of the electronic energy.

class LCAOMethod : public Scine::Utils::SinglePointMethod
#include <LCAOMethod.h>

Base class describing a method based on a LCAO ansatz.

Subclassed by Scine::Utils::SCFMethod

class OverlapCalculator
#include <OverlapCalculator.h>

Abstract class for the calculation of the overlap matrix.

class RepulsionCalculator
#include <RepulsionCalculator.h>

Interface for the calculation of the repulsion energy of LCAO methods.

class ScfConvergenceAccelerator
#include <ScfConvergenceAccelerator.h>

This class sets up the convergence acceleration for a SCF method.

class SCFMethod : public Scine::Utils::LCAOMethod
#include <SCFMethod.h>

Base class for self-consistent-field methods.

class SCFModifier
#include <SCFModifier.h>

Base class for SCF modifiers. SCF modifiers allow to act in between the different steps of a SCF calculations. The functions it defines are anchors and are called after the corresponding steps of a SCF calculation.

See

SCFMethod

Subclassed by Scine::Utils::Charge_Simple, Scine::Utils::EdiisDiisModifier, Scine::Utils::EdiisModifier, Scine::Utils::Fock_Simple, Scine::Utils::FockDiisModifier

class SinglePointMethod
#include <SinglePointMethod.h>

This class is the base class of any method to be used in Real-time quantum chemistry. It contain declarations for variables methods concerning energy and gradients.

Subclassed by Scine::Utils::LCAOMethod

class StructureDependentInitializer
#include <StructureDependentInitializer.h>

Interface for method-specific settings which must be updated f.i. when the molecular structure changes.

struct ScfOptions
#include <ScfOptions.h>

SCF options for the UniversalSettings syntax.

class DiisError
#include <DiisError.h>

Class to handle the error vectors used in DIIS.

class ECQPP
#include <ECQPP.h>

Equality Constrained Quadratic Programming Problem solver for EDIIS. This class solves the optimization problem for an EDIIS iteration. It minimizes the expression E*c - 0.5 c^T*B*c with the constraint that sum(c_i) = 1. To do so, it solves the above equation for all subsets of c, rejects inadmissible solutions (c_i < 0) and keeps the best one.

class Ediis
#include <Ediis.h>

Class performing the calculation of a Fock matrix based on the EDIIS algorithm.

class EdiisCoefficientOptimizer
#include <EdiisCoefficientOptimizer.h>

Calculates optimal coefficients to use in EDIIS. It solves the optimization problem given the B matrix and energies and calculates the coefficients with the Reduced Gradient Method.

class EdiisDiisModifier : public Scine::Utils::SCFModifier
#include <EdiisDiisModifier.h>

Class to use the EDIIS + DIIS combination as in Gaussian. Garza, Scuseria, J Chem Phys 137 (2012) 054110: “In the default EDIIS + DIIS procedure in Gaussian, EDIIS is used when the largest DIIS error (errMax) is greater than 10^-1 a.u. but DIIS is employed when this error goes below 10^-4 a.u. In between these values, the EDIIS and DIIS coefficients are weighted such that c = 10(errorMax) C(EDIIS) + (1-10(errorMax))c(DIIS); however, if the error of the last cycle is 10% greater than the minimum error, pure EDIIS is used.”

class EdiisModifier : public Scine::Utils::SCFModifier
#include <EdiisModifier.h>

SCFModifier for the EDIIS algorithm. This class does not perform directly the EDIIS algorithm, it is like an interface between the SCF method and the class for the actual EDIIS algorithm.

class FockDiis
#include <FockDiis.h>

Class performing the calculation of a Fock matrix based on the direct inversion of the iterative subspace (DIIS) algorithm.

class FockDiisModifier : public Scine::Utils::SCFModifier
#include <FockDiisModifier.h>

SCFModifier for the DIIS algorithm. This class does not perform directly the DIIS algorithm, it is like an interface between the SCF method and the class for the actual DIIS algorithm.

template<class IntegerType>
class UniqueRandomNumbersGenerator
#include <UniqueRandomNumbersGenerator.h>

Class to sample without replacement N numbers from a range bounded by min_ and max_.

The user gives a minimum number and a maximum number as the boundaries of a range of possible numbers. Upon calling the UniqueRandomNumbersGenerator::generate(unsigned N) function, N values are sampled without replacement from the possible values.

class AtomicGTOs
#include <atomicGTOs.h>

Class containing the GTO expansions for one atom.

class AtomsOrbitalsIndexes
#include <AtomsOrbitalsIndexes.h>

Structure containing the information about AO indexes and their corresponding atom indexes. TODO: separate state from construction (addAtom). In principle nextAtom_ and nextAO_ shouldn’t be members of this class.

class DensityMatrix
#include <DensityMatrix.h>

Class defining a density matrix for use in electronic structure calculation methods. Is adequate for both restricted and unrestricted formulations. Allows for an easy transfer of density matrices between different instances of a quantum chemical method. There is no overhead if only the restricted formulation is needed.

Accessors to the matrix elements.

Using template functions to allow perfect forwarding to the Eigen functions.

Arithmetic operations

In the case of functions for pairs of density matrices, the precondition is that they both have the same size and are either RHF-RHF or UHF-UHF.

class DensityMatrixIO
#include <DensityMatrixIO.h>

Class to write density matrices to disk or read them from disk, in binary format. TODO: save space and only write and read triangular matrix (since symmetric)? TODO: Make this class employ EigenMatrixIO? NB: then the memory layout would change

class DipoleMatrix

This class stores the dipole integrals <mu|r|nu> in their three components x, y, z.

class GTF
#include <GTF.h>

Class that contains one of the gaussian functions of a ,f.i., STO-nG expansion

class GTOExpansion
#include <GTOExpansion.h>

The class GTOExpansion is the container for the coefficients of a STO-nG expansion.

class MatrixWithDerivatives
#include <MatrixWithDerivatives.h>

Container for a matrix and the derivatives of its elements.

class MolecularOrbitals
#include <MolecularOrbitals.h>

Class for the coefficient (molecular orbital) matrix. Contains all the eigenfunctions (molecular orbitals) produced when solving the (generalized) eigenvalue problem. The contained matrices are therefore quadratic; their dimension is the number of basis functions. TODO: Take SingleParticleEnergies inside this class?

See

OccupiedMolecularOrbitals

class OccupiedMolecularOrbitals
#include <OccupiedMolecularOrbitals.h>

Class for the occupied part of the coefficient (molecular orbital) matrix. Contains only the occupied eigenfunctions (molecular orbitals). The dimension of the matrices is the number of basis functions (rows) and the number of occupied orbitals (columns) TODO: Take SingleParticleEnergies inside this class?

See

MolecularOrbitals

class SingleParticleEnergies
#include <SingleParticleEnergies.h>

Class handling the single-particle energies obtained as eigenvalues of the (generalized) eigenvalue problem. NB: no check is made if the correct (i.e., Restricted / Unrestricted) function variants are called.

class SpinAdaptedMatrix
#include <SpinAdaptedMatrix.h>

Class defining a matrix and which can be used for both spin-restricted and spin-unrestricted formalisms in electronic structure calculation methods. There is only slight overhead if only the restricted formulation is needed.

Accessors to the matrix elements.

Using template functions to allow perfect forwarding to the Eigen functions.

Setters for the matrix elements.

Using template functions to allow perfect forwarding to the Eigen functions.

class STO_nG
#include <STO_nG.h>

The class STO_nG delivers the coefficients of a STO-nG expansion for a given STO. Expansion coefficients reference: 1s-5g: Robert F. Stewart, Small Gaussian Expansions of Slater‐Type Orbitals, J. Chem. Phys. 52, 431 (1970). >= 6s: Prof. Dr. Jürg Hutter.

class OSUtilsModule : public Module
#include <OSUtilsModule.h>

Module provided by OSUtils.

This provides the ‘default’ implementations of common interfaces.

class MolecularTrajectory
#include <MolecularTrajectory.h>

A container class for Scine::Utils::PositionCollection classes.

Mainly an override of std::vector<PositionCollection>, plus Element Type information. Generally used for storing the positions of several consecutive molecular structures.

class ConvergenceCheck
#include <ConvergenceCheck.h>

The base class for all convergence checks.

Similar to the actual Optimizer the convergence check is abstracted into a class that shall allow both usage via Settings object and as a ‘raw’ class using the public members to adjust the convergence criteria.

To this end the base class defines and enforces the presence of the functions that parse the

Settings into these members and vice versa.

In general, it is adviseable to derive classes that also only require the use of a default constructor, after all the convergence checks are mostly defined by a set of thresholds and any data should be required in the actual optimization cycles and not upon creation of the derived object.

Subclassed by Scine::Utils::GradientBasedCheck

class Bofill : public Scine::Utils::Optimizer
#include <Bofill.h>

An implementation of the Bofill optimization algorithm for saddlepoints.

Minimizes along all coordinates except one, along which it maximizes.

Implemented, as described in: Phys. Chem. Chem. Phys., 2002, 4, 11–15 The original paper by Bofill is the following: J. Comput. Chem., 1994, 15, 1

class GradientBasedCheck : public Scine::Utils::ConvergenceCheck
#include <GradientBasedCheck.h>

A convergence check based on parameter, value and gradient information.

Checks the RMS of the gradient and step, the maximum coefficient of the gradient and step and also the change in value.

For convergence the value and a set number of the other four criteria need to converge.

class LBFGS : public Scine::Utils::Optimizer
#include <LBFGS.h>

An implementation of a limited memory BFGS optimization algorithm.

The implementation includes an optional linesearch using both Wolfe conditions.

class SteepestDescent : public Scine::Utils::Optimizer
#include <SteepestDescent.h>

An implementation of a steepest descent optimization algorithm.

The steepest descent algorithm is on of the most simple gradient based optimization algorithms. In each step the parameters are adjusted by subtracting the scaled negative gradient. The gradient is usually scaled by a factor that is smaller than 1.0 .

class EigenvectorFollowing : public Scine::Utils::Optimizer
#include <EigenvectorFollowing.h>

An implementation of an Eigenvector following optimization algorithm.

This algorithm is intended to find the maximum along one single eigenvector and the minimum along all other eigenvectors of a given system/Hessian.

class NewtonRaphson : public Scine::Utils::Optimizer
#include <NewtonRaphson.h>

An implementation of a Newton-Raphson optimization algorithm.

class Optimizer
#include <Optimizer.h>

The base class for all Optimizers.

Given the fact that optimizers can function in many different ways, this base class is somewhat slim and only unifies the usage of auxiliary functions.

In general the optimizers are intended to be used for internal optimizations where a fixed optimizer and a fixed set of options for that optimizer are used, and also those cases where user settings at runtime shall be used.

To this end, the optimizers should have their options given as public member variables which can alternatively be set using the

Optimizer::applySettings function. In order to then allow the addition of the encoded options into a Settings object a Optimizer::addSettingsDescriptors() function has to be generated in any derived class.

Furthermore all optimizers shall track the number of cycles they run and also call the

Optimizer::triggerObservers() function once in each iteration.

The specific function handling the optimization is dependent on the particular type of optimizer, however it should have the following properties:

  • it should return the final cycle count upon completion

  • it should accept the (in this order):

    1. the initial parameters

    2. a template lambda function for the update of needed variables

    3. an object derived from ConvergenceCheck to signal the end of the optimization

Subclassed by Scine::Utils::Bofill, Scine::Utils::EigenvectorFollowing, Scine::Utils::LBFGS, Scine::Utils::NewtonRaphson, Scine::Utils::SteepestDescent

class Settings : public Scine::Utils::UniversalSettings::ValueCollection
#include <Settings.h>

An interface for settings of all sorts.

This interface allows all other interfaces in Scine::Core to have a central place for Settings, while their implementations maintain their specific settings independently.

Implementation Notes:

  • A derived class of the Settings class has to define a construcor that populates the protected _fields member with a set of possible options and their default/descriptions.

Subclassed by Scine::Utils::AFIROptimizerSettings< OptimizerType, ConvergenceCheckType >, Scine::Utils::GeometryOptimizerSettings< OptimizerType, ConvergenceCheckType >, Scine::Utils::IRCOptimizerSettings< OptimizerType, ConvergenceCheckType >, Scine::Utils::OrcaCalculatorSettings

template<class T>
class Abstract
#include <CloneInterface.h>

Strong typing for the abstract level class.

To provide template parameter resolution at compile time, a strong type is provided to differentiate the leaf class from the abstract, middle class.

template<class Derived, class ...Bases>
class CloneInterface : public Bases
#include <CloneInterface.h>

Class serving the purpose of a generic clone interface.

This class provides the clone method for any interface in core needing it. It uses the curiously recurrent template pattern (CRTP) to acquire infos at compile time about both the interface and the derived class. This allows for the one-time-only generation of code for the clone method, avoiding unnecessary code duplication. Basically, a derived class must just inherit from the CloneInterface with template parameters Derived = derived class and Base = interface class to have a functioning clone method. See https://www.fluentcpp.com/2017/09/12/how-to-return-a-smart-pointer-and-use-covariance/ for further details.

Note that classes deriving from CloneInterface<Derived, Base> do implicitly derive from Base, hence they do not need an additional ‘public Base’ statement.

template<class Derived, class ...Bases>
class CloneInterface<Abstract<Derived>, Bases...> : public Bases
#include <CloneInterface.h>

Specialization for abstract classes.

This is necessary, because otherwise the inheritance pattern would not work. It would not be necessary if the inheritance would just have 2 level: interface class and derived class. From the moment that there are also abstract classes between interface and derived class, as the Core::Calculator in the Sparrow module is the most notorious example, the cloneImpl() method must be virtualized in the abstract class.

class ScopedLocale
#include <ScopedLocale.h>

Introduces a scope to change the locale, the original locale will be set back upon destruction.

template<typename T, typename Parameter>
class StrongType
#include <StrongType.h>

From: https://www.fluentcpp.com/2016/12/08/strong-types-for-strong-interfaces/.

Template Parameters
  • T: The underlying type.

  • Parameter: A “phantom type” serves the purpose of specializing the type.

class UniqueIdentifier
#include <UniqueIdentifier.h>

Class for a unique identifier (handle).

It can f.i. be used to identify an instance unequivocally. Uses pimpl idiom to hide boost::uuid dependency.

class Displacement
#include <Typenames.h>

Another name for an Eigen::RowVector3d .

class DisplacementCollection
#include <Typenames.h>

Another name for a row-major Eigen::MatrixX3d .

class ElementTypeCollection
#include <Typenames.h>

Another name for an std::vector<ElementType> .

class Gradient
#include <Typenames.h>

Another name for an Eigen::RowVector3d .

class GradientCollection
#include <Typenames.h>

Another name for a row-major Eigen::MatrixX3d.

class HessianMatrix
#include <Typenames.h>

Another name for a Eigen::MatrixXd .

class Dipole
#include <Typenames.h>

Another name for an Eigen::RowVector3d .

class DipoleGradient
#include <Typenames.h>

Another name for a row-major Eigen::MatrixX3d.

class Position
#include <Typenames.h>

Another name for an Eigen::RowVector3d .

class PositionCollection
#include <Typenames.h>

Another name for a row-major Eigen::MatrixX3d .

namespace AutomaticDifferentiation

Typedefs

using Value1DType = typename Value1DOrder<o>::ValueType

Templated type for a value in 1 dimension.

using Value3DType = typename Value3DOrder<o>::ValueType

Templated type for a value in 3 dimension.

using DerivativeType = typename Derivative3DImpl<o>::ValueType

Templated type for a derivative.

using DerivativeContainerType = typename DerivativeContainer3DImpl<o>::ValueType

Templated type for a derivative container.

Functions

derivOrder getDerivativeOrderEnum(unsigned order)

Get the enum corresponding to the derivative order.

int getDerivativeOrder(derivOrder order)

Get the integer corresponding to the derivative order.

template<derivOrder o>
Value1DType<o> getFromFull(double v, double firstDer, double secondDer)

Create a value with derivatives in 1 dimension and neglect the unnecessary derivatives.

template<derivOrder o>
Value1DType<o> constant1D(double c)

Transform a double to a ValueWithDerivative in one dimension, with derivatives equal to zero.

template<derivOrder o>
Value1DType<o> variableWithUnitDerivative(double v)

Transform v to a ValueWithDerivative in one dimension, with first derivative 1 and second derivative 0.

template<derivOrder o>
Value3DType<o> constant3D(double c)

Transform a double to a ValueWithDerivative in three dimensions, with derivatives equal to zero.

template<derivOrder o>
Value3DType<o> toX(double x)

Get X as a value with derivatives in three dimensions.

The only non-zero derivative is the first derivative in x-direction.

template<derivOrder o>
Value3DType<o> toY(double y)

Get Y as a value with derivatives in three dimensions.

The only non-zero derivative is the first derivative in y-direction.

template<derivOrder o>
Value3DType<o> toZ(double z)

Get Z as a value with derivatives in three dimensions.

The only non-zero derivative is the first derivative in z-direction.

template<derivOrder o>
Value3DType<o> toRSquared(double x, double y, double z)

Get R-squared as a value with derivatives in three dimensions.

template<derivOrder o>
Value3DType<o> get3Dfrom1D(Value1DType<o> v, const Eigen::Vector3d &R)

Get a value with derivatives in 3 dimensions from the value with derivatives in one dimension, given a vector R.

template<derivOrder o>
Value3DType<o> getValueWithOppositeDerivative(const Value3DType<o> &v)

Get the value with inverted derivatives (useful for pairs of derivatives for two atoms).

The sign of the first derivatives changes.

template<derivOrder o>
double getValue1DAsDouble(const Value1DType<o> &v)

Extract the value with derivatives in 1 dimension as a double.

template<derivOrder o>
double getValue3DAsDouble(const Value3DType<o> &v)

Extract the value with derivatives in 3 dimension as a double.

template<>
double constant1D<derivOrder::zero>(double c)
template<>
double variableWithUnitDerivative<derivOrder::zero>(double v)
template<>
double getFromFull<derivOrder::zero>(double v, double, double)
template<>
double constant3D<derivOrder::zero>(double c)
template<>
double toX<derivOrder::zero>(double x)
template<>
double toY<derivOrder::zero>(double y)
template<>
double toZ<derivOrder::zero>(double z)
template<>
double toRSquared<derivOrder::zero>(double x, double y, double z)
template<>
double get3Dfrom1D<derivOrder::zero>(double v, const Eigen::Vector3d&)
template<>
double getValueWithOppositeDerivative<derivOrder::zero>(const double &v)
template<>
double getValue1DAsDouble<derivOrder::zero>(const double &v)
template<>
double getValue3DAsDouble<derivOrder::zero>(const double &v)
template<>
First1D constant1D<derivOrder::one>(double c)
template<>
First1D variableWithUnitDerivative<derivOrder::one>(double v)
template<>
First1D getFromFull<derivOrder::one>(double v, double firstDer, double)
template<>
First3D constant3D<derivOrder::one>(double c)
template<>
First3D toX<derivOrder::one>(double x)
template<>
First3D toY<derivOrder::one>(double y)
template<>
First3D toZ<derivOrder::one>(double z)
template<>
First3D toRSquared<derivOrder::one>(double x, double y, double z)
template<>
First3D get3Dfrom1D<derivOrder::one>(First1D v, const Eigen::Vector3d &R)
template<>
First3D getValueWithOppositeDerivative<derivOrder::one>(const First3D &v)
template<>
double getValue1DAsDouble<derivOrder::one>(const First1D &v)
template<>
double getValue3DAsDouble<derivOrder::one>(const First3D &v)
template<>
Second1D constant1D<derivOrder::two>(double c)
template<>
Second1D variableWithUnitDerivative<derivOrder::two>(double v)
template<>
Second1D getFromFull<derivOrder::two>(double v, double firstDer, double secondDer)
template<>
Second3D constant3D<derivOrder::two>(double c)
template<>
Second3D toX<derivOrder::two>(double x)
template<>
Second3D toY<derivOrder::two>(double y)
template<>
Second3D toZ<derivOrder::two>(double z)
template<>
Second3D toRSquared<derivOrder::two>(double x, double y, double z)
template<>
Second3D get3Dfrom1D<derivOrder::two>(Second1D v, const Eigen::Vector3d &R)
template<>
Second3D getValueWithOppositeDerivative<derivOrder::two>(const Second3D &v)
template<>
double getValue1DAsDouble<derivOrder::two>(const Second1D &v)
template<>
double getValue3DAsDouble<derivOrder::two>(const Second3D &v)
template<typename DerivativeT, typename Crtp>
Crtp operator+(double v, const FirstBase<DerivativeT, Crtp> &rhs)
template<typename DerivativeT, typename Crtp>
Crtp operator-(double v, const FirstBase<DerivativeT, Crtp> &rhs)
template<typename DerivativeT, typename Crtp>
Crtp operator*(double f, const FirstBase<DerivativeT, Crtp> &rhs)
template<typename DerivativeT, typename Crtp>
Crtp operator/(double f, const FirstBase<DerivativeT, Crtp> &rhs)
template<typename DerivativeT, typename Crtp>
Crtp square(const FirstBase<DerivativeT, Crtp> &value)
template<typename DerivativeT, typename Crtp>
Crtp sqrt(const FirstBase<DerivativeT, Crtp> &value)
template<typename DerivativeT, typename Crtp>
Crtp exp(const FirstBase<DerivativeT, Crtp> &value)
template<derivativeType o>
DerivativeType<o> getOppositeDerivative(const DerivativeType<o> &v)

Get opposed derivatives.

i.e. sign of first derivative changes.

template<derivativeType o>
DerivativeType<o> getDerivativeFromValueWithDerivatives(const Value3DType<UnderlyingOrder<o>> &v)

Extract the derivatives from a value type with derivatives.

template<derivativeType o>
void addDerivativeToContainer(DerivativeContainerType<o> &container, int a, int b, const DerivativeType<o> &v)

Add derivatives to a derivative container.

Parameters
  • v: Derivative for the pair of indexes a-b. Considering the direction a->b.

template<>
DerivativeType<derivativeType::first> getOppositeDerivative<derivativeType::first>(const DerivativeType<derivativeType::first> &v)
template<>
DerivativeType<derivativeType::second_atomic> getOppositeDerivative<derivativeType::second_atomic>(const DerivativeType<derivativeType::second_atomic> &v)
template<>
DerivativeType<derivativeType::second_full> getOppositeDerivative<derivativeType::second_full>(const DerivativeType<derivativeType::second_full> &v)
template<>
DerivativeType<derivativeType::first> getDerivativeFromValueWithDerivatives<derivativeType::first>(const Value3DType<derivOrder::one> &v)
template<>
DerivativeType<derivativeType::second_atomic> getDerivativeFromValueWithDerivatives<derivativeType::second_atomic>(const Value3DType<derivOrder::two> &v)
template<>
DerivativeType<derivativeType::second_full> getDerivativeFromValueWithDerivatives<derivativeType::second_full>(const Value3DType<derivOrder::two> &v)
template<>
void addDerivativeToContainer<derivativeType::first>(DerivativeContainerType<derivativeType::first> &container, int a, int b, const DerivativeType<derivativeType::first> &v)
template<>
void addDerivativeToContainer<derivativeType::second_atomic>(DerivativeContainerType<derivativeType::second_atomic> &container, int a, int b, const DerivativeType<derivativeType::second_atomic> &v)
template<>
void addDerivativeToContainer<derivativeType::second_full>(DerivativeContainerType<derivativeType::second_full> &container, int a, int b, const DerivativeType<derivativeType::second_full> &v)
Second1D operator*(double f, const Second1D &der)
Second1D operator+(double f, const Second1D &h)
Second1D operator+(const Second1D &h, double f)
Second1D operator-(double f, const Second1D &h)
Second1D operator-(const Second1D &h, double f)
Second1D operator/(double v, const Second1D &der)
Second1D sqrt(const Second1D &der)
Second1D exp(const Second1D &der)
Second1D cos(const Second1D &der)
Second3D operator*(double f, const Second3D &der)
Second3D sqrt(const Second3D &der)
Second3D exp(const Second3D &der)
Second3D arccos(const Second3D &der)

Variables

constexpr derivOrder UnderlyingOrder = UnderlyingOrderImpl<O>::o

Template variable for the derivative order corresponding to a derivative type.

template<derivOrder O>
struct Value3DOrder
#include <AutomaticDifferentiationTypesHelper.h>

Values in 3 dimensions.

template<derivOrder O>
struct Value1DOrder
#include <AutomaticDifferentiationTypesHelper.h>

Values in 1 dimension.

class First1D : public Scine::Utils::AutomaticDifferentiation::FirstBase<double, First1D>
#include <First1D.h>

Class representing values in one dimension and allowing for the automatic calculation of first derivatives.

class First3D : public Scine::Utils::AutomaticDifferentiation::FirstBase<Eigen::Vector3d, First3D>
#include <First3D.h>

Class representing values in one dimensions and allowing for the automatic calculation of first derivatives.

template<typename DerivativeT, typename Crtp>
class FirstBase
#include <FirstBase.h>

Base class representing values with some derivative type and allowing for the automatic calculation of first derivatives.

This class makes use of CRTP so that the functions will directly have the correct return type for the derived classes.

class FirstND : public Scine::Utils::AutomaticDifferentiation::FirstBase<Eigen::MatrixXd, FirstND>
#include <FirstND.h>

Class representing values in N dimensions and allowing for the automatic calculation of first derivatives in those N dimensions.

class Second1D
#include <Second1D.h>

Class representing values in one dimensions and allowing for the automatic calculation of first and second derivatives.

class Second3D
#include <Second3D.h>

Class representing values in three dimensions and allowing for the automatic calculation of first and second derivatives.

namespace Constants

A namespace for all constant (hardcoded) data.

This namespace does not include fitted parameters for specific methods, only general constant parameters, such as natural constants and atomic data.

Variables

constexpr double elementaryCharge = 1.6021766208e-19
constexpr double avogadroNumber = 6.022140857e23
constexpr double bohrRadius = 0.52917721067e-10
constexpr double pi = 3.14159265358979323846
constexpr double rad_per_degree = pi / 180
constexpr double degree_per_rad = 180 / pi
constexpr double meter_per_bohr = bohrRadius
constexpr double bohr_per_meter = 1 / meter_per_bohr
constexpr double angstrom_per_meter = 1e10
constexpr double meter_per_angstrom = 1 / angstrom_per_meter
constexpr double angstrom_per_bohr = angstrom_per_meter * meter_per_bohr
constexpr double bohr_per_angstrom = 1. / angstrom_per_bohr
constexpr double hartree_per_ev = 3.674932248e-2
constexpr double ev_per_hartree = 1 / hartree_per_ev
constexpr double joule_per_hartree = 4.359744650e-18
constexpr double hartree_per_joule = 1 / joule_per_hartree
constexpr double joule_per_calorie = 4.184
constexpr double calorie_per_joule = 1 / joule_per_calorie
constexpr double kJPerMol_per_hartree = joule_per_hartree / 1000 * avogadroNumber
constexpr double hartree_per_kJPerMol = 1 / kJPerMol_per_hartree
constexpr double kCalPerMol_per_hartree = joule_per_hartree * calorie_per_joule / 1000 * avogadroNumber
constexpr double hartree_per_kCalPerMol = 1 / kCalPerMol_per_hartree
class ElementDataSingleton
#include <ElementData.h>

Provides a mapping of ElementType to ElementData and accessing Elements by a symbol string.

This class is a singleton to avoid multiple instances of the same hardcoded data in memory.

Features:

  • Fast lookup for subscript operator [ElementType]. Throws std::out_of_range.

  • Slow lookup for subscript operator [std::string]. Throws ElementData::DataNotAvailable.

namespace ExternalQC
class Exception : public runtime_error
#include <Exceptions.h>

Base exception.

Subclassed by Scine::Utils::ExternalQC::OutputFileParsingError, Scine::Utils::ExternalQC::UnsuccessfulSystemCommand

class UnsuccessfulSystemCommand : public Scine::Utils::ExternalQC::Exception
#include <Exceptions.h>

Exception thrown when there is an error while executing the command for the external program.

class OutputFileParsingError : public Scine::Utils::ExternalQC::Exception
#include <Exceptions.h>

Exception thrown for errors during output file parsing.

class ExternalProgram
#include <ExternalProgram.h>

This class allows for running external programs through SCINE.

This class is used by the ExternalQC calculators.

class OrcaHessianOutputParser
#include <OrcaHessianOutputParser.h>

This class parses information out of the ORCA hessian output file.

class InputFileCreator

This class creates ORCA input files.

class OrcaMainOutputParser
#include <OrcaMainOutputParser.h>

This class parses information out of the main ORCA output file.

class StateNotAvailableException : public exception
#include <OrcaState.h>

Exception for the case that a state is requested which is not available.

struct OrcaState : public Scine::Utils::State
#include <OrcaState.h>

Definition of a calculation state for ORCA calculations.

The calculation state is defined as a unique identifier. Only a string state is saved here.

namespace Geometry

Functionalities working with an entire geometry (PositionCollection).

Functions

PositionCollection translatePositions(const PositionCollection &positions, const Displacement &translation)

Translates a set of postions by a given displacement.

Return

PositionCollection Returns the translated positions.

Parameters
  • positions: The original positions.

  • translation: The displacement to be added.

void translatePositions(PositionCollection &positions, const Displacement &translation)

Translates a set of postions by a given displacement.

The translatrion happend in-place.

Parameters
  • positions: The initial positions, will be transformed in-place.

  • translation: The displacement to be added.

unsigned int getIndexOfClosestAtom(const PositionCollection &positions, const Position &targetPosition)

Get the index of closest position (atom) to a given position in space.

Return

unsigned int The index of the closest atom to the given target.

Parameters
  • positions: A set of positions to be traversed.

  • targetPosition: The target position.

Eigen::MatrixXd positionVectorToMatrix(const Eigen::VectorXd &v)

Transforms a 3N-dimensional vector {x0, y0, z0, x1, y1, z1, …} to a Nx3 matrix.

Return

Eigen::MatrixXd Returns the final Nx3 matrix.

Parameters
  • v: The positions in vector form.

Eigen::VectorXd positionMatrixToVector(const Eigen::MatrixXd &m)

Transforms a Nx3 matrix to a 3N-dimensional vector {x0, y0, z0, x1, y1, z1, …}.

Return

Eigen::VectorXd Returns the final vector.

Parameters
  • m: The positions in matrix form.

void alignPositions(const PositionCollection &reference, PositionCollection &positions)

Rotate and translate positions so that it is as close as possible to referencePositions.

Parameters
  • reference: The reference positions.

  • positions: The positions to be aligned, will be transformed in place.

std::vector<double> getMasses(const ElementTypeCollection &elements)

Get a vector of all masses (in a.u.).

Return

std::vector<double> Returns the masses listed in a vector.

Parameters
  • elements: A collection of elements.

Position getCenterOfMass(const PositionCollection &positions, const std::vector<double> &masses)

Get the center of mass.

Return

Position Returns the center of mass (COM).

Parameters
  • positions: The positions.

  • masses: The masses (sorted according to the positions).

Position getCenterOfMass(const AtomCollection &structure)

Get the center of mass.

Return

Position Returns the center of mass (COM).

Parameters
  • structure: The structure (positions and masses are relevant).

Position getAveragePosition(const PositionCollection &positions)

Get the average Position.

(The average Position is identical to the center of mass if all masses are identical)

Return

Position Returns the average position.

Parameters
  • positions: The positions.

Eigen::Matrix3d calculateInertiaTensor(const PositionCollection &positions, const std::vector<double> &masses, const Position &centerOfMass)

Calculates the inertia tensor.

Return

Eigen::Matrix3d Returns the inertia tensor.

Parameters
  • positions: The positions.

  • masses: The masses (sorted according to the positions).

  • centerOfMass: The center of mass.

PrincipalMomentsOfInertia calculatePrincipalMoments(const PositionCollection &positions, const std::vector<double> &masses, const Position &centerOfMass)

Calculates the principal moments of inertia.

Return

PrincipalMomentsOfInertia Returns the principal moments of inertia.

Parameters
  • positions: The positions.

  • masses: The masses (sorted according to the positions).

  • centerOfMass: The center of mass.

Eigen::MatrixXd calculateTranslationAndRotationModes(const PositionCollection &positions, const ElementTypeCollection &elements)

Calculated the cartesion modes correspoding to translations and rotations of the entire system.

Return

Eigen::MatrixXd The rotation and translation modes. Translation modes (x,y,z) first-third column, rotation modes 3-final column (depending on the geometry)

Parameters
  • positions: The positions of all atoms.

  • elements: The ElemenetTypes of all atoms.

Eigen::MatrixXd calculateRotTransFreeTransformMatrix(const PositionCollection &positions, const ElementTypeCollection &elements)

Generates the matrix removing rotation and translation modes from the given geometry if applied.

Return

Eigen::MatrixXd The transformation matrix (applied as X^T*H*X to the Hessian).

Parameters
  • positions: The positions of all atoms.

  • elements: The ElemenetTypes of all atoms.

class PrincipalMomentsOfInertia
#include <Geometry.h>

The principal moments of inertia stored.

namespace LcaoUtil

Functions

void getNumberUnrestrictedElectrons(int &nAlpha, int &nBeta, int nElectrons, int spinMultiplicity)

Calculate the numbers of alpha and beta electrons from the total number of electrons and the spin multiplicity.

void solveRestrictedEigenvalueProblem(const SpinAdaptedMatrix &fockMatrix, MolecularOrbitals &coefficientMatrix, SingleParticleEnergies &singleParticleEnergies)
void solveRestrictedGeneralizedEigenvalueProblem(const SpinAdaptedMatrix &fockMatrix, const Eigen::MatrixXd &overlapMatrix, MolecularOrbitals &coefficientMatrix, SingleParticleEnergies &singleParticleEnergies)
void solveUnrestrictedEigenvalueProblem(const SpinAdaptedMatrix &fockMatrix, MolecularOrbitals &coefficientMatrix, SingleParticleEnergies &singleParticleEnergies)
void solveUnrestrictedGeneralizedEigenvalueProblem(const SpinAdaptedMatrix &fockMatrix, const Eigen::MatrixXd &overlapMatrix, MolecularOrbitals &coefficientMatrix, SingleParticleEnergies &singleParticleEnergies)
void calculateRestrictedDensityMatrix(DensityMatrix &densityMatrix, const MolecularOrbitals &coefficientMatrix, int nElectrons)
void calculateUnrestrictedDensityMatrices(DensityMatrix &densityMatrix, const MolecularOrbitals &coefficientMatrix, int nElectrons, int spinMultiplicity)
void calculateRestrictedEnergyWeightedDensityMatrix(Eigen::MatrixXd &energyWeightedDensityMatrix, const MolecularOrbitals &coefficientMatrix, const SingleParticleEnergies &singleParticleEnergies, int nElectrons)
void calculateUnrestrictedEnergyWeightedDensityMatrix(Eigen::MatrixXd &energyWeightedDensityMatrix, const MolecularOrbitals &coefficientMatrix, const SingleParticleEnergies &singleParticleEnergies, int nElectrons, int spinMultiplicity)
void calculateBondOrderMatrix(Utils::BondOrderCollection &bondOrderMatrix, const DensityMatrix &densityMatrix, const Eigen::MatrixXd &overlapMatrix, const AtomsOrbitalsIndexes &aoIndexes)

Computes the lower-triangular bond-order matrix for an non-orthogonal basis.

void calculateOrthonormalBondOrderMatrix(Utils::BondOrderCollection &bondOrderMatrix, const DensityMatrix &densityMatrix, const AtomsOrbitalsIndexes &aoIndexes)

Computes the lower-triangular bond-order matrix for an orthogonal basis.

void calculateOrthonormalAtomicCharges(std::vector<double> &mullikenCharges, const std::vector<double> &coreCharges, const DensityMatrix &densityMatrix, const AtomsOrbitalsIndexes &aoIndexes)
void calculateMullikenAtomicCharges(std::vector<double> &mullikenCharges, const std::vector<double> &coreCharges, const DensityMatrix &densityMatrix, const Eigen::MatrixXd &overlapMatrix, const AtomsOrbitalsIndexes &aoIndexes)
class DensityMatrixBuilder
#include <DensityMatrixBuilder.h>

Class to generate density matrices from coefficient matrices. Different possibilities to do so:

  • Specify the number of electrons

  • Specify which orbitals to consider

  • Specify the number of electrons and some number of swaps TODO: Use the class MolecularOrbitalsManipulation instead of reimplementing the mixes and swaps

class DensityMatrixGenerator
#include <DensityMatrixGenerator.h>

Class to convert an electronic occupation and the corresponding molecular orbitals to a density matrix.

See

DensityMatrixBuilder

class ElectronicOccupation
#include <ElectronicOccupation.h>

Class to hold information about which molecular orbitals are occupied. TODO: Implement fractional occupation ?

class ElectronicOccupationGenerator
#include <ElectronicOccupationGenerator.h>

This interface generates an ElectronicOccupation instance, from which the density matrix can be generated (in combination with the molecular orbitals).

Subclassed by Scine::Utils::LcaoUtil::AufbauPrincipleOccupationGenerator

class EnergyWeightedDensityMatrixBuilder
#include <EnergyWeightedDensityMatrixBuilder.h>

Class to generate energy-weighted density matrices for given occupations

class HFWaveFunctionOverlap
#include <HFWaveFunctionOverlap.h>

Class to calculate the overlap between two HartreeFock-like wave functions.

class HomoLumoGapCalculator
#include <HomoLumoGapCalculator.h>

This class calculates a Homo-Lumo gap, given the single-particle energies and the electronic occupation. In the unrestricted case, puts returns energy difference between the orbitals irrespectively of whether they are alpha- or beta-polarized.

class MolecularOrbitalsManipulation
#include <MolecularOrbitalsManipulation.h>

Class for transformations on molecular orbitals.

namespace Methods
class InitializationException : public runtime_error

Subclassed by Scine::Utils::Methods::ParameterFileCannotBeOpenedException, Scine::Utils::Methods::ParameterFileIsInvalidException, Scine::Utils::Methods::ParametersDoNotExistForElementException, Scine::Utils::Methods::ParametersDoNotExistForElementPairException

namespace MultipleScfSolutions
class RandomOrbitalMixer
#include <RandomOrbitalMixer.h>

Class to randomly mix randomly chosen molecular orbitals. In the restricted case, the alphaHomo_ and betaHomo_ are both holding the value for restrictedHomo_.

namespace Regex

Helpers for std::regex.

The functions deliver regular expression building blocks to be used when composing regular expressions. Functions starting with “capturing” add a parenthesis to capture the value.

Functions

std::string lineBegin()

Regex code for the beginning of a line.

Return

std::string Returns Regex code for the beginning of a line.

std::string lineEnd()

Regex code for the end of a line.

Return

std::string Returns Regex code for the end of a line.

std::string floatingPointNumber()

Regex code for a floating point number.

Return

std::string Returns Regex code for a floating point number.

std::string capturingFloatingPointNumber()

Regex code for a floating point number including capture.

Return

std::string Returns Regex code for a floating point number including capture.

std::string integerNumber()

Regex code for an integer number.

Return

std::string Returns Regex code for an integer number.

std::string capturingIntegerNumber()

Regex code for an integer number including capture.

Return

std::string Returns Regex code for an integer number including capture.

std::string addCaptureParenthesis(const std::string &s)

Adds capture paranteses to any string.

Return

std::string Returns the string s with capture parenthesis added.

Parameters
  • s: The Regex string describing the type to be captured.

namespace SettingsNames

Variables

constexpr const char *orcaMethod = "orca_method"
constexpr const char *orcaNumProcs = "orca_nprocs"
constexpr const char *orcaBinaryPath = "orca_binary_path"
constexpr const char *orcaFilenameBase = "orca_filename_base"
constexpr const char *baseWorkingDirectory = "base_working_directory"
constexpr const char *molecularCharge = "molecular_charge"

SettingsNames provides a consistent list of settings names throughout the whole module.

constexpr const char *spinMultiplicity = "spin_multiplicity"
constexpr const char *unrestrictedCalculation = "unrestricted_calculation"
constexpr const char *selfConsistanceCriterion = "self_consistence_criterion"
constexpr const char *maxIterations = "max_scf_iterations"
constexpr const char *mixer = "scf_mixer"
constexpr const char *loggerVerbosity = "log"
constexpr const char *NDDODipoleApproximation = "nddo_dipole"
constexpr const char *parameterFile = "parameter_file"
constexpr const char *parameterRootDirectory = "parameter_root"
struct SCFMixers
#include <SettingsNames.h>

Struct to contain the name of the mixers available.

struct LogLevels
#include <SettingsNames.h>

Struct to contain the name of the level of verbosity for the logger.

namespace UniversalSettings

Functions

ValueCollection createDefaultValueCollection(const DescriptorCollection &descriptors)
class BoolDescriptor : public Scine::Utils::UniversalSettings::SettingDescriptor
#include <BoolDescriptor.h>

Setting descriptor-derived boolean value.

SettingDescriptor for a bool value. This is a class that contains a string that explains what the setting is for and also a default value.

class CollectionListDescriptor : public Scine::Utils::UniversalSettings::SettingDescriptor
#include <CollectionListDescriptor.h>

SettingDescriptor for multiple SettingsDescriptors.

class DescriptorCollection : public Scine::Utils::UniversalSettings::SettingDescriptor
#include <DescriptorCollection.h>

Setting that wraps a vector<pair<str, GenericDescriptor>>.

A setting (string descriptor, generic value pair) whose value is a list of string, GenericDescriptor pairs, which altogether define a nested list of descriptors.

Information

Checks whether a particular descriptor field exists in the configuration

class DirectoryDescriptor : public Scine::Utils::UniversalSettings::SettingDescriptor
#include <DirectoryDescriptor.h>

SettingDescriptor for a directory path.

class DoubleDescriptor : public Scine::Utils::UniversalSettings::SettingDescriptor
#include <DoubleDescriptor.h>

SettingDescriptor for a floating-point value.

Can set boundaries on valid values and a default value for this kind of setting.

class Exception : public runtime_error
#include <Exceptions.h>

Base class for UniversalSettings exceptions.

Subclassed by Scine::Utils::UniversalSettings::AlreadyExistingDescriptorException, Scine::Utils::UniversalSettings::AlreadyExistingValueException, Scine::Utils::UniversalSettings::EmptyOptionListException, Scine::Utils::UniversalSettings::InexistingDescriptorException, Scine::Utils::UniversalSettings::InexistingValueException, Scine::Utils::UniversalSettings::InvalidDescriptorConversionException, Scine::Utils::UniversalSettings::InvalidValueConversionException, Scine::Utils::UniversalSettings::OptionAlreadyExistsException, Scine::Utils::UniversalSettings::OptionDoesNotExistException, Scine::Utils::UniversalSettings::ValueHasDifferentTypeException

class FileDescriptor : public Scine::Utils::UniversalSettings::SettingDescriptor
#include <FileDescriptor.h>

SettingDescriptor for a file path.

class GenericDescriptor
#include <GenericDescriptor.h>

Wrapper around SettingDescriptor that hides the setting type.

Makes it possible to handle setting descriptors as objects (but is in principle similar to a pointer to the base class SettingDescriptor).

template<typename BaseT>
class GenericInstanceEditor
#include <GenericInstanceEditor.h>

Same as GenericInstanceEditor, defined below, without create function for the case that the constructor needs more parameters, which can then be specified in the derived class.

template<typename Base>
class GenericInstanceEditorWithDefaultConstructor : public Scine::Utils::UniversalSettings::GenericInstanceEditor<Base>
#include <GenericInstanceEditor.h>

Template for a class being able to create and modify instances of some polymorphic type from settings specified in the UniversalSettings syntax.

template<typename BaseEditor, typename T>
class GenericInstanceEditorImpl : public BaseEditor
#include <GenericInstanceEditor.h>

Specification of GenericInstanceEditorWithoutCreateFunction with the type T, to generate automatically some of the virtual functions.

Template Parameters
  • BaseEditor: must be some GenericInstanceEditor

  • T: class to instantiate, must be a * derived class of the template parameter of BaseEditor.

Subclassed by Scine::Utils::UniversalSettings::GenericInstanceEditorWithDefaultConstructorImpl< BaseEditor, T >

template<typename BaseEditor, typename T>
class GenericInstanceEditorWithDefaultConstructorImpl : public Scine::Utils::UniversalSettings::GenericInstanceEditorImpl<BaseEditor, T>
#include <GenericInstanceEditor.h>

Specification of GenericInstanceEditor with the type T, to generate automatically some of the virtual functions.

Template Parameters

class GenericValue
#include <GenericValue.h>

Class that uniformly stores multiple types of values.

Wrapper around some parameter value that hides the actual value. Makes it possible to handle different types as one object.

class InformationOutput
#include <InformationOutput.h>

This class prints information about setting descriptors in a human-readable format.

class IntDescriptor : public Scine::Utils::UniversalSettings::SettingDescriptor
#include <IntDescriptor.h>

SettingDescriptor for an integer value.

class IntListDescriptor : public Scine::Utils::UniversalSettings::SettingDescriptor
#include <IntListDescriptor.h>

SettingDescriptor for a list of integer values.

class OptionListDescriptor : public Scine::Utils::UniversalSettings::SettingDescriptor
#include <OptionListDescriptor.h>

SettingDescriptor for a list of options, of which one must be chosen.

class ParametrizedOptionListDescriptor : public Scine::Utils::UniversalSettings::SettingDescriptor
#include <ParametrizedOptionListDescriptor.h>

SettingDescriptor for a list of options, of which one must be chosen, with corresponding settings (that depend on the exact option).

This is for example useful if one of several algorithms must be chosen for some task and each algorithm has specific settings.

struct ParametrizedOptionValue
#include <ParametrizedOptionValue.h>

Value struct for a ParametrizedOptionList NB: settingsKey is relevant for the generation of input files, as the option and the settings would be saved in different nodes in YAML, for instance. The key for the selectedOption is not registered, it is taken to be the one saved in SettingValueCollection.

class SettingDescriptor
#include <SettingDescriptor.h>

Class to restrict a particular setting’s type and values.

Contains a string explaining what the setting is for, possibly a default value and a function that checks if the contained value is valid for the setting described.

This is an abstract base class for type-specific descriptor implementations. Other classes derive from it not because of the need for polymorphism, but to collect common functionality (i.e. property description).

Subclassed by Scine::Utils::UniversalSettings::BoolDescriptor, Scine::Utils::UniversalSettings::CollectionListDescriptor, Scine::Utils::UniversalSettings::DescriptorCollection, Scine::Utils::UniversalSettings::DirectoryDescriptor, Scine::Utils::UniversalSettings::DoubleDescriptor, Scine::Utils::UniversalSettings::FileDescriptor, Scine::Utils::UniversalSettings::IntDescriptor, Scine::Utils::UniversalSettings::IntListDescriptor, Scine::Utils::UniversalSettings::OptionListDescriptor, Scine::Utils::UniversalSettings::ParametrizedOptionListDescriptor, Scine::Utils::UniversalSettings::StringDescriptor, Scine::Utils::UniversalSettings::StringListDescriptor

class SettingPopulator
#include <SettingPopulator.h>

This class populates the common settings of many calculators.

These settings include molecular charge, spin multiplicity, restricted/unrestricted formalism, and SCF options. It populates a Utils::UniversalSettings::DescriptorCollection with default values.

class StringDescriptor : public Scine::Utils::UniversalSettings::SettingDescriptor
#include <StringDescriptor.h>

SettingDescriptor for a string value.

class StringListDescriptor : public Scine::Utils::UniversalSettings::SettingDescriptor
#include <StringListDescriptor.h>

SettingDescriptor for a list of strings.

class ValueCollection
#include <ValueCollection.h>

Wrapper around vector<pair<string, GenericValue>>

Class holding values corresponding to some SettingDescriptorCollection. The functions addXXX() throw a AlreadyExistingValueException if the given name already exists. The functions getXXX() throw a InexistingValueException if the given value does not exist, and a InvalidValueConversionException if the incorrect type is requested. The functions modifyXXX() assume the key is already existing and throw a InexistingValueException if not, and a InvalidValueConversionException if the type is not the same.

Subclassed by Scine::Utils::Settings