File Global2c2eMatrix.h

Copyright

This code is licensed under the 3-clause BSD license.

Copyright ETH Zurich, Laboratory for Physical Chemistry, Reiher Group.

See LICENSE.txt for details.

namespace Scine
namespace Sparrow
namespace nddo
namespace multipole

Functions

template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::first> getDerivative<Utils::derivativeType::first>(orbPair_index_t op1, orbPair_index_t op2) const
template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::second_atomic> getDerivative<Utils::derivativeType::second_atomic>(orbPair_index_t op1, orbPair_index_t op2) const
template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::second_full> getDerivative<Utils::derivativeType::second_full>(orbPair_index_t op1, orbPair_index_t op2) const
class Global2c2eMatrix
#include <Global2c2eMatrix.h>

This class calculates the two-center two-electron integrals in the global coordinate system.

Public Types

using orb_index_t = int
using orbPair_index_t = int

Public Functions

Global2c2eMatrix(int l1, int l2, const ChargeSeparationParameter &D1, const ChargeSeparationParameter &D2, const KlopmanParameter &r1, const KlopmanParameter &r2)

constructor of the Global2c2eMatrix, constaining the 2 centers 2 electrons integrals between two atoms.

Parameters
  • l1: orbital angular momentum quantum number of the first orbital

  • l2: orbital angular momentum quantum number of the second orbital

  • D1:

  • D2:

  • r1:

  • r2:

void setSymmetric(bool sym)

set whether the class refers to integrals within two elements of the same type.

Parameters
  • sym: true if the elements are the same

template<Utils::derivOrder O>
void calculate(const Eigen::Vector3d &Rab)
template<Utils::derivativeType O>
Utils::AutomaticDifferentiation::DerivativeType<O> getDerivative(orb_index_t o1, orb_index_t o2, orb_index_t o3, orb_index_t o4) const
template<Utils::derivativeType O>
Utils::AutomaticDifferentiation::DerivativeType<O> getDerivative(orbPair_index_t op1, orbPair_index_t op2) const
double get(orb_index_t o1, orb_index_t o2, orb_index_t o3, orb_index_t o4) const
double get(orbPair_index_t op1, orbPair_index_t op2) const
orbPair_index_t getPairIndex(orb_index_t o1, orb_index_t o2) const
void output() const

Private Functions

Eigen::RowVector3d getFirstDerivative(orbPair_index_t op1, orbPair_index_t op2) const
Utils::AutomaticDifferentiation::Second3D getSecondDerivative(orbPair_index_t op1, orbPair_index_t op2) const
template<Utils::derivOrder O>
void evaluate()
template<Utils::derivOrder O>
Utils::AutomaticDifferentiation::Value3DType<O> evaluateMatrixElement(orbPair_index_t op1, orbPair_index_t op2)

Private Members

const int d1_
const int d2_
Global2c2eTerms terms_
const ChargeSeparationParameter &dist1
const ChargeSeparationParameter &dist2
const KlopmanParameter &rho1
const KlopmanParameter &rho2
Local2c2eMatrix<Utils::derivOrder::zero> localZero_
Local2c2eMatrix<Utils::derivOrder::one> localOne_
Local2c2eMatrix<Utils::derivOrder::two> localTwo_
OrbitalRotation<Utils::derivOrder::zero> rotationsZero_
OrbitalRotation<Utils::derivOrder::one> rotationsOne_
OrbitalRotation<Utils::derivOrder::two> rotationsTwo_
double R_
double RNormX_
double RNormY_
double RNormZ_
Utils::AutomaticDifferentiation::First3D nullDeriv
Utils::AutomaticDifferentiation::First3D oneDeriv
Eigen::MatrixXd globalMatrix_
Eigen::Matrix<Utils::AutomaticDifferentiation::First3D, Eigen::Dynamic, Eigen::Dynamic> globalMatrixOne_
Eigen::Matrix<Utils::AutomaticDifferentiation::Second3D, Eigen::Dynamic, Eigen::Dynamic> globalMatrixTwo_
bool sameElement_
TwoElectronIntegralIndexes pairIndexes_