DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
dftBase.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 dftBase_H_
19#define dftBase_H_
20
21#include <vector>
22#include <tuple>
23#include <deal.II/base/tensor_function.h>
24#include "dftParameters.h"
25
26namespace dftfe
27{
28 /**
29 * @brief abstract base class for dft
30 *
31 * @author Sambit Das
32 */
33 class dftBase
34 {
35 public:
36 /**
37 * @brief This is required to correctly delete the derived class object
38 * through the base class ptr
39 */
40 virtual ~dftBase(){};
41
42 virtual void
43 set() = 0;
44
45 virtual void
46 init() = 0;
47
48 virtual void
49 run() = 0;
50
51 virtual void
52 writeMesh() = 0;
53
54 virtual void
56 const std::vector<dealii::Tensor<1, 3, double>> &globalAtomsDisplacements,
57 const double maxJacobianRatioFactor = 1.25,
58 const bool useSingleAtomSolutionsOverride = false) = 0;
59
60 /**
61 *@brief Deforms the domain by the given deformation gradient and reinitializes the
62 * dftClass datastructures.
63 */
64 virtual void
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;
69
70
71 virtual std::tuple<bool, double>
72 solve(const bool computeForces = true,
73 const bool computeStress = true,
74 const bool isRestartGroundStateCalcFromChk = false) = 0;
75
76 virtual void
78
79 virtual void
81
82 virtual double
83 getInternalEnergy() const = 0;
84
85 virtual double
86 getEntropicEnergy() const = 0;
87
88 virtual double
89 getFreeEnergy() const = 0;
90
91 virtual const distributedCPUVec<double> &
92 getRhoNodalOut() const = 0;
93
94 virtual const distributedCPUVec<double> &
96
97 virtual double
99
100 virtual void
102
103 virtual void
105
106 /**
107 * @brief Gets the current atom Locations in cartesian form
108 * (origin at center of domain) from dftClass
109 */
110 virtual const std::vector<std::vector<double>> &
112
113 /**
114 * @brief Gets the current image atom Locations in cartesian form
115 * (origin at center of domain) from dftClass
116 */
117 virtual const std::vector<std::vector<double>> &
119
120 /**
121 * @brief Gets the current image atom ids from dftClass
122 */
123 virtual const std::vector<int> &
124 getImageAtomIDs() const = 0;
125
126
127 /**
128 * @brief Gets the current atom Locations in fractional form
129 * from dftClass (only applicable for periodic and semi-periodic BCs)
130 */
131 virtual const std::vector<std::vector<double>> &
133
134
135
136 /**
137 * @brief Gets the current cell lattice vectors
138 *
139 * @return std::vector<std::vector<double>> 3 \times 3 matrix,lattice[i][j]
140 * corresponds to jth component of ith lattice vector
141 */
142 virtual const std::vector<std::vector<double>> &
143 getCell() const = 0;
144
145 /**
146 * @brief Gets the current cell volume
147 *
148 */
149 virtual double
150 getCellVolume() const = 0;
151
152 /**
153 * @brief Gets the current atom types from dftClass
154 */
155 virtual const std::set<unsigned int> &
156 getAtomTypes() const = 0;
157
158 /**
159 * @brief Gets the current atomic forces (configurational forces) from dftClass
160 */
161 virtual const std::vector<double> &
162 getForceonAtoms() const = 0;
163
164
165 /**
166 * @brief Gets the current cell stress from dftClass
167 */
168 virtual const dealii::Tensor<2, 3, double> &
169 getCellStress() const = 0;
170
171 /**
172 * @brief Get reference to dftParameters object
173 */
174 virtual dftParameters &
176
177 /**
178 * @brief writes the current domain bounding vectors and atom coordinates to files, which are required for
179 * geometry relaxation restart
180
181 */
182 virtual void
184
185
186 /**
187 * @brief writes the current domain bounding vectors and atom coordinates to files for
188 * structural optimization and dynamics restarts. The coordinates are stored
189 * as: 1. fractional for semi-periodic/periodic 2. Cartesian for
190 * non-periodic.
191 * @param[in] Path The folder path to store the atom coordinates required
192 * during restart.
193 */
194 virtual void
195 writeDomainAndAtomCoordinates(const std::string Path) const = 0;
196
197
198 /**
199 * @brief writes atomistics data for subsequent post-processing. Related to
200 * WRITE STRUCTURE ENERGY FORCES DATA POST PROCESS input parameter.
201 * @param[in] Path The folder path to store the atomistics data.
202 */
203 virtual void
204 writeStructureEnergyForcesDataPostProcess(const std::string Path) const = 0;
205
206 /**
207 * @brief writes quadrature grid information and associated spin-up
208 * and spin-down electron-density for post-processing
209 * @param[in] Path The folder path to store the atomistics data.
210 */
211 virtual void
212 writeGSElectronDensity(const std::string Path) const = 0;
213
214
215 virtual const MPI_Comm &
216 getMPIDomain() const = 0;
217
218 virtual const MPI_Comm &
219 getMPIParent() const = 0;
220
221 virtual const MPI_Comm &
222 getMPIInterPool() const = 0;
223
224 virtual const MPI_Comm &
225 getMPIInterBand() const = 0;
226 };
227
228} // namespace dftfe
229
230#endif
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 run()=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 void init()=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 void set()=0
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