DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
dftfeWrapper.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (c) 2017-2025 The Regents of the University of Michigan and DFT-FE
4// authors.
5//
6// This file is part of the DFT-FE code.
7//
8// The DFT-FE code is free software; you can use it, redistribute
9// it, and/or modify it under the terms of the GNU Lesser General
10// Public License as published by the Free Software Foundation; either
11// version 2.1 of the License, or (at your option) any later version.
12// The full text of the license can be found in the file LICENSE at
13// the top level of the DFT-FE distribution.
14//
15// ---------------------------------------------------------------------
16//
17
18#ifndef dftfeWrapper_H_
19#define dftfeWrapper_H_
20
21#include <mpi.h>
22#include <string>
23#include <vector>
24
25namespace dftfe
26{
27 class dftBase;
28 class dftParameters;
29 /**
30 * @brief wrapper class for dftfe
31 *
32 * @author Sambit Das
33 */
35 {
36 public:
37 /**
38 * @brief must be called only once at start of program from all processors
39 * after calling MPI_Init
40 */
41 static void
42 globalHandlesInitialize(const MPI_Comm &mpi_comm_world);
43
44 /**
45 * @brief must be called only once at end of program from all processors
46 * but before calling MPI_Finalize
47 */
48 static void
50
51
52 /**
53 * @brief empty constructor
54 */
56
57 /**
58 * @brief constructor based on input parameter_file
59 */
60 dftfeWrapper(const std::string parameter_file,
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);
68
69 /**
70 * @brief constructor based on input parameter_file and restart
71 * coordinates and domain vectors file paths
72 */
73 dftfeWrapper(const std::string parameter_file,
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);
84
85
86 /**
87 * @brief constructor based on input list of atomic coordinates,
88 * list of atomic numbers,cell, boundary conditions,
89 * Monkhorst-Pack k-point grid, and other optional parameters.
90 * This constructor currently only sets up GGA PBE pseudopotential
91 * DFT calculations using ONCV pseudopotentials in .upf format
92 * (read from DFTFE_PSP_PATH folder provided as an environment
93 * variable). The pseudpotential directory must contain files in the
94 * format: AtomicSymbol.upf
95 *
96 * @param[in] mpi_comm_parent mpi communicator to be used by the
97 * dftfeWrapper.
98 * @param[in] useDevice toggle use of Device accelerated DFT-FE
99 * @param[in] atomicPositionsCart vector of atomic positions for
100 * each atom (in Bohr units), Origin is at cell corner
101 * @param[in] atomicNumbers vector of atomic numbers
102 * @param[in] cell 3 \times 3 matrix in Bohr units, cell[i] denotes the ith
103 * cell vector. DFT-FE requires the cell vectors to form a
104 * right-handed coordinate system i.e.
105 * dotProduct(crossProduct(cell[0],cell[1]),cell[2])>0
106 * @param[in] pbc vector of bools denoting periodic boundary conditions
107 * along the three cell vectors, false denotes non-periodic and true is
108 * periodic
109 * @param[in] mpgrid vector of Monkhorst-Pack grid points along the
110 * reciprocal lattice vector directions for sampling the Brillouin zone
111 * along periodic directions. Default value is a Gamma point.
112 * @param[in] mpgridShift vector of bools where false denotes no shift and
113 * true denotes shift by half the Monkhost-Pack grid spacing. Default value
114 * is no shift.
115 * @param[in] spinPolarizedDFT toggles spin-polarized DFT calculations.
116 * Default value is false
117 * @param[in] startMagnetization Starting magnetization to be used for
118 * spin-polarized DFT calculations (must be between -0.5 and +0.5).
119 * Corresponding magnetization per simulation domain will be
120 * (2 x START MAGNETIZATION x Number of electrons) in Bohr magneton units.
121 * @param[in] fermiDiracSmearingTemp Fermi-Dirac smearing temperature in
122 * Kelvin. Default value is 500 K.
123 * @param[in] npkpt Number of groups of MPI tasks across which the work load
124 * of the irreducible k-points is parallelised. npkpt must be a divisor of
125 * total number of MPI tasks. Default value of 0 internally sets npkt to an
126 * heuristically determined value.
127 * @param[in] meshSize Finite-element mesh size around the atoms in Bohr
128 * units. The default value of 0.8 is sufficient to achieve chemical
129 * accuracy in energy (0.1 mHa/atom discretization error) and forces (0.1
130 * mHa/Bohr discretization error) for the ONCV pseudo-dojo
131 * pseudopotentials. Note that this function assumes a sixth order
132 * finite-element interpolating polynomial
133 * @param[in] scfMixingParameter mixing paramter for SCF fixed point
134 * iteration. Currently the Anderson mixing strategy is used.
135 * @param[in] verbosity printing verbosity. Default value is -1: no printing
136 * @param[in] setDeviceToMPITaskBindingInternally This option is only valid
137 * for Device runs. If set to true Device to MPI task binding is set inside
138 * the DFT-FE code. Default behaviour is false which assumes the binding has
139 * been externally set.
140 */
141 dftfeWrapper(const MPI_Comm & mpi_comm_parent,
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,
150 false,
151 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);
160
161
163
164 /**
165 * @brief clear and reinitialize based on input parameter_file
166 */
167 void
168 reinit(const std::string parameter_file,
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);
176
177 /**
178 * @brief clear and reinitialize based on input parameter_file and restart
179 * coordinates and domain vectors file paths
180 */
181 void
182 reinit(const std::string parameter_file,
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);
193
194 void
195 reinit(const MPI_Comm & mpi_comm_parent,
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,
204 false,
205 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);
214
215 void
217
218 /**
219 * @brief Legacy function (to be deprecated)
220 */
221 void
223
224 /**
225 * @brief Calls dftBasepointer public function writeMesh(). Here the inital density and mesh are stored in a file. Useful for visluiing meshes without running solve.
226 */
227 void
229
230
231 /**
232 * @brief solve ground-state and return DFT free energy which is sum of internal
233 * energy and negative of electronic entropic energy (in Hartree units)
234 *
235 * @return tuple of ground-state energy, boolean flag on whether scf converged,
236 * and L2 norm of residual electron-density of the last SCF iteration
237 */
238 std::tuple<double, bool, double>
239 computeDFTFreeEnergy(const bool computeIonForces = true,
240 const bool computeCellStress = false);
241
242 void
244
245 /**
246 * @brief Get DFT free energy (in Hartree units). This function can
247 * only be called after calling computeDFTFreeEnergy
248 */
249 double
251
252 /**
253 * @brief Get electronic entropic energy (in Hartree units). This function can
254 * only be called after calling computeDFTFreeEnergy
255 */
256 double
258
259 /**
260 * @brief Get ionic forces: negative of gradient of DFT free energy with
261 * respect to ionic positions (in Hartree/Bohr units). This function can
262 * only be called after calling computeDFTFreeEnergy
263 *
264 * @return vector of forces on each atom
265 */
266 std::vector<std::vector<double>>
268
269 /**
270 * @brief Get cell stress: negative of gradient of DFT free energy
271 * with respect to affine strain components scaled by volume
272 * (Hartree/Bohr^3) units. This function can only
273 * be called after calling computeDFTFreeEnergy
274 *
275 * @return cell stress 3 \times 3 matrix given by
276 * sigma[i][j]=\frac{1}{\Omega}\frac{\partial E}{\partial \epsilon_{ij}}
277 */
278 std::vector<std::vector<double>>
280
281 /**
282 * @brief update atom positions and reinitialize all related data-structures
283 *
284 * @param[in] atomsDisplacements vector of displacements for
285 * each atom (in Bohr units)
286 */
287 void
289 const std::vector<std::vector<double>> atomsDisplacements);
290
291
292 /**
293 *@brief Deforms the cell by applying the given affine deformation gradient and
294 * reinitializes the underlying data-structures.
295 *
296 *@param[in] deformationGradient deformation gradient
297 * matrix given by F[i][j]=\frac{\partial x_i}{\partial X_j}
298 */
299 void
300 deformCell(const std::vector<std::vector<double>> deformationGradient);
301
302 /**
303 * @brief Gets the current atom Positions in cartesian form (in Bohr units)
304 * (origin at corner of cell against which the cell vectors are defined)
305 *
306 * @return array of coords for each atom
307 */
308 std::vector<std::vector<double>>
310
311 /**
312 * @brief Gets the current atom Positions in fractional form
313 * (only applicable for periodic and semi-periodic BCs).
314 * CAUTION: during relaxation and MD fractional coordinates may have negaive
315 * values
316 *
317 * @return array of coords for each atom
318 */
319 std::vector<std::vector<double>>
321
322
323
324 /**
325 * @brief Gets the current cell vectors
326 *
327 * @return 3 \times 3 matrix, cell[i][j] corresponds to jth component of
328 * ith cell vector (in Bohr units)
329 */
330 std::vector<std::vector<double>>
331 getCell() const;
332
333
334 /**
335 * @brief Gets the boundary conditions for each cell vector direction
336 *
337 * @return vector of bools, false denotes non-periodic BC and true denotes periodic BC
338 */
339 std::vector<bool>
340 getPBC() const;
341
342 /**
343 * @brief Gets the atomic numbers vector
344 *
345 * @return vector of atomic numbers
346 */
347 std::vector<int>
349
350
351 /**
352 * @brief Gets the number of valence electrons for each atom
353 *
354 * @return array of number of valence for each atom
355 */
356 std::vector<int>
358
359
360 /**
361 * @brief writes the current domain bounding vectors and atom coordinates to files for
362 * structural optimization and dynamics restarts. The coordinates are stored
363 * as: 1. fractional for semi-periodic/periodic 2. Cartesian for
364 * non-periodic.
365 * @param[in] Path The folder path to store the atom coordinates required
366 * during restart.
367 */
368 void
369 writeDomainAndAtomCoordinates(const std::string Path) const;
370
371
372 dftBase *
374
377
378 private:
379 void
381
382 void
383 initialize(const bool setDeviceToMPITaskBindingInternally,
384 const bool useDevice);
385
391 };
392} // namespace dftfe
393#endif
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 ...
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