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#include <TypeConfig.h>
25
26namespace dftfe
27{
28 class dftBase;
29 class dftParameters;
30 /**
31 * @brief wrapper class for dftfe
32 *
33 * @author Sambit Das
34 */
36 {
37 public:
38 /**
39 * @brief must be called only once at start of program from all processors
40 * after calling MPI_Init
41 */
42 static void
43 globalHandlesInitialize(const MPI_Comm &mpi_comm_world);
44
45 /**
46 * @brief must be called only once at end of program from all processors
47 * but before calling MPI_Finalize
48 */
49 static void
51
52
53 /**
54 * @brief empty constructor
55 */
57
58 /**
59 * @brief constructor based on input parameter_file
60 */
61 dftfeWrapper(const std::string parameter_file,
62 const MPI_Comm &mpi_comm_parent,
63 const bool printParams = false,
64 const bool setDeviceToMPITaskBindingInternally = false,
65 const std::string solverMode = "GS",
66 const std::string restartFilesPath = ".",
67 const dftfe::Int _verbosity = 1,
68 const bool useDevice = false);
69
70 /**
71 * @brief constructor based on input parameter_file and restart
72 * coordinates and domain vectors file paths
73 */
74 dftfeWrapper(const std::string parameter_file,
75 const std::string restartCoordsFile,
76 const std::string restartDomainVectorsFile,
77 const MPI_Comm &mpi_comm_parent,
78 const bool printParams = false,
79 const bool setDeviceToMPITaskBindingInternally = false,
80 const std::string solverMode = "GS",
81 const std::string restartFilesPath = ".",
82 const dftfe::Int _verbosity = 1,
83 const bool useDevice = false,
84 const bool isScfRestart = true);
85
86
87 /**
88 * @brief constructor based on input list of atomic coordinates,
89 * list of atomic numbers,cell, boundary conditions,
90 * Monkhorst-Pack k-point grid, and other optional parameters.
91 * This constructor currently only sets up GGA PBE pseudopotential
92 * DFT calculations using ONCV pseudopotentials in .upf format
93 * (read from DFTFE_PSP_PATH folder provided as an environment
94 * variable). The pseudpotential directory must contain files in the
95 * format: AtomicSymbol.upf
96 *
97 * @param[in] mpi_comm_parent mpi communicator to be used by the
98 * dftfeWrapper.
99 * @param[in] useDevice toggle use of Device accelerated DFT-FE
100 * @param[in] atomicPositionsCart vector of atomic positions for
101 * each atom (in Bohr units), Origin is at cell corner
102 * @param[in] atomicNumbers vector of atomic numbers
103 * @param[in] cell 3 \times 3 matrix in Bohr units, cell[i] denotes the ith
104 * cell vector. DFT-FE requires the cell vectors to form a
105 * right-handed coordinate system i.e.
106 * dotProduct(crossProduct(cell[0],cell[1]),cell[2])>0
107 * @param[in] pbc vector of bools denoting periodic boundary conditions
108 * along the three cell vectors, false denotes non-periodic and true is
109 * periodic
110 * @param[in] mpgrid vector of Monkhorst-Pack grid points along the
111 * reciprocal lattice vector directions for sampling the Brillouin zone
112 * along periodic directions. Default value is a Gamma point.
113 * @param[in] mpgridShift vector of bools where false denotes no shift and
114 * true denotes shift by half the Monkhost-Pack grid spacing. Default value
115 * is no shift.
116 * @param[in] spinPolarizedDFT toggles spin-polarized DFT calculations.
117 * Default value is false
118 * @param[in] startMagnetization Starting magnetization to be used for
119 * spin-polarized DFT calculations (must be between -0.5 and +0.5).
120 * Corresponding magnetization per simulation domain will be
121 * (2 x START MAGNETIZATION x Number of electrons) in Bohr magneton units.
122 * @param[in] fermiDiracSmearingTemp Fermi-Dirac smearing temperature in
123 * Kelvin. Default value is 500 K.
124 * @param[in] npkpt Number of groups of MPI tasks across which the work load
125 * of the irreducible k-points is parallelised. npkpt must be a divisor of
126 * total number of MPI tasks. Default value of 0 internally sets npkt to an
127 * heuristically determined value.
128 * @param[in] meshSize Finite-element mesh size around the atoms in Bohr
129 * units. The default value of 0.8 is sufficient to achieve chemical
130 * accuracy in energy (0.1 mHa/atom discretization error) and forces (0.1
131 * mHa/Bohr discretization error) for the ONCV pseudo-dojo
132 * pseudopotentials. Note that this function assumes a sixth order
133 * finite-element interpolating polynomial
134 * @param[in] scfMixingParameter mixing paramter for SCF fixed point
135 * iteration. Currently the Anderson mixing strategy is used.
136 * @param[in] verbosity printing verbosity. Default value is -1: no printing
137 * @param[in] setDeviceToMPITaskBindingInternally This option is only valid
138 * for Device runs. If set to true Device to MPI task binding is set inside
139 * the DFT-FE code. Default behaviour is false which assumes the binding has
140 * been externally set.
141 */
143 const MPI_Comm &mpi_comm_parent,
144 const bool useDevice,
145 const std::vector<std::vector<double>> atomicPositionsCart,
146 const std::vector<dftfe::uInt> atomicNumbers,
147 const std::vector<std::vector<double>> cell,
148 const std::vector<bool> pbc,
149 const std::vector<dftfe::uInt> mpGrid = std::vector<dftfe::uInt>{1, 1, 1},
150 const std::vector<bool> mpGridShift = std::vector<bool>{false,
151 false,
152 false},
153 const bool spinPolarizedDFT = false,
154 const double startMagnetization = 0.0,
155 const double fermiDiracSmearingTemp = 500.0,
156 const dftfe::uInt npkpt = 0,
157 const double meshSize = 0.8,
158 const double scfMixingParameter = 0.2,
159 const dftfe::Int verbosity = -1,
160 const bool setDeviceToMPITaskBindingInternally = false);
161
162
164
165 /**
166 * @brief clear and reinitialize based on input parameter_file
167 */
168 void
169 reinit(const std::string parameter_file,
170 const MPI_Comm &mpi_comm_parent,
171 const bool printParams = false,
172 const bool setDeviceToMPITaskBindingInternally = false,
173 const std::string solverMode = "GS",
174 const std::string restartFilesPath = ".",
175 const dftfe::Int _verbosity = 1,
176 const bool useDevice = false);
177
178 /**
179 * @brief clear and reinitialize based on input parameter_file and restart
180 * coordinates and domain vectors file paths
181 */
182 void
183 reinit(const std::string parameter_file,
184 const std::string restartCoordsFile,
185 const std::string restartDomainVectorsFile,
186 const MPI_Comm &mpi_comm_parent,
187 const bool printParams = false,
188 const bool setDeviceToMPITaskBindingInternally = false,
189 const std::string solverMode = "GS",
190 const std::string restartFilesPath = ".",
191 const dftfe::Int _verbosity = 1,
192 const bool useDevice = false,
193 const bool isScfRestart = true);
194
195 void
196 reinit(const MPI_Comm &mpi_comm_parent,
197 const bool useDevice,
198 const std::vector<std::vector<double>> atomicPositionsCart,
199 const std::vector<dftfe::uInt> atomicNumbers,
200 const std::vector<std::vector<double>> cell,
201 const std::vector<bool> pbc,
202 const std::vector<dftfe::uInt> mpGrid = std::vector<dftfe::uInt>{1,
203 1,
204 1},
205 const std::vector<bool> mpGridShift = std::vector<bool>{false,
206 false,
207 false},
208 const bool spinPolarizedDFT = false,
209 const double startMagnetization = 0.0,
210 const double fermiDiracSmearingTemp = 500.0,
211 const dftfe::uInt npkpt = 0,
212 const double meshSize = 0.8,
213 const double scfMixingParameter = 0.2,
214 const dftfe::Int verbosity = -1,
215 const bool setDeviceToMPITaskBindingInternally = false);
216
217 void
219
220 /**
221 * @brief Legacy function (to be deprecated)
222 */
223 void
225
226 /**
227 * @brief Calls dftBasepointer public function writeMesh(). Here the inital density and mesh are stored in a file. Useful for visluiing meshes without running solve.
228 */
229 void
231
232
233 /**
234 * @brief solve ground-state and return DFT free energy which is sum of internal
235 * energy and negative of electronic entropic energy (in Hartree units)
236 *
237 * @return tuple of ground-state energy, boolean flag on whether scf converged,
238 * and L2 norm of residual electron-density of the last SCF iteration
239 */
240 std::tuple<double, bool, double>
241 computeDFTFreeEnergy(const bool computeIonForces = true,
242 const bool computeCellStress = false);
243
244 void
246
247 /**
248 * @brief Get DFT free energy (in Hartree units). This function can
249 * only be called after calling computeDFTFreeEnergy
250 */
251 double
253
254 /**
255 * @brief Get electronic entropic energy (in Hartree units). This function can
256 * only be called after calling computeDFTFreeEnergy
257 */
258 double
260
261 /**
262 * @brief Get ionic forces: negative of gradient of DFT free energy with
263 * respect to ionic positions (in Hartree/Bohr units). This function can
264 * only be called after calling computeDFTFreeEnergy
265 *
266 * @return vector of forces on each atom
267 */
268 std::vector<std::vector<double>>
270
271 /**
272 * @brief Get cell stress: negative of gradient of DFT free energy
273 * with respect to affine strain components scaled by volume
274 * (Hartree/Bohr^3) units. This function can only
275 * be called after calling computeDFTFreeEnergy
276 *
277 * @return cell stress 3 \times 3 matrix given by
278 * sigma[i][j]=\frac{1}{\Omega}\frac{\partial E}{\partial \epsilon_{ij}}
279 */
280 std::vector<std::vector<double>>
282
283 /**
284 * @brief update atom positions and reinitialize all related data-structures
285 *
286 * @param[in] atomsDisplacements vector of displacements for
287 * each atom (in Bohr units)
288 */
289 void
291 const std::vector<std::vector<double>> atomsDisplacements);
292
293
294 /**
295 *@brief Deforms the cell by applying the given affine deformation gradient and
296 * reinitializes the underlying data-structures.
297 *
298 *@param[in] deformationGradient deformation gradient
299 * matrix given by F[i][j]=\frac{\partial x_i}{\partial X_j}
300 */
301 void
302 deformCell(const std::vector<std::vector<double>> deformationGradient);
303
304 /**
305 * @brief Gets the current atom Positions in cartesian form (in Bohr units)
306 * (origin at corner of cell against which the cell vectors are defined)
307 *
308 * @return array of coords for each atom
309 */
310 std::vector<std::vector<double>>
312
313 /**
314 * @brief Gets the current atom Positions in fractional form
315 * (only applicable for periodic and semi-periodic BCs).
316 * CAUTION: during relaxation and MD fractional coordinates may have negaive
317 * values
318 *
319 * @return array of coords for each atom
320 */
321 std::vector<std::vector<double>>
323
324
325
326 /**
327 * @brief Gets the current cell vectors
328 *
329 * @return 3 \times 3 matrix, cell[i][j] corresponds to jth component of
330 * ith cell vector (in Bohr units)
331 */
332 std::vector<std::vector<double>>
333 getCell() const;
334
335
336 /**
337 * @brief Gets the boundary conditions for each cell vector direction
338 *
339 * @return vector of bools, false denotes non-periodic BC and true denotes periodic BC
340 */
341 std::vector<bool>
342 getPBC() const;
343
344 /**
345 * @brief Gets the atomic numbers vector
346 *
347 * @return vector of atomic numbers
348 */
349 std::vector<dftfe::Int>
351
352
353 /**
354 * @brief Gets the number of valence electrons for each atom
355 *
356 * @return array of number of valence for each atom
357 */
358 std::vector<dftfe::Int>
360
361
362 /**
363 * @brief writes the current domain bounding vectors and atom coordinates to files for
364 * structural optimization and dynamics restarts. The coordinates are stored
365 * as: 1. fractional for semi-periodic/periodic 2. Cartesian for
366 * non-periodic.
367 * @param[in] Path The folder path to store the atom coordinates required
368 * during restart.
369 */
370 void
371 writeDomainAndAtomCoordinates(const std::string Path) const;
372
373
374 dftBase *
376
379
380 private:
381 void
383
384 void
385 initialize(const bool setDeviceToMPITaskBindingInternally,
386 const bool useDevice);
387
393 };
394} // namespace dftfe
395#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:36
void reinit(const MPI_Comm &mpi_comm_parent, const bool useDevice, const std::vector< std::vector< double > > atomicPositionsCart, const std::vector< dftfe::uInt > atomicNumbers, const std::vector< std::vector< double > > cell, const std::vector< bool > pbc, const std::vector< dftfe::uInt > mpGrid=std::vector< dftfe::uInt >{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 dftfe::uInt npkpt=0, const double meshSize=0.8, const double scfMixingParameter=0.2, const dftfe::Int verbosity=-1, const bool setDeviceToMPITaskBindingInternally=false)
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 dftfe::Int _verbosity=1, const bool useDevice=false)
clear and reinitialize based on input parameter_file
std::vector< dftfe::Int > getValenceElectronNumbers() const
Gets the number of valence electrons for each atom.
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 dftfe::Int _verbosity=1, const bool useDevice=false)
constructor based on input parameter_file
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
std::vector< std::vector< double > > getCellStress() const
Get cell stress: negative of gradient of DFT free energy with respect to affine strain components sca...
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.
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 dftfe::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
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< dftfe::Int > getAtomicNumbers() const
Gets the atomic numbers vector.
std::string d_scratchFolderName
Definition dftfeWrapper.h:391
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...
dftfeWrapper(const MPI_Comm &mpi_comm_parent, const bool useDevice, const std::vector< std::vector< double > > atomicPositionsCart, const std::vector< dftfe::uInt > atomicNumbers, const std::vector< std::vector< double > > cell, const std::vector< bool > pbc, const std::vector< dftfe::uInt > mpGrid=std::vector< dftfe::uInt >{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 dftfe::uInt npkpt=0, const double meshSize=0.8, const double scfMixingParameter=0.2, const dftfe::Int verbosity=-1, const bool setDeviceToMPITaskBindingInternally=false)
constructor based on input list of atomic coordinates, list of atomic numbers,cell,...
dftParameters * d_dftfeParamsPtr
Definition dftfeWrapper.h:390
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:389
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...
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 dftfe::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 ...
MPI_Comm d_mpi_comm_parent
Definition dftfeWrapper.h:388
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:392
Definition pseudoPotentialToDftfeConverter.cc:34
std::uint32_t uInt
Definition TypeConfig.h:10
std::int32_t Int
Definition TypeConfig.h:11