23#include <deal.II/base/tensor_function.h>
56 const std::vector<dealii::Tensor<1, 3, double>> &globalAtomsDisplacements,
57 const double maxJacobianRatioFactor = 1.25,
58 const bool useSingleAtomSolutionsOverride =
false) = 0;
65 deformDomain(
const dealii::Tensor<2, 3, double> &deformationGradient,
66 const bool vselfPerturbationUpdateForStress =
false,
67 const bool useSingleAtomSolutionsInitialGuess =
false,
68 const bool print =
true) = 0;
71 virtual std::tuple<bool, double>
72 solve(
const bool computeForces =
true,
74 const bool isRestartGroundStateCalcFromChk =
false) = 0;
110 virtual const std::vector<std::vector<double>> &
117 virtual const std::vector<std::vector<double>> &
123 virtual const std::vector<int> &
131 virtual const std::vector<std::vector<double>> &
142 virtual const std::vector<std::vector<double>> &
155 virtual const std::set<unsigned int> &
161 virtual const std::vector<double> &
168 virtual const dealii::Tensor<2, 3, double> &
215 virtual const MPI_Comm &
218 virtual const MPI_Comm &
221 virtual const MPI_Comm &
224 virtual const MPI_Comm &
abstract base class for dft
Definition dftBase.h:34
virtual double getInternalEnergy() const =0
virtual dftParameters & getParametersObject() const =0
Get reference to dftParameters object.
virtual void updateAtomPositionsAndMoveMesh(const std::vector< dealii::Tensor< 1, 3, double > > &globalAtomsDisplacements, const double maxJacobianRatioFactor=1.25, const bool useSingleAtomSolutionsOverride=false)=0
virtual const std::vector< double > & getForceonAtoms() const =0
Gets the current atomic forces (configurational forces) from dftClass.
virtual const std::vector< std::vector< double > > & getImageAtomLocationsCart() const =0
Gets the current image atom Locations in cartesian form (origin at center of domain) from dftClass.
virtual void deformDomain(const dealii::Tensor< 2, 3, double > &deformationGradient, const bool vselfPerturbationUpdateForStress=false, const bool useSingleAtomSolutionsInitialGuess=false, const bool print=true)=0
Deforms the domain by the given deformation gradient and reinitializes the dftClass datastructures.
virtual ~dftBase()
This is required to correctly delete the derived class object through the base class ptr.
Definition dftBase.h:40
virtual const MPI_Comm & getMPIDomain() const =0
virtual void resetRhoNodalSplitIn(distributedCPUVec< double > &OutDensity)=0
virtual void computeStress()=0
virtual void writeDomainAndAtomCoordinates()=0
writes the current domain bounding vectors and atom coordinates to files, which are required for geom...
virtual const MPI_Comm & getMPIParent() const =0
virtual const std::vector< int > & getImageAtomIDs() const =0
Gets the current image atom ids from dftClass.
virtual void writeGSElectronDensity(const std::string Path) const =0
writes quadrature grid information and associated spin-up and spin-down electron-density for post-pro...
virtual const MPI_Comm & getMPIInterPool() const =0
virtual const distributedCPUVec< double > & getRhoNodalSplitOut() const =0
virtual const distributedCPUVec< double > & getRhoNodalOut() const =0
virtual double getEntropicEnergy() const =0
virtual double getCellVolume() const =0
Gets the current cell volume.
virtual const std::vector< std::vector< double > > & getAtomLocationsFrac() const =0
Gets the current atom Locations in fractional form from dftClass (only applicable for periodic and se...
virtual const MPI_Comm & getMPIInterBand() const =0
virtual void trivialSolveForStress()=0
virtual void writeStructureEnergyForcesDataPostProcess(const std::string Path) const =0
writes atomistics data for subsequent post-processing. Related to WRITE STRUCTURE ENERGY FORCES DATA ...
virtual std::tuple< bool, double > solve(const bool computeForces=true, const bool computeStress=true, const bool isRestartGroundStateCalcFromChk=false)=0
virtual double getTotalChargeforRhoSplit()=0
virtual void writeDomainAndAtomCoordinates(const std::string Path) const =0
writes the current domain bounding vectors and atom coordinates to files for structural optimization ...
virtual double getFreeEnergy() const =0
virtual const std::vector< std::vector< double > > & getAtomLocationsCart() const =0
Gets the current atom Locations in cartesian form (origin at center of domain) from dftClass.
virtual const std::vector< std::vector< double > > & getCell() const =0
Gets the current cell lattice vectors.
virtual const dealii::Tensor< 2, 3, double > & getCellStress() const =0
Gets the current cell stress from dftClass.
virtual void writeMesh()=0
virtual const std::set< unsigned int > & getAtomTypes() const =0
Gets the current atom types from dftClass.
virtual void resetRhoNodalIn(distributedCPUVec< double > &OutDensity)=0
Namespace which declares the input parameters and the functions to parse them from the input paramete...
Definition dftParameters.h:35
Definition pseudoPotentialToDftfeConverter.cc:34
dealii::LinearAlgebra::distributed::Vector< elem_type, dealii::MemorySpace::Host > distributedCPUVec
Definition headers.h:92