18#ifndef dftfeWrapper_H_
19#define dftfeWrapper_H_
61 const MPI_Comm & mpi_comm_parent,
62 const bool printParams =
false,
63 const bool setDeviceToMPITaskBindingInternally =
false,
64 const std::string solverMode =
"GS",
65 const std::string restartFilesPath =
".",
66 const int _verbosity = 1,
67 const bool useDevice =
false);
74 const std::string restartCoordsFile,
75 const std::string restartDomainVectorsFile,
76 const MPI_Comm & mpi_comm_parent,
77 const bool printParams =
false,
78 const bool setDeviceToMPITaskBindingInternally =
false,
79 const std::string solverMode =
"GS",
80 const std::string restartFilesPath =
".",
81 const int _verbosity = 1,
82 const bool useDevice =
false,
83 const bool isScfRestart =
true);
142 const bool useDevice,
143 const std::vector<std::vector<double>> atomicPositionsCart,
144 const std::vector<unsigned int> atomicNumbers,
145 const std::vector<std::vector<double>> cell,
146 const std::vector<bool> pbc,
147 const std::vector<unsigned int> mpGrid =
148 std::vector<unsigned int>{1, 1, 1},
149 const std::vector<bool> mpGridShift = std::vector<bool>{
false,
152 const bool spinPolarizedDFT =
false,
153 const double startMagnetization = 0.0,
154 const double fermiDiracSmearingTemp = 500.0,
155 const unsigned int npkpt = 0,
156 const double meshSize = 0.8,
157 const double scfMixingParameter = 0.2,
158 const int verbosity = -1,
159 const bool setDeviceToMPITaskBindingInternally =
false);
169 const MPI_Comm & mpi_comm_parent,
170 const bool printParams =
false,
171 const bool setDeviceToMPITaskBindingInternally =
false,
172 const std::string solverMode =
"GS",
173 const std::string restartFilesPath =
".",
174 const int _verbosity = 1,
175 const bool useDevice =
false);
183 const std::string restartCoordsFile,
184 const std::string restartDomainVectorsFile,
185 const MPI_Comm & mpi_comm_parent,
186 const bool printParams =
false,
187 const bool setDeviceToMPITaskBindingInternally =
false,
188 const std::string solverMode =
"GS",
189 const std::string restartFilesPath =
".",
190 const int _verbosity = 1,
191 const bool useDevice =
false,
192 const bool isScfRestart =
true);
196 const bool useDevice,
197 const std::vector<std::vector<double>> atomicPositionsCart,
198 const std::vector<unsigned int> atomicNumbers,
199 const std::vector<std::vector<double>> cell,
200 const std::vector<bool> pbc,
201 const std::vector<unsigned int> mpGrid =
202 std::vector<unsigned int>{1, 1, 1},
203 const std::vector<bool> mpGridShift = std::vector<bool>{
false,
206 const bool spinPolarizedDFT =
false,
207 const double startMagnetization = 0.0,
208 const double fermiDiracSmearingTemp = 500.0,
209 const unsigned int npkpt = 0,
210 const double meshSize = 0.8,
211 const double scfMixingParameter = 0.2,
212 const int verbosity = -1,
213 const bool setDeviceToMPITaskBindingInternally =
false);
238 std::tuple<double, bool, double>
240 const bool computeCellStress =
false);
266 std::vector<std::vector<double>>
278 std::vector<std::vector<double>>
289 const std::vector<std::vector<double>> atomsDisplacements);
300 deformCell(
const std::vector<std::vector<double>> deformationGradient);
308 std::vector<std::vector<double>>
319 std::vector<std::vector<double>>
330 std::vector<std::vector<double>>
384 const bool useDevice);
abstract base class for dft
Definition dftBase.h:34
Namespace which declares the input parameters and the functions to parse them from the input paramete...
Definition dftParameters.h:35
std::vector< std::vector< double > > getAtomPositionsFrac() const
Gets the current atom Positions in fractional form (only applicable for periodic and semi-periodic BC...
dftParameters * getDftfeParamsPtr()
void reinit(const std::string parameter_file, const MPI_Comm &mpi_comm_parent, const bool printParams=false, const bool setDeviceToMPITaskBindingInternally=false, const std::string solverMode="GS", const std::string restartFilesPath=".", const int _verbosity=1, const bool useDevice=false)
clear and reinitialize based on input parameter_file
void reinit(const MPI_Comm &mpi_comm_parent, const bool useDevice, const std::vector< std::vector< double > > atomicPositionsCart, const std::vector< unsigned int > atomicNumbers, const std::vector< std::vector< double > > cell, const std::vector< bool > pbc, const std::vector< unsigned int > mpGrid=std::vector< unsigned int >{1, 1, 1}, const std::vector< bool > mpGridShift=std::vector< bool >{false, false, false}, const bool spinPolarizedDFT=false, const double startMagnetization=0.0, const double fermiDiracSmearingTemp=500.0, const unsigned int npkpt=0, const double meshSize=0.8, const double scfMixingParameter=0.2, const int verbosity=-1, const bool setDeviceToMPITaskBindingInternally=false)
void updateAtomPositions(const std::vector< std::vector< double > > atomsDisplacements)
update atom positions and reinitialize all related data-structures
void initialize(const bool setDeviceToMPITaskBindingInternally, const bool useDevice)
void run()
Legacy function (to be deprecated)
static void globalHandlesInitialize(const MPI_Comm &mpi_comm_world)
must be called only once at start of program from all processors after calling MPI_Init
dftfeWrapper(const std::string parameter_file, const std::string restartCoordsFile, const std::string restartDomainVectorsFile, const MPI_Comm &mpi_comm_parent, const bool printParams=false, const bool setDeviceToMPITaskBindingInternally=false, const std::string solverMode="GS", const std::string restartFilesPath=".", const int _verbosity=1, const bool useDevice=false, const bool isScfRestart=true)
constructor based on input parameter_file and restart coordinates and domain vectors file paths
std::vector< std::vector< double > > getCellStress() const
Get cell stress: negative of gradient of DFT free energy with respect to affine strain components sca...
std::vector< int > getAtomicNumbers() const
Gets the atomic numbers vector.
void writeMesh()
Calls dftBasepointer public function writeMesh(). Here the inital density and mesh are stored in a fi...
dftBase * getDftfeBasePtr()
std::vector< bool > getPBC() const
Gets the boundary conditions for each cell vector direction.
static void globalHandlesFinalize()
must be called only once at end of program from all processors but before calling MPI_Finalize
dftfeWrapper()
empty constructor
std::vector< int > getValenceElectronNumbers() const
Gets the number of valence electrons for each atom.
void reinit(const std::string parameter_file, const std::string restartCoordsFile, const std::string restartDomainVectorsFile, const MPI_Comm &mpi_comm_parent, const bool printParams=false, const bool setDeviceToMPITaskBindingInternally=false, const std::string solverMode="GS", const std::string restartFilesPath=".", const int _verbosity=1, const bool useDevice=false, const bool isScfRestart=true)
clear and reinitialize based on input parameter_file and restart coordinates and domain vectors file ...
void createScratchFolder()
dftfeWrapper(const std::string parameter_file, const MPI_Comm &mpi_comm_parent, const bool printParams=false, const bool setDeviceToMPITaskBindingInternally=false, const std::string solverMode="GS", const std::string restartFilesPath=".", const int _verbosity=1, const bool useDevice=false)
constructor based on input parameter_file
dftfeWrapper(const MPI_Comm &mpi_comm_parent, const bool useDevice, const std::vector< std::vector< double > > atomicPositionsCart, const std::vector< unsigned int > atomicNumbers, const std::vector< std::vector< double > > cell, const std::vector< bool > pbc, const std::vector< unsigned int > mpGrid=std::vector< unsigned int >{1, 1, 1}, const std::vector< bool > mpGridShift=std::vector< bool >{false, false, false}, const bool spinPolarizedDFT=false, const double startMagnetization=0.0, const double fermiDiracSmearingTemp=500.0, const unsigned int npkpt=0, const double meshSize=0.8, const double scfMixingParameter=0.2, const int verbosity=-1, const bool setDeviceToMPITaskBindingInternally=false)
constructor based on input list of atomic coordinates, list of atomic numbers,cell,...
std::string d_scratchFolderName
Definition dftfeWrapper.h:389
void deformCell(const std::vector< std::vector< double > > deformationGradient)
Deforms the cell by applying the given affine deformation gradient and reinitializes the underlying d...
double getDFTFreeEnergy() const
Get DFT free energy (in Hartree units). This function can only be called after calling computeDFTFree...
std::vector< std::vector< double > > getForcesAtoms() const
Get ionic forces: negative of gradient of DFT free energy with respect to ionic positions (in Hartree...
dftParameters * d_dftfeParamsPtr
Definition dftfeWrapper.h:388
void writeDomainAndAtomCoordinates(const std::string Path) const
writes the current domain bounding vectors and atom coordinates to files for structural optimization ...
dftBase * d_dftfeBasePtr
Definition dftfeWrapper.h:387
std::vector< std::vector< double > > getCell() const
Gets the current cell vectors.
std::vector< std::vector< double > > getAtomPositionsCart() const
Gets the current atom Positions in cartesian form (in Bohr units) (origin at corner of cell against w...
MPI_Comm d_mpi_comm_parent
Definition dftfeWrapper.h:386
double getElectronicEntropicEnergy() const
Get electronic entropic energy (in Hartree units). This function can only be called after calling com...
std::tuple< double, bool, double > computeDFTFreeEnergy(const bool computeIonForces=true, const bool computeCellStress=false)
solve ground-state and return DFT free energy which is sum of internal energy and negative of electro...
bool d_isDeviceToMPITaskBindingSetInternally
Definition dftfeWrapper.h:390
Definition pseudoPotentialToDftfeConverter.cc:34