Pseudopotential for a species of ions, and atom positions and other properties for that species. More...
#include <SpeciesInfo.h>
Classes | |
| struct | Constraint |
| Contains the information on the constraints of motion for each ion. More... | |
Public Types | |
| enum | PseudopotentialFormat { Fhi , Uspp , UPF } |
Public Member Functions | |
| void | sync_atpos () |
| update changes in atpos; call whenever atpos is changed (this will update atposManaged and invalidate cached projectors, if any) | |
| void | setup (const Everything &) |
| void | print (FILE *fp) const |
| print ionic positions from current species | |
| void | populationAnalysis (const std::vector< matrix > &RhoAll) const |
| print population analysis given the density matrix in the Lowdin basis | |
| bool | isRelativistic () const |
| whether pseudopotential is relativistic | |
| bool | isUltrasoft () const |
| whether pseudopotential is ultrasoft | |
| PseudopotentialFormat | getPSPFormat () |
| Returns the pseudopotential format. | |
| std::shared_ptr< ColumnBundle > | getV (const ColumnBundle &Cq, const vector3<> *derivDir=0, const int stressDir=-1) const |
| ColumnBundle | getVatom (const ColumnBundle &Cq, int atom) const |
| Single atom version. | |
| int | nProjectors () const |
| total number of projectors for all atoms in this species (number of columns in result of getV) | |
| double | EnlAndGrad (const QuantumNumber &qnum, const diagMatrix &Fq, const matrix &VdagCq, matrix &HVdagCq, int atom=-1) const |
| void | augmentOverlap (const ColumnBundle &Cq, ColumnBundle &OCq, matrix *VdagCq=0) const |
| Accumulate pseudopotential contribution to the overlap in OCq. | |
| void | augmentOverlapDeriv (const ColumnBundle &Cq, ColumnBundle &OCq, ColumnBundle &V, ColumnBundle &dV) const |
| Derivative of augmented overlap operator w.r.t. atomic position. | |
| void | augmentOverlapSecondDeriv (const ColumnBundle &Cq, ColumnBundle &OCq, ColumnBundle &V, ColumnBundle &dV_A, ColumnBundle &dV_B, ColumnBundle &dsqV) const |
| Second derivative of augmented overlap operator w.r.t. atomic position. | |
| void | augmentSpinOverlap (const ColumnBundle &Cq, vector3< matrix > &Sq) const |
| Accumulate pseudopotential contribution to spin overlap of a columnbundle. | |
| void | augmentDensityInit (int atom=-1) |
| Clear internal data and prepare for density augmentation (call before a loop over augmentDensitySpherical per k-point) | |
| void | augmentDensitySpherical (const QuantumNumber &qnum, const diagMatrix &Fq, const matrix &VdagCq, const matrix *VdagdCqL=0, const matrix *VdagdCqR=0, int atom=-1) |
| Accumulate the pseudopotential dependent contribution to the density in the spherical functions nAug (call once per k-point) | |
| void | augmentDensityGrid (ScalarFieldArray &n, int atom=-1, const vector3<> *atposDeriv=0) const |
| Accumulate the spherical augmentation functions nAug to the grid electron density (call only once, after augmentDensitySpherical on all k-points) | |
| void | augmentDensityGridGrad (const ScalarFieldArray &E_n, std::vector< vector3<> > *forces=0, matrix3<> *Eaug_RRT=0) |
| Gradient propagation corresponding to augmentDensityGrid (stores intermediate spherical function results to E_nAug; call only once). Optionally collect forces and stress contributions. | |
| void | augmentDensityGridGradDeriv (const ScalarFieldArray &E_n, int atom, const vector3<> *atposDeriv) |
| void | augmentDensitySphericalGrad (const QuantumNumber &qnum, const matrix &VdagCq, matrix &HVdagCq, int atom=-1) const |
| Gradient propagation corresponding to augmentDensitySpherical (uses intermediate spherical function results from E_nAug; call once per k-point after augmentDensityGridGrad) | |
| size_t | rhoAtom_nMatrices () const |
| void | rhoAtom_initZero (matrix *rhoAtomPtr) const |
| void | rhoAtom_calc (const std::vector< diagMatrix > &F, const std::vector< ColumnBundle > &C, matrix *rhoAtomPtr) const |
| double | rhoAtom_computeU (const matrix *rhoAtomPtr, matrix *U_rhoAtomPtr) const |
| void | rhoAtom_grad (const ColumnBundle &Cq, const matrix *U_rhoAtomPtr, ColumnBundle &HCq) const |
| void | rhoAtom_forces (const std::vector< diagMatrix > &F, const std::vector< ColumnBundle > &C, const matrix *U_rhoAtomPtr, std::vector< vector3<> > &forces, matrix3<> *EU_RRT) const |
| void | rhoAtom_getV (const ColumnBundle &Cq, const matrix *U_rhoAtomPtr, ColumnBundle &Opsi, matrix &M, const vector3<> *derivDir=0, const int stressDir=-1) const |
| void | accumulateAtomicDensity (ScalarFieldTildeArray &nTilde) const |
| Accumulate atomic density from this species. | |
| void | accumulateAtomicPotential (ScalarFieldTilde &dTilde) const |
| Accumulate electrostatic potential of neutral atoms from this species. | |
| void | setAtomicOrbitals (ColumnBundle &Y, bool applyO, int colOffset=0, const vector3<> *derivDir=0, const int stressDir=-1) const |
| Calculate atomic orbitals with/without O preapplied (store in Y with an optional column offset, or calculate derivatives instead) | |
| void | setAtomicOrbitals (ColumnBundle &Y, bool applyO, unsigned n, int l, int colOffset=0, int atomColStride=0, const vector3<> *derivDir=0, const int stressDir=-1) const |
| If non-zero, atomColStride overrides the number of columns between the same orbital of multiple atoms (default = number of orbitals at current n and l) More... | |
| int | nAtomicOrbitals () const |
| return number of atomic orbitals in this species (all atoms) | |
| int | lMaxAtomicOrbitals () const |
| return maximum angular momentum in available atomic orbitals | |
| int | nAtomicOrbitals (int l) const |
| return number of (pseudo-)principal quantum numbers for atomic orbitals of given l | |
| int | atomicOrbitalOffset (unsigned iAtom, unsigned n, int l, int m, int s) const |
| s is 0/1 for up/dn spinors in non-relativistic case, s=0/1 is for j=l+/-0.5 and mj=m+/-0.5 in relativistic case More... | |
| void | updateLocal (ScalarFieldTilde &Vlocps, ScalarFieldTilde &rhoIon, ScalarFieldTilde &nChargeball, ScalarFieldTilde &nCore, ScalarFieldTilde &tauCore, int atom=-1) const |
| Add contributions from this species to Vlocps, rhoIon, nChargeball and nCore/tauCore (if any) | |
| std::vector< vector3<> > | getLocalForces (const ScalarFieldTilde &ccgrad_Vlocps, const ScalarFieldTilde &ccgrad_rhoIon, const ScalarFieldTilde &ccgrad_nChargeball, const ScalarFieldTilde &ccgrad_nCore, const ScalarFieldTilde &ccgrad_tauCore) const |
| Return the local forces (due to Vlocps, rhoIon, nChargeball and nCore/tauCore) | |
| matrix3 | getLocalStress (const ScalarFieldTilde &ccgrad_Vlocps, const ScalarFieldTilde &ccgrad_rhoIon, const ScalarFieldTilde &ccgrad_nChargeball, const ScalarFieldTilde &ccgrad_nCore, const ScalarFieldTilde &ccgrad_tauCore) const |
| Return the local forces (due to Vlocps, rhoIon, nChargeball and nCore/tauCore) | |
| void | accumNonlocalForces (const ColumnBundle &Cq, const matrix &VdagC, const matrix &E_VdagC, const matrix &grad_CdagOCq, std::vector< vector3<> > &forces, matrix3<> *Enl_RRT) const |
Static Public Member Functions | |
| static matrix | getYlmToSpinAngleMatrix (int l, int j2) |
| Spin-angle helper functions: More... | |
| static matrix | getYlmOverlapMatrix (int l, int j2) |
| Get the ((2l+1)*2)x((2l+1)*2) overlap matrix of the spin-spherical harmonics for total angular momentum j (note j2=2*j) | |
Public Attributes | |
| double | Z |
| Valence charge of the species (prefactor to 1/r in long-range part of pseudopotential) | |
| int | atomicNumber |
| Atomic number of the species (0 if unavailable) | |
| string | name |
| Identifier. | |
| string | potfilename |
| pseudopotential filename | |
| string | pulayfilename |
| pulay correction filename | |
| bool | fromWildcard |
| whether this pseudopotential was automatically added using a wildcard (for command printing purposes only) | |
| std::vector< vector3<> > | atpos |
| array of atomic positions of this species | |
| std::vector< vector3<> > | velocities |
| array of atomic velocities (NAN unless running MD) in lattice coordinates | |
| ManagedArray< vector3<> > | atposManaged |
| managed copy of atpos accessed from operator code (for auto cpu/gpu transfers) | |
| double | dE_dnG |
| Derivative of [total energy per atom] w.r.t [nPlanewaves per unit volume] (for Pulay corrections) | |
| double | mass |
| ionic mass (currently unused) | |
| double | coreRadius |
| maximum pseudopotential core radius (used for core overlap checks during ionic/lattice relaxation) | |
| double | ZfullCore |
| number of electrons in full-core correction (atomicNumber - Z - integral(nCore)) | |
| std::vector< Constraint > | constraints |
| List of all constraints on ions of this species. | |
| std::vector< vector3<> > | initialMagneticMoments |
| Initial magnetic moments of each atom (used only for LCAO and symmetries) (x and y magnetizations only used in noncollinear calculations) | |
| double | initialOxidationState |
| Initial oxidation state of this species (only affects LCAO) | |
| matrix | MnlAll |
| block matrix containing Mnl for all l,m | |
| matrix | E_nAug |
| Gradient w.r.t nAug (in the same layout) | |
Pseudopotential for a species of ions, and atom positions and other properties for that species.
| void SpeciesInfo::accumNonlocalForces | ( | const ColumnBundle & | Cq, |
| const matrix & | VdagC, | ||
| const matrix & | E_VdagC, | ||
| const matrix & | grad_CdagOCq, | ||
| std::vector< vector3<> > & | forces, | ||
| matrix3<> * | Enl_RRT | ||
| ) | const |
Propagate gradient with respect to atomic projections (in E_VdagC, along with additional overlap contributions from grad_CdagOC) to forces Additionally accumulate nonlocal stresses if Enl_RRT is non-null
| int SpeciesInfo::atomicOrbitalOffset | ( | unsigned | iAtom, |
| unsigned | n, | ||
| int | l, | ||
| int | m, | ||
| int | s | ||
| ) | const |
s is 0/1 for up/dn spinors in non-relativistic case, s=0/1 is for j=l+/-0.5 and mj=m+/-0.5 in relativistic case
offset of specified atomic orbital in output of current species (when not using the fixed n and l version)
| double SpeciesInfo::EnlAndGrad | ( | const QuantumNumber & | qnum, |
| const diagMatrix & | Fq, | ||
| const matrix & | VdagCq, | ||
| matrix & | HVdagCq, | ||
| int | atom = -1 |
||
| ) | const |
Return non-local energy for this species and quantum number q and optionally accumulate projected electronic gradient in HVdagCq (if non-null)
| std::shared_ptr<ColumnBundle> SpeciesInfo::getV | ( | const ColumnBundle & | Cq, |
| const vector3<> * | derivDir = 0, |
||
| const int | stressDir = -1 |
||
| ) | const |
Get projectors with qnum and basis matching Cq (optionally cached). If derivDir is non-null, return the derivative with respct to Cartesian k direction *derivDir instead (never cached). If stressDir is >=0, then calculate (i,j) component of dVnl/dR . RT where stressDir = 3*i+j
|
static |
Spin-angle helper functions:
Get the ((2l+1)*2)x(j2+1) matrix that transforms the Ylm+spin to the spin-angle functions, where j2=2*j with j = l+/-0.5
| void SpeciesInfo::setAtomicOrbitals | ( | ColumnBundle & | Y, |
| bool | applyO, | ||
| unsigned | n, | ||
| int | l, | ||
| int | colOffset = 0, |
||
| int | atomColStride = 0, |
||
| const vector3<> * | derivDir = 0, |
||
| const int | stressDir = -1 |
||
| ) | const |
If non-zero, atomColStride overrides the number of columns between the same orbital of multiple atoms (default = number of orbitals at current n and l)
Same as above, but for specific n and l.