DFT-FE 1.3.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
poissonSolverProblemDevice.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/**
19 * @author Gourab Panigrahi
20 *
21 */
22
23#if defined(DFTFE_WITH_DEVICE)
24# ifndef poissonSolverProblemDevice_H_
25# define poissonSolverProblemDevice_H_
26
28# include <constraintMatrixInfo.h>
29# include <constants.h>
30# include <dftUtils.h>
31# include <headers.h>
32# include "FEBasisOperations.h"
33# include "BLASWrapper.h"
34# include "MatrixFreeWrapper.h"
35# include <DeviceAPICalls.h>
36
37namespace dftfe
38{
39 /**
40 * @brief poisson solver problem device class template. template parameter FEOrderElectro
41 * is the finite element polynomial order. The class should not be used with
42 * FLOATING NUCLEAR CHARGES = false or POINT WISE DIRICHLET CONSTRAINT = true
43 *
44 * @author Gourab Panigrahi
45 */
46 template <dftfe::uInt FEOrderElectro>
47 class poissonSolverProblemDevice : public linearSolverProblemDevice
48 {
49 public:
50 /// Constructor
51 poissonSolverProblemDevice(const MPI_Comm &mpi_comm);
52
53 /**
54 * @brief clears all datamembers and reset to original state.
55 *
56 *
57 */
58 void
59 clear();
60
61 /**
62 * @brief reinitialize data structures for total electrostatic potential solve.
63 *
64 * For Hartree electrostatic potential solve give an empty map to the atoms
65 * parameter.
66 *
67 */
68 void
69 reinit(
70 const std::shared_ptr<
71 dftfe::basis::
72 FEBasisOperations<double, double, dftfe::utils::MemorySpace::HOST>>
73 &basisOperationsPtr,
74 distributedCPUVec<double> &x,
75 const dealii::AffineConstraints<double> &constraintMatrix,
76 const dftfe::uInt matrixFreeVectorComponent,
77 const dftfe::uInt matrixFreeQuadratureComponentRhsDensity,
78 const dftfe::uInt matrixFreeQuadratureComponentAX,
79 const std::map<dealii::types::global_dof_index, double> &atoms,
80 const std::map<dealii::CellId, std::vector<double>> &smearedChargeValues,
81 const dftfe::uInt smearedChargeQuadratureId,
82 const dftfe::utils::MemoryStorage<double, dftfe::utils::MemorySpace::HOST>
83 &rhoValues,
84 const std::shared_ptr<
85 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
86 BLASWrapperPtr,
87 const bool isComputeDiagonalA = true,
88 const bool isComputeMeanValueConstraints = false,
89 const bool smearedNuclearCharges = false,
90 const bool isRhoValues = true,
91 const bool isGradSmearedChargeRhs = false,
92 const dftfe::uInt smearedChargeGradientComponentId = 0,
93 const bool storeSmearedChargeRhs = false,
94 const bool reuseSmearedChargeRhs = false,
95 const bool reinitializeFastConstraints = false);
96
97 /**
98 * @brief Compute A matrix multipled by x.
99 *
100 */
101 void
102 computeAX(distributedDeviceVec<double> &Ax,
103 distributedDeviceVec<double> &x);
104
105 /**
106 * @brief Compute right hand side vector for the problem Ax = rhs.
107 *
108 * @param rhs vector for the right hand side values
109 */
110 void
111 computeRhs(distributedCPUVec<double> &rhs);
112
113 /**
114 * @brief get the reference to x field
115 *
116 * @return reference to x field. Assumes x field data structure is already initialized
117 */
118 distributedDeviceVec<double> &
119 getX();
120
121 /**
122 * @brief get the reference to Preconditioner
123 *
124 * @return reference to Preconditioner
125 */
126 distributedDeviceVec<double> &
127 getPreconditioner();
128
129 /**
130 * @brief Copies x from Device to Host
131 *
132 */
133 void
134 copyXfromDeviceToHost();
135
136 /**
137 * @brief distribute x to the constrained nodes.
138 *
139 */
140 void
141 distributeX();
142
143
144 void
145 setX();
146
147
148 private:
149 /**
150 * @brief Sets up the constraints matrix
151 *
152 */
153 void
154 setupConstraints();
155
156 /**
157 * @brief Compute the diagonal of A.
158 *
159 */
160 void
161 computeDiagonalA();
162
163 /**
164 * @brief Compute mean value constraint which is required in case of fully periodic
165 * boundary conditions.
166 *
167 */
168 void
169 computeMeanValueConstraint();
170
171 /**
172 * @brief Mean value constraint distibute
173 *
174 */
175 void
176 meanValueConstraintDistribute(distributedDeviceVec<double> &vec) const;
177
178 /**
179 * @brief Mean value constraint distibute slave to master
180 *
181 */
182 void
183 meanValueConstraintDistributeSlaveToMaster(
184 distributedDeviceVec<double> &vec) const;
185
186 void
187 meanValueConstraintDistributeSlaveToMaster(
188 distributedCPUVec<double> &vec) const;
189
190 /**
191 * @brief Mean value constraint set zero
192 *
193 */
194 void
195 meanValueConstraintSetZero(distributedCPUVec<double> &vec) const;
196
197 /// storage for diagonal of the A matrix
198 distributedCPUVec<double> d_diagonalA;
199 distributedDeviceVec<double> d_diagonalAdevice;
200
201 /// storage for smeared charge rhs in case of total potential solve (doesn't
202 /// change every scf)
203 distributedCPUVec<double> d_rhsSmearedCharge;
204
205 /// pointer to dealii MatrixFree object
206 const dealii::MatrixFree<3, double> *d_matrixFreeDataPtr;
207
208 /// pointer to the x vector being solved for
209 distributedCPUVec<double> *d_xPtr;
210 distributedDeviceVec<double> d_xDevice;
211
212 // number of cells local to each mpi task, number of degrees of freedom
213 // locally owned and total degrees of freedom including ghost
214 dftfe::Int d_nLocalCells, d_xLocalDof, d_xLen;
215
216 // Matrix free wrapper object
217 std::unique_ptr<
218 dftfe::MatrixFreeWrapperClass<double,
221 false>>
222 d_matrixFreeWrapperDevice;
223
224 // constraints
225 dftUtils::constraintMatrixInfo<dftfe::utils::MemorySpace::DEVICE>
226 d_inhomogenousConstraintsTotalPotentialInfo;
227
228 /// pointer to dealii dealii::AffineConstraints<double> object
229 const dealii::AffineConstraints<double> *d_constraintMatrixPtr;
230
231 /// matrix free index required to access the DofHandler and
232 /// dealii::AffineConstraints<double> objects corresponding to the problem
233 dftfe::uInt d_matrixFreeVectorComponent;
234
235 /// matrix free quadrature index
236 dftfe::uInt d_matrixFreeQuadratureComponentRhsDensity;
237
238 /// matrix free quadrature index
239 dftfe::uInt d_matrixFreeQuadratureComponentAX;
240
241 /// pointer to electron density cell quadrature data
242 const dftfe::utils::MemoryStorage<double, dftfe::utils::MemorySpace::HOST>
243 *d_rhoValuesPtr;
244
245 /// pointer to smeared charge cell quadrature data
246 const std::map<dealii::CellId, std::vector<double>>
247 *d_smearedChargeValuesPtr;
248
249 ///
250 dftfe::uInt d_smearedChargeQuadratureId;
251
252 /// pointer to map between global dof index in current processor and the
253 /// atomic charge on that dof
254 const std::map<dealii::types::global_dof_index, double> *d_atomsPtr;
255
256 /// shape function gradient integral storage
257 std::vector<double> d_cellShapeFunctionGradientIntegralFlattened;
258
259 /// storage for mean value constraint vector
260 distributedCPUVec<double> d_meanValueConstraintVec;
261
262 /// storage for mean value constraint device vector
263 distributedDeviceVec<double> d_meanValueConstraintDeviceVec;
264
265 /// boolean flag to query if mean value constraint datastructures are
266 /// precomputed
267 bool d_isMeanValueConstraintComputed;
268
269 ///
270 bool d_isGradSmearedChargeRhs;
271
272 ///
273 bool d_isStoreSmearedChargeRhs;
274
275 ///
276 bool d_isReuseSmearedChargeRhs;
277
278 ///
279 dftfe::uInt d_smearedChargeGradientComponentId;
280
281 /// mean value constraints: mean value constrained node
282 dealii::types::global_dof_index d_meanValueConstraintNodeId;
283
284 /// mean value constrained node local id
285 dealii::types::global_dof_index d_meanValueConstraintNodeIdLocal;
286
287 /// mean value constraints: constrained proc id containing the mean value
288 /// constrained node
289 dftfe::uInt d_meanValueConstraintProcId;
290
291 /// duplicate constraints object with flattened maps for faster access
292 dftUtils::constraintMatrixInfo<dftfe::utils::MemorySpace::HOST>
293 d_constraintsInfo;
294 std::shared_ptr<
295 dftfe::basis::
296 FEBasisOperations<double, double, dftfe::utils::MemorySpace::HOST>>
297 d_basisOperationsPtr;
298 ///
299 std::shared_ptr<
300 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
301 d_BLASWrapperPtr;
302 bool d_isFastConstraintsInitialized;
303 bool d_isHomogenousConstraintsInitialized;
304
305 const MPI_Comm mpi_communicator;
306 const dftfe::uInt n_mpi_processes;
307 const dftfe::uInt this_mpi_process;
308 dealii::ConditionalOStream pcout;
309 };
310
311} // namespace dftfe
312# endif // poissonSolverProblemDevice_H_
313#endif
@ DEVICE
Definition MemorySpaceType.h:36
Definition pseudoPotentialToDftfeConverter.cc:34
std::uint32_t uInt
Definition TypeConfig.h:10
@ Laplace
Definition MatrixFreeDevice.h:38
std::int32_t Int
Definition TypeConfig.h:11