DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
force.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 force_H_
19#define force_H_
20#include "vselfBinsManager.h"
21#include "dftParameters.h"
22
23#include "constants.h"
24#include "headers.h"
26#include <dftd.h>
27#include <oncvClass.h>
29
30
31namespace dftfe
32{
33 // forward declaration
34 template <dftfe::uInt T1, dftfe::uInt T2, dftfe::utils::MemorySpace memory>
35 class dftClass;
36
37 /**
38 * @brief computes configurational forces in KSDFT
39 *
40 * This class computes and stores the configurational forces corresponding to
41 * geometry optimization. It uses the formulation in the paper by Motamarri
42 * et.al. (https://link.aps.org/doi/10.1103/PhysRevB.97.165132) which provides
43 * an unified approach to atomic forces corresponding to internal atomic
44 * relaxation and cell stress corresponding to cell relaxation.
45 *
46 * @author Sambit Das
47 */
48 template <dftfe::uInt FEOrder,
49 dftfe::uInt FEOrderElectro,
50 dftfe::utils::MemorySpace memorySpace>
52 {
53 friend class dftClass<FEOrder, FEOrderElectro, memorySpace>;
54
55 public:
56 /** @brief Constructor.
57 *
58 * @param _dftPtr pointer to dftClass
59 * @param mpi_comm_parent parent mpi_communicator
60 * @param mpi_comm_domain domain decomposition mpi_communicator
61 */
63 const MPI_Comm &mpi_comm_parent,
64 const MPI_Comm &mpi_comm_domain,
65 const dftParameters &dftParams);
66
67 /** @brief initializes data structures inside forceClass assuming unmoved triangulation.
68 *
69 * initUnmoved is the first step of the initialization/reinitialization of
70 * force class when starting from a new unmoved triangulation. It creates
71 * the dofHandler with linear finite elements and three components
72 * corresponding to the three force components. It also creates the
73 * corresponding constraint matrices which is why an unmoved triangulation
74 * is necessary. Finally this function also initializes the gaussianMovePar
75 * data member.
76 *
77 * @param triangulation reference to unmoved triangulation where the mesh nodes have not
78 * been manually moved.
79 * @param isElectrostaticsMesh boolean parameter specifying whether this triangulatio is to be used for
80 * for the electrostatics part of the configurational force.
81 * @return void.
82 */
83 void
84 initUnmoved(const dealii::Triangulation<3, 3> &triangulation,
85 const dealii::Triangulation<3, 3> &serialTriangulation,
86 const std::vector<std::vector<double>> &domainBoundingVectors,
87 const bool isElectrostaticsMesh);
88
89 /** @brief initializes data structures inside forceClass which depend on the moved mesh.
90 *
91 * initMoved is the second step (first step call initUnmoved) of the
92 * initialization/reinitialization of force class when starting from a new
93 * mesh, and the first step when recomputing forces on a perturbed mesh.
94 * initMoved assumes that the triangulation whose reference was passed to
95 * the forceClass object in the initUnmoved call now has its nodes moved
96 * such that all atomic positions lie on nodes.
97 *
98 * @return void.
99 */
100 void
102 std::vector<const dealii::DoFHandler<3> *> &dofHandlerVectorMatrixFree,
103 std::vector<const dealii::AffineConstraints<double> *>
104 &constraintsVectorMatrixFree,
105 const bool isElectrostaticsMesh);
106
107 /** @brief initializes and precomputes pseudopotential related data structuers required for configurational force
108 * and stress computation.
109 *
110 * This function is only activated for pseudopotential calculations and is
111 * currently called when initializing/reinitializing the dftClass object.
112 * This function initializes and precomputes the pseudopotential
113 * datastructures for local and non-local parts. Separate internal function
114 * calls are made for KB and ONCV projectors.
115 *
116 * @return void.
117 */
118 void
120
121 /** @brief computes the configurational force on all atoms corresponding to a Gaussian generator,
122 * which represents perturbation of the underlying space.
123 *
124 * The Gaussian generator is taken to be exp(-d_gaussianConstant*r^2), r
125 * being the distance from the atom. Currently d_gaussianConstant is
126 * hardcoded to be 4.0. To get the computed atomic forces use getAtomsForces
127 *
128 * @return void.
129 */
130 void
132 const dealii::MatrixFree<3, double> &matrixFreeData,
133 const dispersionCorrection &dispersionCorr,
134 const dftfe::uInt eigenDofHandlerIndex,
135 const dftfe::uInt smearedChargeQuadratureId,
136 const dftfe::uInt lpspQuadratureIdElectro,
137 const dealii::MatrixFree<3, double> &matrixFreeDataElectro,
138 const dftfe::uInt phiTotDofHandlerIndexElectro,
139 const distributedCPUVec<double> &phiTotRhoOutElectro,
140 const std::vector<
142 &rhoOutValues,
143 const std::vector<
145 &gradRhoOutValues,
147 &rhoTotalOutValuesLpsp,
149 &gradRhoTotalOutValuesLpsp,
150 const std::map<dealii::CellId, std::vector<double>> &rhoCoreValues,
151 const std::map<dealii::CellId, std::vector<double>> &gradRhoCoreValues,
152 const std::map<dealii::CellId, std::vector<double>> &hessianRhoCoreValues,
153 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
154 &gradRhoCoreAtoms,
155 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
156 &hessianRhoCoreAtoms,
157 const std::map<dealii::CellId, std::vector<double>> &pseudoVLocElectro,
158 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
159 &pseudoVLocAtomsElectro,
160 const dealii::AffineConstraints<double> &hangingPlusPBCConstraintsElectro,
161 const vselfBinsManager<FEOrder, FEOrderElectro> &vselfBinsManagerElectro);
162
163 /** @brief returns a copy of the configurational force on all global atoms.
164 *
165 * computeAtomsForces must be called prior to this function call.
166 *
167 * @return std::vector<double> flattened array of the configurational force on all atoms,
168 * the three force components on each atom being the leading dimension.
169 * Units- Hartree/Bohr
170 */
171 std::vector<double> &
173
174 /** @brief prints the currently stored configurational forces on atoms and the Gaussian generator constant
175 * used to compute them.
176 *
177 * @return void.
178 */
179 void
181
182 /** @brief Update force generator Gaussian constant.
183 *
184 * @return void.
185 */
186 // void updateGaussianConstant(const double newGaussianConstant);
187
188 /** @brief computes the configurational stress on the domain corresponding to
189 * affine deformation of the periodic cell.
190 *
191 * This function cannot be called for fully non-periodic calculations.
192 *
193 * @return void.
194 */
195 void
197 const dealii::MatrixFree<3, double> &matrixFreeData,
198 const dispersionCorrection &dispersionCorr,
199 const dftfe::uInt eigenDofHandlerIndex,
200 const dftfe::uInt smearedChargeQuadratureId,
201 const dftfe::uInt lpspQuadratureIdElectro,
202 const dealii::MatrixFree<3, double> &matrixFreeDataElectro,
203 const dftfe::uInt phiTotDofHandlerIndexElectro,
204 const distributedCPUVec<double> &phiTotRhoOutElectro,
205 const std::vector<
207 &rhoOutValues,
208 const std::vector<
210 &gradRhoOutValues,
212 &rhoTotalOutValuesLpsp,
214 &gradRhoTotalOutValuesLpsp,
215 const std::map<dealii::CellId, std::vector<double>> &pseudoVLocElectro,
216 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
217 &pseudoVLocAtomsElectro,
218 const std::map<dealii::CellId, std::vector<double>> &rhoCoreValues,
219 const std::map<dealii::CellId, std::vector<double>> &gradRhoCoreValues,
220 const std::map<dealii::CellId, std::vector<double>> &hessianRhoCoreValues,
221 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
222 &gradRhoCoreAtoms,
223 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
224 &hessianRhoCoreAtoms,
225 const dealii::AffineConstraints<double> &hangingPlusPBCConstraintsElectro,
226 const vselfBinsManager<FEOrder, FEOrderElectro> &vselfBinsManagerElectro);
227
228 /** @brief prints the currently stored configurational stress tensor.
229 *
230 * @return void.
231 */
232 void
234
235 /** @brief returns a copy of the current stress tensor value.
236 *
237 * computeStress must be call prior to this function call.
238 *
239 * @return dealii::Tensor<2,3,double> second order stress Tensor in Hartree/Bohr^3
240 */
241 dealii::Tensor<2, 3, double> &
243
244 /** @brief get the value of Gaussian generator parameter (d_gaussianConstant).
245 * Gaussian generator: Gamma(r)= exp(-d_gaussianConstant*r^2).
246 *
247 */
248 // double getGaussianGeneratorParameter() const;
249
250 private:
251 /** @brief Locates and stores the global dof indices of d_dofHandlerForce whose cooridinates match
252 * with the atomic positions.
253 *
254 * @return void.
255 */
256 void
257 locateAtomCoreNodesForce(const dealii::DoFHandler<3> &dofHandlerForce,
258 const dealii::IndexSet &locally_owned_dofsForce,
259 std::map<std::pair<dftfe::uInt, dftfe::uInt>,
260 dftfe::uInt> &atomsForceDofs);
261
262 void
264 const dealii::DoFHandler<3> &dofHandler,
265 const dealii::DoFHandler<3> &dofHandlerForce,
266 const dealii::AffineConstraints<double> &hangingPlusPBCConstraints,
268 std::vector<std::vector<dealii::DoFHandler<3>::active_cell_iterator>>
269 &cellsVselfBallsDofHandler,
270 std::vector<std::vector<dealii::DoFHandler<3>::active_cell_iterator>>
271 &cellsVselfBallsDofHandlerForce,
272 std::vector<std::map<dealii::CellId, dftfe::uInt>>
273 &cellsVselfBallsClosestAtomIdDofHandler,
274 std::map<dftfe::uInt, dftfe::uInt> &AtomIdBinIdLocalDofHandler,
275 std::vector<std::map<dealii::DoFHandler<3>::active_cell_iterator,
276 std::vector<dftfe::uInt>>>
277 &cellFacesVselfBallSurfacesDofHandler,
278 std::vector<std::map<dealii::DoFHandler<3>::active_cell_iterator,
279 std::vector<dftfe::uInt>>>
280 &cellFacesVselfBallSurfacesDofHandlerForce);
281
282 void
284 const dealii::MatrixFree<3, double> &matrixFreeData,
285 const dealii::MatrixFree<3, double> &matrixFreeDataElectro);
286
287 void
289
290 void
292 const dealii::MatrixFree<3, double> &matrixFreeData,
293 const dftfe::uInt eigenDofHandlerIndex,
294 const dftfe::uInt smearedChargeQuadratureId,
295 const dftfe::uInt lpspQuadratureIdElectro,
296 const dealii::MatrixFree<3, double> &matrixFreeDataElectro,
297 const dftfe::uInt phiTotDofHandlerIndexElectro,
298 const distributedCPUVec<double> &phiTotRhoOutElectro,
299 const std::vector<
301 &rhoOutValues,
302 const std::vector<
304 &gradRhoOutValues,
306 &rhoTotalOutValuesLpsp,
308 &gradRhoTotalOutValuesLpsp,
309 const std::map<dealii::CellId, std::vector<double>> &rhoCoreValues,
310 const std::map<dealii::CellId, std::vector<double>> &gradRhoCoreValues,
311 const std::map<dealii::CellId, std::vector<double>> &hessianRhoCoreValues,
312 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
313 &gradRhoCoreAtoms,
314 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
315 &hessianRhoCoreAtoms,
316 const std::map<dealii::CellId, std::vector<double>> &pseudoVLocElectro,
317 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
318 &pseudoVLocAtomsElectro,
319 const vselfBinsManager<FEOrder, FEOrderElectro> &vselfBinsManagerElectro);
320
321 void
323 const dealii::MatrixFree<3, double> &matrixFreeDataElectro,
324 const dftfe::uInt phiTotDofHandlerIndexElectro,
325 const dftfe::uInt smearedChargeQuadratureId,
326 const dftfe::uInt lpspQuadratureIdElectro,
327 const distributedCPUVec<double> &phiTotRhoOutElectro,
329 &rhoTotalOutValues,
331 &rhoTotalOutValuesLpsp,
333 &gradRhoTotalOutValues,
335 &gradRhoTotalOutValuesLpsp,
336 const std::map<dealii::CellId, std::vector<double>> &pseudoVLocElectro,
337 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
338 &pseudoVLocAtomsElectro,
339 const vselfBinsManager<FEOrder, FEOrderElectro> &vselfBinsManagerElectro);
340
341 void
343
344 void
346 const dealii::DoFHandler<3> &dofHandlerElectro,
347 const vselfBinsManager<FEOrder, FEOrderElectro> &vselfBinsManagerElectro,
348 const dealii::MatrixFree<3, double> &matrixFreeDataElectro,
349 const dftfe::uInt smearedChargeQuadratureId);
350
351 void
353
354 void
356 const dealii::MatrixFree<3, double> &matrixFreeData,
357 const dftfe::uInt eigenDofHandlerIndex,
358 const dftfe::uInt smearedChargeQuadratureId,
359 const dftfe::uInt lpspQuadratureIdElectro,
360 const dealii::MatrixFree<3, double> &matrixFreeDataElectro,
361 const dftfe::uInt phiTotDofHandlerIndexElectro,
362 const distributedCPUVec<double> &phiTotRhoOutElectro,
363 const std::vector<
365 &rhoOutValues,
366 const std::vector<
368 &gradRhoOutValues,
370 &rhoTotalOutValuesLpsp,
372 &gradRhoTotalOutValuesLpsp,
373 const std::map<dealii::CellId, std::vector<double>> &rhoCoreValues,
374 const std::map<dealii::CellId, std::vector<double>> &gradRhoCoreValues,
375 const std::map<dealii::CellId, std::vector<double>> &hessianRhoCoreValues,
376 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
377 &gradRhoCoreAtoms,
378 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
379 &hessianRhoCoreAtoms,
380 const std::map<dealii::CellId, std::vector<double>> &pseudoVLocElectro,
381 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
382 &pseudoVLocAtomsElectro,
383 const vselfBinsManager<FEOrder, FEOrderElectro> &vselfBinsManagerElectro);
384
385 void
387 std::map<dftfe::uInt, std::vector<double>>
388 &forceContributionFPSPLocalGammaAtoms,
389 dealii::FEValues<3> &feValues,
390 dealii::FEFaceValues<3> &feFaceValues,
391 dealii::FEEvaluation<3,
392 1,
394 3> &forceEval,
395 const dealii::MatrixFree<3, double> &matrixFreeData,
396 const dftfe::uInt phiTotDofHandlerIndexElectro,
397 const dftfe::uInt cell,
398 const dealii::AlignedVector<dealii::VectorizedArray<double>> &rhoQuads,
399 const dealii::AlignedVector<
400 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>> &gradRhoQuads,
401 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
402 &pseudoVLocAtoms,
404 const std::vector<std::map<dealii::CellId, dftfe::uInt>>
405 &cellsVselfBallsClosestAtomIdDofHandler);
406
407 void
409 std::map<dftfe::uInt, std::vector<double>>
410 &forceContributionSmearedChargesGammaAtoms,
411 dealii::FEEvaluation<3, -1, 1, 3> &forceEval,
412 const dealii::MatrixFree<3, double> &matrixFreeData,
413 const dftfe::uInt cell,
414 const dealii::AlignedVector<
415 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>> &gradPhiTotQuads,
416 const std::vector<dftfe::uInt> &nonTrivialAtomIdsMacroCell,
417 const std::map<dealii::CellId, std::vector<dftfe::Int>>
418 &bQuadAtomIdsAllAtoms,
419 const dealii::AlignedVector<dealii::VectorizedArray<double>>
420 &smearedbQuads);
421
422 void
424 std::map<dftfe::uInt, std::vector<double>>
425 &forceContributionSmearedChargesGammaAtoms,
426 dealii::FEEvaluation<3, -1, 1, 3> &forceEval,
427 const dealii::MatrixFree<3, double> &matrixFreeData,
428 const dftfe::uInt cell,
429 const dealii::AlignedVector<
430 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>>
431 &gradVselfBinQuads,
432 const std::vector<dftfe::uInt> &nonTrivialAtomIdsMacroCell,
433 const std::map<dealii::CellId, std::vector<dftfe::Int>>
434 &bQuadAtomIdsAllAtoms,
435 const dealii::AlignedVector<dealii::VectorizedArray<double>>
436 &smearedbQuads);
437
438 void
440 std::map<dftfe::uInt, std::vector<double>>
441 &forceContributionFNonlinearCoreCorrectionGammaAtoms,
442 dealii::FEEvaluation<
443 3,
444 1,
446 3> &forceEval,
447 const dealii::MatrixFree<3, double> &matrixFreeData,
448 const dftfe::uInt cell,
449 const dealii::AlignedVector<dealii::VectorizedArray<double>> &vxcQuads,
450 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
451 &gradRhoCoreAtoms);
452
453
454 void
456 std::map<dftfe::uInt, std::vector<double>>
457 &forceContributionFNonlinearCoreCorrectionGammaAtoms,
458 dealii::FEEvaluation<
459 3,
460 1,
462 3> &forceEval,
463 const dealii::MatrixFree<3, double> &matrixFreeData,
464 const dftfe::uInt cell,
465 const dealii::AlignedVector<
466 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>> &derExcGradRho,
467 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
468 &hessianRhoCoreAtoms);
469
470 void
472 std::map<dftfe::uInt, std::vector<double>>
473 &forceContributionFNonlinearCoreCorrectionGammaAtoms,
474 dealii::FEEvaluation<
475 3,
476 1,
478 3> &forceEval,
479 const dealii::MatrixFree<3, double> &matrixFreeData,
480 const dftfe::uInt cell,
481 const dealii::AlignedVector<dealii::VectorizedArray<double>>
482 &vxcQuadsSpin0,
483 const dealii::AlignedVector<dealii::VectorizedArray<double>>
484 &vxcQuadsSpin1,
485 const dealii::AlignedVector<
486 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>>
487 &derExcGradRhoSpin0,
488 const dealii::AlignedVector<
489 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>>
490 &derExcGradRhoSpin1,
491 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
492 &gradRhoCoreAtoms,
493 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
494 &hessianRhoCoreAtoms,
495 const bool isXCGGA = false);
496
497 void
499 const std::map<dftfe::uInt, std::vector<double>>
500 &forceContributionFPSPLocalGammaAtoms,
501 const std::map<std::pair<dftfe::uInt, dftfe::uInt>, dftfe::uInt>
502 &atomsForceDofs,
503 const dealii::AffineConstraints<double> &constraintsNoneForce,
504 distributedCPUVec<double> &configForceVectorLinFE);
505
506 void
508 const std::map<dftfe::uInt, std::vector<double>>
509 &forceContributionLocalGammaAtoms,
510 std::vector<double> &accumForcesVector);
511
512
513 void
515 std::map<dftfe::uInt, std::vector<double>>
516 &forceContributionFnlGammaAtoms,
517 const dealii::MatrixFree<3, double> &matrixFreeData,
518 dealii::FEEvaluation<3,
519 1,
522 3> &forceEvalNLP,
523 const std::shared_ptr<
525 nonLocalOp,
526 dftfe::uInt numNonLocalAtomsCurrentProcess,
527 const std::vector<dftfe::Int> &globalChargeIdNonLocalAtoms,
528 const std::vector<dftfe::uInt> &numberPseudoWaveFunctionsPerAtom,
529 const dftfe::uInt cell,
530 const std::map<dealii::CellId, dftfe::uInt> &cellIdToCellNumberMap,
531#ifdef USE_COMPLEX
532 const std::vector<dataTypes::number>
533 &projectorKetTimesPsiTimesVTimesPartOccContractionPsiQuadsFlattened,
534#endif
535 const std::vector<dataTypes::number> &zetaDeltaVQuadsFlattened,
536 const std::vector<dataTypes::number> &
537 projectorKetTimesPsiTimesVTimesPartOccContractionGradPsiQuadsFlattened);
538
539
540 void
542 dealii::AlignedVector<
543 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>> &FVectQuads,
544 const dealii::MatrixFree<3, double> &matrixFreeData,
545 const dftfe::uInt numQuadPoints,
546 const std::shared_ptr<
548 nonLocalOp,
549 const dftfe::uInt numNonLocalAtomsCurrentProcess,
550 const std::vector<dftfe::Int> &globalChargeIdNonLocalAtoms,
551 const std::vector<dftfe::uInt> &numberPseudoWaveFunctionsPerAtom,
552 const dftfe::uInt cell,
553 const std::map<dealii::CellId, dftfe::uInt> &cellIdToCellNumberMap,
554 const std::vector<dataTypes::number> &zetaDeltaVQuadsFlattened,
555 const std::vector<dataTypes::number> &
556 projectorKetTimesPsiTimesVTimesPartOccContractionGradPsiQuadsFlattened);
557
558 void
560 const std::map<dftfe::uInt, std::vector<double>>
561 &forceContributionFnlGammaAtoms);
562
563 void
565 dealii::Tensor<2, 3, double> &stressContribution,
566 const dealii::MatrixFree<3, double> &matrixFreeData,
567 const dftfe::uInt numQuadPoints,
568 const std::vector<double> &jxwQuadsSubCells,
569 const dftfe::uInt cell,
570 const dftfe::uInt numNonLocalAtomsCurrentProcess,
571 const std::shared_ptr<
573 nonLocalOp,
574 const std::vector<dftfe::uInt> &numberPseudoWaveFunctionsPerAtom,
575 const std::map<dealii::CellId, dftfe::uInt> &cellIdToCellNumberMap,
576 const std::vector<dataTypes::number> &zetalmDeltaVlProductDistImageAtoms,
577#ifdef USE_COMPLEX
578 const std::vector<dataTypes::number>
579 &projectorKetTimesPsiTimesVTimesPartOccContractionPsiQuadsFlattened,
580#endif
581 const std::vector<dataTypes::number>
582 &projectorKetTimesPsiTimesVTimesPartOccContractionGradPsiQuadsFlattened,
583 const bool isSpinPolarized);
584
585 void
587 bool allowGaussianOverlapOnAtoms = false);
588
589 void
591
592 void
594 const dealii::DoFHandler<3> &dofHandlerElectro,
595 const vselfBinsManager<FEOrder, FEOrderElectro> &vselfBinsManagerElectro,
596 const dealii::MatrixFree<3, double> &matrixFreeDataElectro,
597 const dftfe::uInt smearedChargeQuadratureId);
598
599 void
601 const dealii::MatrixFree<3, double> &matrixFreeData,
602 const dftfe::uInt eigenDofHandlerIndex,
603 const dftfe::uInt smearedChargeQuadratureId,
604 const dftfe::uInt lpspQuadratureIdElectro,
605 const dealii::MatrixFree<3, double> &matrixFreeDataElectro,
606 const dftfe::uInt phiTotDofHandlerIndexElectro,
607 const distributedCPUVec<double> &phiTotRhoOutElectro,
608 const std::vector<
610 &rhoOutValues,
611 const std::vector<
613 &gradRhoOutValues,
615 &rhoTotalOutValuesLpsp,
617 &gradRhoTotalOutValuesLpsp,
618 const std::map<dealii::CellId, std::vector<double>> &pseudoVLocElectro,
619 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
620 &pseudoVLocAtomsElectro,
621 const std::map<dealii::CellId, std::vector<double>> &rhoCoreValues,
622 const std::map<dealii::CellId, std::vector<double>> &gradRhoCoreValues,
623 const std::map<dealii::CellId, std::vector<double>> &hessianRhoCoreValues,
624 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
625 &gradRhoCoreAtoms,
626 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
627 &hessianRhoCoreAtoms,
628 const vselfBinsManager<FEOrder, FEOrderElectro> &vselfBinsManagerElectro);
629
630 void
632 const dealii::MatrixFree<3, double> &matrixFreeDataElectro,
633 const dftfe::uInt phiTotDofHandlerIndexElectro,
634 const dftfe::uInt smearedChargeQuadratureId,
635 const dftfe::uInt lpspQuadratureIdElectro,
636 const distributedCPUVec<double> &phiTotRhoOutElectro,
638 &rhoTotalOutValues,
640 &rhoTotalOutValuesLpsp,
642 &gradRhoTotalOutValuesElectro,
644 &gradRhoTotalOutValuesElectroLpsp,
645 const std::map<dealii::CellId, std::vector<double>> &pseudoVLocElectro,
646 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
647 &pseudoVLocAtomsElectro,
648 const vselfBinsManager<FEOrder, FEOrderElectro> &vselfBinsManagerElectro);
649
650 void
652 dealii::FEValues<3> &feValues,
653 dealii::FEFaceValues<3> &feFaceValues,
654 dealii::FEEvaluation<3,
655 1,
657 3> &forceEval,
658 const dealii::MatrixFree<3, double> &matrixFreeData,
659 const dftfe::uInt phiTotDofHandlerIndexElectro,
660 const dftfe::uInt cell,
661 const dealii::AlignedVector<dealii::VectorizedArray<double>> &rhoQuads,
662 const dealii::AlignedVector<
663 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>> &gradRhoQuads,
664 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
665 &pseudoVLocAtoms,
667 const std::vector<std::map<dealii::CellId, dftfe::uInt>>
668 &cellsVselfBallsClosestAtomIdDofHandler);
669
670 void
672 dealii::FEEvaluation<
673 3,
674 1,
676 3> &forceEval,
677 const dealii::MatrixFree<3, double> &matrixFreeData,
678 const dftfe::uInt cell,
679 const dealii::AlignedVector<dealii::VectorizedArray<double>> &vxcQuads,
680 const dealii::AlignedVector<
681 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>> &derExcGradRho,
682 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
683 &gradRhoCoreAtoms,
684 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
685 &hessianRhoCoreAtoms);
686
687 void
689 dealii::FEEvaluation<
690 3,
691 1,
693 3> &forceEval,
694 const dealii::MatrixFree<3, double> &matrixFreeData,
695 const dftfe::uInt cell,
696 const dealii::AlignedVector<dealii::VectorizedArray<double>>
697 &vxcQuadsSpin0,
698 const dealii::AlignedVector<dealii::VectorizedArray<double>>
699 &vxcQuadsSpin1,
700 const dealii::AlignedVector<
701 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>>
702 &derExcGradRhoSpin0,
703 const dealii::AlignedVector<
704 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>>
705 &derExcGradRhoSpin1,
706 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
707 &gradRhoCoreAtoms,
708 const std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
709 &hessianRhoCoreAtoms,
710 const bool isXCGGA = false);
711
712 void
714 dealii::FEEvaluation<3, -1, 1, 3> &forceEval,
715 const dealii::MatrixFree<3, double> &matrixFreeData,
716 const dftfe::uInt cell,
717 const dealii::AlignedVector<
718 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>> &gradPhiTotQuads,
719 const std::vector<dftfe::uInt> &nonTrivialAtomImageIdsMacroCell,
720 const std::map<dealii::CellId, std::vector<dftfe::Int>>
721 &bQuadAtomIdsAllAtomsImages,
722 const dealii::AlignedVector<dealii::VectorizedArray<double>>
723 &smearedbQuads);
724
725 void
727 dealii::FEEvaluation<3, -1, 1, 3> &forceEval,
728 const dealii::MatrixFree<3, double> &matrixFreeData,
729 const dftfe::uInt cell,
730 const dealii::AlignedVector<
731 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>> &gradVselfQuads,
732 const std::vector<dftfe::uInt> &nonTrivialAtomImageIdsMacroCell,
733 const std::map<dealii::CellId, std::vector<dftfe::Int>>
734 &bQuadAtomIdsAllAtomsImages,
735 const dealii::AlignedVector<dealii::VectorizedArray<double>>
736 &smearedbQuads);
737
738 void
740
741
742 /// Parallel distributed vector field which stores the configurational force
743 /// for each fem node corresponding to linear shape function generator (see
744 /// equations 52-53 in
745 /// (https://link.aps.org/doi/10.1103/PhysRevB.97.165132)). This vector
746 /// doesn't contain contribution from terms which have sums over k points.
748
749 /// Parallel distributed vector field which stores the configurational force
750 /// for each fem node corresponding to linear shape function generator (see
751 /// equations 52-53 in
752 /// (https://link.aps.org/doi/10.1103/PhysRevB.97.165132)). This vector only
753 /// containts contribution from the electrostatic part.
755
756#ifdef USE_COMPLEX
757 /// Parallel distributed vector field which stores the configurational force
758 /// for each fem node corresponding to linear shape function generator (see
759 /// equations 52-53 in
760 /// (https://link.aps.org/doi/10.1103/PhysRevB.97.165132)). This vector only
761 /// contains contribution from terms which have sums over k points.
762 distributedCPUVec<double> d_configForceVectorLinFEKPoints;
763#endif
764
765
766 std::vector<double> d_forceAtomsFloating;
767
768#ifdef USE_COMPLEX
769 std::vector<double> d_forceAtomsFloatingKPoints;
770#endif
771
772
773
774 /// Gaussian generator constant. Gaussian generator: Gamma(r)=
775 /// exp(-d_gaussianConstant*r^2)
776 /// FIXME: Until the hanging nodes surface integral issue is fixed use a
777 /// value >=4.0
778 // double d_gaussianConstant;
779
780 /// Storage for configurational force on all global atoms.
781 std::vector<double> d_globalAtomsForces;
782
783 /// Storage for configurational stress tensor
784 dealii::Tensor<2, 3, double> d_stress;
785
786
787
788 /* Part of the stress tensor which is summed over k points.
789 * It is a temporary data structure required for stress evaluation
790 * (d_stress) when parallization over k points is on.
791 */
792 dealii::Tensor<2, 3, double> d_stressKPoints;
793
794 /* Dont use true except for debugging forces only without mesh movement, as
795 * gaussian ovelap on atoms for move mesh is by default set to false
796 */
798
799 /// pointer to dft class
801
802 /// Finite element object for configurational force computation. Linear
803 /// finite elements with three force field components are used.
804 dealii::FESystem<3> FEForce;
805
806 /* DofHandler on which we define the configurational force field. Each
807 * geometric fem node has three dofs corresponding the the three force
808 * components. The generator for the configurational force on the fem node
809 * is the linear shape function attached to it. This DofHandler is based on
810 * the same triangulation on which we solve the dft problem.
811 */
812 dealii::DoFHandler<3> d_dofHandlerForce;
813
814 /* DofHandler on which we define the configurational force field from
815 * electrostatic part (without psp). Each geometric fem node has three dofs
816 * corresponding the the three force components. The generator for the
817 * configurational force on the fem node is the linear shape function
818 * attached to it. This DofHandler is based on the same triangulation on
819 * which we solve the dft problem.
820 */
821 dealii::DoFHandler<3> d_dofHandlerForceElectro;
822
823 /// Index of the d_dofHandlerForce in the MatrixFree object stored in
824 /// dftClass. This is required to correctly use FEEvaluation class.
826
827 /// Index of the d_dofHandlerForceElectro in the MatrixFree object stored in
828 /// dftClass. This is required to correctly use FEEvaluation class.
830
831 /// dealii::IndexSet of locally owned dofs of in d_dofHandlerForce the
832 /// current processor
834
835 /// dealii::IndexSet of locally owned dofs of in d_dofHandlerForceElectro
836 /// the current processor
838
839 /// dealii::IndexSet of locally relevant dofs of in d_dofHandlerForce the
840 /// current processor
842
843 /// dealii::IndexSet of locally relevant dofs of in d_dofHandlerForceElectro
844 /// the current processor
846
847 /// Constraint matrix for hanging node and periodic constaints on
848 /// d_dofHandlerForce.
849 dealii::AffineConstraints<double> d_constraintsNoneForce;
850
851 /// Constraint matrix for hanging node and periodic constaints on
852 /// d_dofHandlerForceElectro.
853 dealii::AffineConstraints<double> d_constraintsNoneForceElectro;
854
855 /// Internal data: map < <atomId,force component>, globaldof in
856 /// d_dofHandlerForce>
857 std::map<std::pair<dftfe::uInt, dftfe::uInt>, dftfe::uInt> d_atomsForceDofs;
858
859 /// Internal data: map < <atomId,force component>, globaldof in
860 /// d_dofHandlerForceElectro>
861 std::map<std::pair<dftfe::uInt, dftfe::uInt>, dftfe::uInt>
863
864 /// Internal data: stores cell iterators of all cells in
865 /// dftPtr->d_dofHandler which are part of the vself ball. Outer vector is
866 /// over vself bins.
867 std::vector<std::vector<dealii::DoFHandler<3>::active_cell_iterator>>
869
870 /// Internal data: stores cell iterators of all cells in d_dofHandlerForce
871 /// which are part of the vself ball. Outer vector is over vself bins.
872 std::vector<std::vector<dealii::DoFHandler<3>::active_cell_iterator>>
874
875 /// Internal data: stores map of vself ball cell Id to the closest atom Id
876 /// of that cell. Outer vector over vself bins.
877 std::vector<std::map<dealii::CellId, dftfe::uInt>>
879
880 /// Internal data: stores the map of atom Id (only in the local processor)
881 /// to the vself bin Id.
882 std::map<dftfe::uInt, dftfe::uInt> d_AtomIdBinIdLocalDofHandler;
883
884 /* Internal data: stores the face ids of dftPtr->d_dofHandler (single
885 * component field) on which to evaluate the vself ball surface integral in
886 * the configurational force expression. Outer vector is over the vself
887 * bins. Inner map is between the cell iterator and the vector of face ids
888 * to integrate on for that cell iterator.
889 */
890 std::vector<std::map<dealii::DoFHandler<3>::active_cell_iterator,
891 std::vector<dftfe::uInt>>>
893
894 /* Internal data: stores the face ids of d_dofHandlerForce (three component
895 * field) on which to evaluate the vself ball surface integral in the
896 * configurational force expression. Outer vector is over the vself bins.
897 * Inner map is between the cell iterator and the vector of face ids to
898 * integrate on for that cell iterator.
899 */
900 std::vector<std::map<dealii::DoFHandler<3>::active_cell_iterator,
901 std::vector<dftfe::uInt>>>
903
904 /// Internal data: stores cell iterators of all cells in
905 /// dftPtr->d_dofHandler which are part of the vself ball. Outer vector is
906 /// over vself bins.
907 std::vector<std::vector<dealii::DoFHandler<3>::active_cell_iterator>>
909
910 /// Internal data: stores cell iterators of all cells in d_dofHandlerForce
911 /// which are part of the vself ball. Outer vector is over vself bins.
912 std::vector<std::vector<dealii::DoFHandler<3>::active_cell_iterator>>
914
915 /// Internal data: stores map of vself ball cell Id to the closest atom Id
916 /// of that cell. Outer vector over vself bins.
917 std::vector<std::map<dealii::CellId, dftfe::uInt>>
919
920 /// Internal data: stores the map of atom Id (only in the local processor)
921 /// to the vself bin Id.
922 std::map<dftfe::uInt, dftfe::uInt> d_AtomIdBinIdLocalDofHandlerElectro;
923
924 /* Internal data: stores the face ids of dftPtr->d_dofHandler (single
925 * component field) on which to evaluate the vself ball surface integral in
926 * the configurational force expression. Outer vector is over the vself
927 * bins. Inner map is between the cell iterator and the vector of face ids
928 * to integrate on for that cell iterator.
929 */
930 std::vector<std::map<dealii::DoFHandler<3>::active_cell_iterator,
931 std::vector<dftfe::uInt>>>
933
934 /* Internal data: stores the face ids of d_dofHandlerForce (three component
935 * field) on which to evaluate the vself ball surface integral in the
936 * configurational force expression. Outer vector is over the vself bins.
937 * Inner map is between the cell iterator and the vector of face ids to
938 * integrate on for that cell iterator.
939 */
940 std::vector<std::map<dealii::DoFHandler<3>::active_cell_iterator,
941 std::vector<dftfe::uInt>>>
943
944 std::map<dealii::CellId, dealii::DoFHandler<3>::active_cell_iterator>
946
947 std::vector<distributedCPUVec<double>> d_gaussianWeightsVecAtoms;
948
949
950
951 /// domain decomposition mpi_communicator
952 const MPI_Comm d_mpiCommParent;
953
954 /// domain decomposition mpi_communicator
955 const MPI_Comm mpi_communicator;
956
958
959 /// number of mpi processes in the current pool
961
962 /// current mpi process id in the current pool
964
965 /// conditional stream object to enable printing only on root processor
966 /// across all pools
967 dealii::ConditionalOStream pcout;
968 };
969
970} // namespace dftfe
971#endif
Definition AtomicCenteredNonLocalOperator.h:58
This class is the primary interface location of all other parts of the DFT-FE code for all steps invo...
Definition dft.h:116
Namespace which declares the input parameters and the functions to parse them from the input paramete...
Definition dftParameters.h:36
Calculates dispersion correction to energy, force and stress.
Definition dftd.h:37
void createBinObjectsForce(const dealii::DoFHandler< 3 > &dofHandler, const dealii::DoFHandler< 3 > &dofHandlerForce, const dealii::AffineConstraints< double > &hangingPlusPBCConstraints, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManager, std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > &cellsVselfBallsDofHandler, std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > &cellsVselfBallsDofHandlerForce, std::vector< std::map< dealii::CellId, dftfe::uInt > > &cellsVselfBallsClosestAtomIdDofHandler, std::map< dftfe::uInt, dftfe::uInt > &AtomIdBinIdLocalDofHandler, std::vector< std::map< dealii::DoFHandler< 3 >::active_cell_iterator, std::vector< dftfe::uInt > > > &cellFacesVselfBallSurfacesDofHandler, std::vector< std::map< dealii::DoFHandler< 3 >::active_cell_iterator, std::vector< dftfe::uInt > > > &cellFacesVselfBallSurfacesDofHandlerForce)
std::map< std::pair< dftfe::uInt, dftfe::uInt >, dftfe::uInt > d_atomsForceDofsElectro
Definition force.h:862
const dftfe::uInt n_mpi_processes
number of mpi processes in the current pool
Definition force.h:960
std::vector< std::map< dealii::DoFHandler< 3 >::active_cell_iterator, std::vector< dftfe::uInt > > > d_cellFacesVselfBallSurfacesDofHandlerElectro
Definition force.h:932
std::map< std::pair< dftfe::uInt, dftfe::uInt >, dftfe::uInt > d_atomsForceDofs
Definition force.h:857
void computeAtomsForcesGaussianGenerator(bool allowGaussianOverlapOnAtoms=false)
dealii::DoFHandler< 3 > d_dofHandlerForceElectro
Definition force.h:821
std::vector< std::map< dealii::CellId, dftfe::uInt > > d_cellsVselfBallsClosestAtomIdDofHandlerElectro
Definition force.h:918
forceClass(dftClass< FEOrder, FEOrderElectro, memorySpace > *_dftPtr, const MPI_Comm &mpi_comm_parent, const MPI_Comm &mpi_comm_domain, const dftParameters &dftParams)
Constructor.
std::map< dftfe::uInt, dftfe::uInt > d_AtomIdBinIdLocalDofHandler
Definition force.h:882
const dftParameters & d_dftParams
Definition force.h:957
dealii::Tensor< 2, 3, double > & getStress()
returns a copy of the current stress tensor value.
std::vector< double > d_globalAtomsForces
Storage for configurational force on all global atoms.
Definition force.h:781
void FPSPLocalGammaAtomsElementalContribution(std::map< dftfe::uInt, std::vector< double > > &forceContributionFPSPLocalGammaAtoms, dealii::FEValues< 3 > &feValues, dealii::FEFaceValues< 3 > &feFaceValues, dealii::FEEvaluation< 3, 1, C_num1DQuadLPSP< FEOrder >() *C_numCopies1DQuadLPSP(), 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const dftfe::uInt phiTotDofHandlerIndexElectro, const dftfe::uInt cell, const dealii::AlignedVector< dealii::VectorizedArray< double > > &rhoQuads, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &gradRhoQuads, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtoms, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManager, const std::vector< std::map< dealii::CellId, dftfe::uInt > > &cellsVselfBallsClosestAtomIdDofHandler)
dftfe::uInt d_forceDofHandlerIndexElectro
Definition force.h:829
dealii::Tensor< 2, 3, double > d_stress
Storage for configurational stress tensor.
Definition force.h:784
std::vector< double > & getAtomsForces()
returns a copy of the configurational force on all global atoms.
void configForceLinFEInit(const dealii::MatrixFree< 3, double > &matrixFreeData, const dealii::MatrixFree< 3, double > &matrixFreeDataElectro)
void computeAtomsForces(const dealii::MatrixFree< 3, double > &matrixFreeData, const dispersionCorrection &dispersionCorr, const dftfe::uInt eigenDofHandlerIndex, const dftfe::uInt smearedChargeQuadratureId, const dftfe::uInt lpspQuadratureIdElectro, const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const dftfe::uInt phiTotDofHandlerIndexElectro, const distributedCPUVec< double > &phiTotRhoOutElectro, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &rhoOutValues, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &gradRhoOutValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValuesLpsp, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValuesLpsp, const std::map< dealii::CellId, std::vector< double > > &rhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &gradRhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &hessianRhoCoreValues, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms, const std::map< dealii::CellId, std::vector< double > > &pseudoVLocElectro, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtomsElectro, const dealii::AffineConstraints< double > &hangingPlusPBCConstraintsElectro, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro)
computes the configurational force on all atoms corresponding to a Gaussian generator,...
const bool d_allowGaussianOverlapOnAtoms
Definition force.h:797
void initUnmoved(const dealii::Triangulation< 3, 3 > &triangulation, const dealii::Triangulation< 3, 3 > &serialTriangulation, const std::vector< std::vector< double > > &domainBoundingVectors, const bool isElectrostaticsMesh)
initializes data structures inside forceClass assuming unmoved triangulation.
void accumulateForceContributionGammaAtomsFloating(const std::map< dftfe::uInt, std::vector< double > > &forceContributionLocalGammaAtoms, std::vector< double > &accumForcesVector)
void FNonlinearCoreCorrectionGammaAtomsElementalContribution(std::map< dftfe::uInt, std::vector< double > > &forceContributionFNonlinearCoreCorrectionGammaAtoms, dealii::FEEvaluation< 3, 1, C_num1DQuad< C_rhoNodalPolyOrder< FEOrder, FEOrderElectro >()>(), 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const dftfe::uInt cell, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &derExcGradRho, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms)
dftfe::uInt d_forceDofHandlerIndex
Definition force.h:825
dealii::ConditionalOStream pcout
Definition force.h:967
void addENonlinearCoreCorrectionStressContributionSpinPolarized(dealii::FEEvaluation< 3, 1, C_num1DQuad< C_rhoNodalPolyOrder< FEOrder, FEOrderElectro >()>(), 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const dftfe::uInt cell, const dealii::AlignedVector< dealii::VectorizedArray< double > > &vxcQuadsSpin0, const dealii::AlignedVector< dealii::VectorizedArray< double > > &vxcQuadsSpin1, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &derExcGradRhoSpin0, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &derExcGradRhoSpin1, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms, const bool isXCGGA=false)
void computeStressEEshelbyEPSPEnlEk(const dealii::MatrixFree< 3, double > &matrixFreeData, const dftfe::uInt eigenDofHandlerIndex, const dftfe::uInt smearedChargeQuadratureId, const dftfe::uInt lpspQuadratureIdElectro, const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const dftfe::uInt phiTotDofHandlerIndexElectro, const distributedCPUVec< double > &phiTotRhoOutElectro, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &rhoOutValues, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &gradRhoOutValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValuesLpsp, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValuesLpsp, const std::map< dealii::CellId, std::vector< double > > &pseudoVLocElectro, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtomsElectro, const std::map< dealii::CellId, std::vector< double > > &rhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &gradRhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &hessianRhoCoreValues, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro)
dealii::AffineConstraints< double > d_constraintsNoneForceElectro
Definition force.h:853
void printStress()
prints the currently stored configurational stress tensor.
void computeConfigurationalForceTotalLinFE(const dealii::MatrixFree< 3, double > &matrixFreeData, const dftfe::uInt eigenDofHandlerIndex, const dftfe::uInt smearedChargeQuadratureId, const dftfe::uInt lpspQuadratureIdElectro, const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const dftfe::uInt phiTotDofHandlerIndexElectro, const distributedCPUVec< double > &phiTotRhoOutElectro, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &rhoOutValues, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &gradRhoOutValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValuesLpsp, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValuesLpsp, const std::map< dealii::CellId, std::vector< double > > &rhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &gradRhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &hessianRhoCoreValues, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms, const std::map< dealii::CellId, std::vector< double > > &pseudoVLocElectro, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtomsElectro, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro)
void computeConfigurationalForceEselfNoSurfaceLinFE()
void computeConfigurationalForceEEshelbyTensorFPSPFnlLinFE(const dealii::MatrixFree< 3, double > &matrixFreeData, const dftfe::uInt eigenDofHandlerIndex, const dftfe::uInt smearedChargeQuadratureId, const dftfe::uInt lpspQuadratureIdElectro, const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const dftfe::uInt phiTotDofHandlerIndexElectro, const distributedCPUVec< double > &phiTotRhoOutElectro, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &rhoOutValues, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &gradRhoOutValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValuesLpsp, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValuesLpsp, const std::map< dealii::CellId, std::vector< double > > &rhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &gradRhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &hessianRhoCoreValues, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms, const std::map< dealii::CellId, std::vector< double > > &pseudoVLocElectro, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtomsElectro, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro)
const MPI_Comm mpi_communicator
domain decomposition mpi_communicator
Definition force.h:955
std::vector< std::map< dealii::DoFHandler< 3 >::active_cell_iterator, std::vector< dftfe::uInt > > > d_cellFacesVselfBallSurfacesDofHandlerForceElectro
Definition force.h:942
std::vector< std::map< dealii::DoFHandler< 3 >::active_cell_iterator, std::vector< dftfe::uInt > > > d_cellFacesVselfBallSurfacesDofHandler
Definition force.h:892
void FnlGammaxElementalContribution(dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &FVectQuads, const dealii::MatrixFree< 3, double > &matrixFreeData, const dftfe::uInt numQuadPoints, const std::shared_ptr< AtomicCenteredNonLocalOperator< dataTypes::number, memorySpace > > nonLocalOp, const dftfe::uInt numNonLocalAtomsCurrentProcess, const std::vector< dftfe::Int > &globalChargeIdNonLocalAtoms, const std::vector< dftfe::uInt > &numberPseudoWaveFunctionsPerAtom, const dftfe::uInt cell, const std::map< dealii::CellId, dftfe::uInt > &cellIdToCellNumberMap, const std::vector< dataTypes::number > &zetaDeltaVQuadsFlattened, const std::vector< dataTypes::number > &projectorKetTimesPsiTimesVTimesPartOccContractionGradPsiQuadsFlattened)
void computeConfigurationalForceEselfLinFE(const dealii::DoFHandler< 3 > &dofHandlerElectro, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro, const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const dftfe::uInt smearedChargeQuadratureId)
void addEVselfSmearedStressContribution(dealii::FEEvaluation< 3, -1, 1, 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const dftfe::uInt cell, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &gradVselfQuads, const std::vector< dftfe::uInt > &nonTrivialAtomImageIdsMacroCell, const std::map< dealii::CellId, std::vector< dftfe::Int > > &bQuadAtomIdsAllAtomsImages, const dealii::AlignedVector< dealii::VectorizedArray< double > > &smearedbQuads)
void FVselfSmearedChargesGammaAtomsElementalContribution(std::map< dftfe::uInt, std::vector< double > > &forceContributionSmearedChargesGammaAtoms, dealii::FEEvaluation< 3, -1, 1, 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const dftfe::uInt cell, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &gradVselfBinQuads, const std::vector< dftfe::uInt > &nonTrivialAtomIdsMacroCell, const std::map< dealii::CellId, std::vector< dftfe::Int > > &bQuadAtomIdsAllAtoms, const dealii::AlignedVector< dealii::VectorizedArray< double > > &smearedbQuads)
std::map< dealii::CellId, dealii::DoFHandler< 3 >::active_cell_iterator > d_cellIdToActiveCellIteratorMapDofHandlerRhoNodalElectro
Definition force.h:945
const dftfe::uInt this_mpi_process
current mpi process id in the current pool
Definition force.h:963
distributedCPUVec< double > d_configForceVectorLinFEElectro
Definition force.h:754
const MPI_Comm d_mpiCommParent
domain decomposition mpi_communicator
Definition force.h:952
std::map< dftfe::uInt, dftfe::uInt > d_AtomIdBinIdLocalDofHandlerElectro
Definition force.h:922
void computeFloatingAtomsForces()
dealii::AffineConstraints< double > d_constraintsNoneForce
Definition force.h:849
std::vector< double > d_forceAtomsFloating
Definition force.h:766
dealii::FESystem< 3 > FEForce
Definition force.h:804
void distributeForceContributionFnlGammaAtoms(const std::map< dftfe::uInt, std::vector< double > > &forceContributionFnlGammaAtoms)
void FNonlinearCoreCorrectionGammaAtomsElementalContribution(std::map< dftfe::uInt, std::vector< double > > &forceContributionFNonlinearCoreCorrectionGammaAtoms, dealii::FEEvaluation< 3, 1, C_num1DQuad< C_rhoNodalPolyOrder< FEOrder, FEOrderElectro >()>(), 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const dftfe::uInt cell, const dealii::AlignedVector< dealii::VectorizedArray< double > > &vxcQuads, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms)
void stressEnlElementalContribution(dealii::Tensor< 2, 3, double > &stressContribution, const dealii::MatrixFree< 3, double > &matrixFreeData, const dftfe::uInt numQuadPoints, const std::vector< double > &jxwQuadsSubCells, const dftfe::uInt cell, const dftfe::uInt numNonLocalAtomsCurrentProcess, const std::shared_ptr< AtomicCenteredNonLocalOperator< dataTypes::number, memorySpace > > nonLocalOp, const std::vector< dftfe::uInt > &numberPseudoWaveFunctionsPerAtom, const std::map< dealii::CellId, dftfe::uInt > &cellIdToCellNumberMap, const std::vector< dataTypes::number > &zetalmDeltaVlProductDistImageAtoms, const std::vector< dataTypes::number > &projectorKetTimesPsiTimesVTimesPartOccContractionGradPsiQuadsFlattened, const bool isSpinPolarized)
void computeStress(const dealii::MatrixFree< 3, double > &matrixFreeData, const dispersionCorrection &dispersionCorr, const dftfe::uInt eigenDofHandlerIndex, const dftfe::uInt smearedChargeQuadratureId, const dftfe::uInt lpspQuadratureIdElectro, const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const dftfe::uInt phiTotDofHandlerIndexElectro, const distributedCPUVec< double > &phiTotRhoOutElectro, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &rhoOutValues, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &gradRhoOutValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValuesLpsp, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValuesLpsp, const std::map< dealii::CellId, std::vector< double > > &pseudoVLocElectro, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtomsElectro, const std::map< dealii::CellId, std::vector< double > > &rhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &gradRhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &hessianRhoCoreValues, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms, const dealii::AffineConstraints< double > &hangingPlusPBCConstraintsElectro, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro)
Update force generator Gaussian constant.
dealii::IndexSet d_locally_owned_dofsForceElectro
Definition force.h:837
void FnlGammaAtomsElementalContribution(std::map< dftfe::uInt, std::vector< double > > &forceContributionFnlGammaAtoms, const dealii::MatrixFree< 3, double > &matrixFreeData, dealii::FEEvaluation< 3, 1, C_num1DQuadNLPSP< FEOrder >() *C_numCopies1DQuadNLPSP(), 3 > &forceEvalNLP, const std::shared_ptr< AtomicCenteredNonLocalOperator< dataTypes::number, memorySpace > > nonLocalOp, dftfe::uInt numNonLocalAtomsCurrentProcess, const std::vector< dftfe::Int > &globalChargeIdNonLocalAtoms, const std::vector< dftfe::uInt > &numberPseudoWaveFunctionsPerAtom, const dftfe::uInt cell, const std::map< dealii::CellId, dftfe::uInt > &cellIdToCellNumberMap, const std::vector< dataTypes::number > &zetaDeltaVQuadsFlattened, const std::vector< dataTypes::number > &projectorKetTimesPsiTimesVTimesPartOccContractionGradPsiQuadsFlattened)
void FPhiTotSmearedChargesGammaAtomsElementalContribution(std::map< dftfe::uInt, std::vector< double > > &forceContributionSmearedChargesGammaAtoms, dealii::FEEvaluation< 3, -1, 1, 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const dftfe::uInt cell, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &gradPhiTotQuads, const std::vector< dftfe::uInt > &nonTrivialAtomIdsMacroCell, const std::map< dealii::CellId, std::vector< dftfe::Int > > &bQuadAtomIdsAllAtoms, const dealii::AlignedVector< dealii::VectorizedArray< double > > &smearedbQuads)
void computeConfigurationalForceEEshelbyEElectroPhiTot(const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const dftfe::uInt phiTotDofHandlerIndexElectro, const dftfe::uInt smearedChargeQuadratureId, const dftfe::uInt lpspQuadratureIdElectro, const distributedCPUVec< double > &phiTotRhoOutElectro, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValuesLpsp, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValuesLpsp, const std::map< dealii::CellId, std::vector< double > > &pseudoVLocElectro, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtomsElectro, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro)
dealii::DoFHandler< 3 > d_dofHandlerForce
Definition force.h:812
void addEPSPStressContribution(dealii::FEValues< 3 > &feValues, dealii::FEFaceValues< 3 > &feFaceValues, dealii::FEEvaluation< 3, 1, C_num1DQuadLPSP< FEOrder >() *C_numCopies1DQuadLPSP(), 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const dftfe::uInt phiTotDofHandlerIndexElectro, const dftfe::uInt cell, const dealii::AlignedVector< dealii::VectorizedArray< double > > &rhoQuads, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &gradRhoQuads, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtoms, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManager, const std::vector< std::map< dealii::CellId, dftfe::uInt > > &cellsVselfBallsClosestAtomIdDofHandler)
dealii::IndexSet d_locally_owned_dofsForce
Definition force.h:833
std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > d_cellsVselfBallsDofHandlerForceElectro
Definition force.h:913
void computeElementalNonLocalPseudoOVDataForce()
void FNonlinearCoreCorrectionGammaAtomsElementalContributionSpinPolarized(std::map< dftfe::uInt, std::vector< double > > &forceContributionFNonlinearCoreCorrectionGammaAtoms, dealii::FEEvaluation< 3, 1, C_num1DQuad< C_rhoNodalPolyOrder< FEOrder, FEOrderElectro >()>(), 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const dftfe::uInt cell, const dealii::AlignedVector< dealii::VectorizedArray< double > > &vxcQuadsSpin0, const dealii::AlignedVector< dealii::VectorizedArray< double > > &vxcQuadsSpin1, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &derExcGradRhoSpin0, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &derExcGradRhoSpin1, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms, const bool isXCGGA=false)
void computeConfigurationalForcePhiExtLinFE()
std::vector< std::map< dealii::CellId, dftfe::uInt > > d_cellsVselfBallsClosestAtomIdDofHandler
Definition force.h:878
void distributeForceContributionFPSPLocalGammaAtoms(const std::map< dftfe::uInt, std::vector< double > > &forceContributionFPSPLocalGammaAtoms, const std::map< std::pair< dftfe::uInt, dftfe::uInt >, dftfe::uInt > &atomsForceDofs, const dealii::AffineConstraints< double > &constraintsNoneForce, distributedCPUVec< double > &configForceVectorLinFE)
std::vector< distributedCPUVec< double > > d_gaussianWeightsVecAtoms
Definition force.h:947
dftClass< FEOrder, FEOrderElectro, memorySpace > * dftPtr
pointer to dft class
Definition force.h:800
void addEPhiTotSmearedStressContribution(dealii::FEEvaluation< 3, -1, 1, 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const dftfe::uInt cell, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &gradPhiTotQuads, const std::vector< dftfe::uInt > &nonTrivialAtomImageIdsMacroCell, const std::map< dealii::CellId, std::vector< dftfe::Int > > &bQuadAtomIdsAllAtomsImages, const dealii::AlignedVector< dealii::VectorizedArray< double > > &smearedbQuads)
void locateAtomCoreNodesForce(const dealii::DoFHandler< 3 > &dofHandlerForce, const dealii::IndexSet &locally_owned_dofsForce, std::map< std::pair< dftfe::uInt, dftfe::uInt >, dftfe::uInt > &atomsForceDofs)
get the value of Gaussian generator parameter (d_gaussianConstant). Gaussian generator: Gamma(r)= exp...
std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > d_cellsVselfBallsDofHandler
Definition force.h:868
dealii::IndexSet d_locally_relevant_dofsForceElectro
Definition force.h:845
void initMoved(std::vector< const dealii::DoFHandler< 3 > * > &dofHandlerVectorMatrixFree, std::vector< const dealii::AffineConstraints< double > * > &constraintsVectorMatrixFree, const bool isElectrostaticsMesh)
initializes data structures inside forceClass which depend on the moved mesh.
void computeStressEself(const dealii::DoFHandler< 3 > &dofHandlerElectro, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro, const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const dftfe::uInt smearedChargeQuadratureId)
std::vector< std::map< dealii::DoFHandler< 3 >::active_cell_iterator, std::vector< dftfe::uInt > > > d_cellFacesVselfBallSurfacesDofHandlerForce
Definition force.h:902
void initPseudoData()
initializes and precomputes pseudopotential related data structuers required for configurational forc...
std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > d_cellsVselfBallsDofHandlerElectro
Definition force.h:908
dealii::IndexSet d_locally_relevant_dofsForce
Definition force.h:841
void addENonlinearCoreCorrectionStressContribution(dealii::FEEvaluation< 3, 1, C_num1DQuad< C_rhoNodalPolyOrder< FEOrder, FEOrderElectro >()>(), 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const dftfe::uInt cell, const dealii::AlignedVector< dealii::VectorizedArray< double > > &vxcQuads, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &derExcGradRho, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms)
void computeStressEEshelbyEElectroPhiTot(const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const dftfe::uInt phiTotDofHandlerIndexElectro, const dftfe::uInt smearedChargeQuadratureId, const dftfe::uInt lpspQuadratureIdElectro, const distributedCPUVec< double > &phiTotRhoOutElectro, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValuesLpsp, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValuesElectro, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValuesElectroLpsp, const std::map< dealii::CellId, std::vector< double > > &pseudoVLocElectro, const std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtomsElectro, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro)
distributedCPUVec< double > d_configForceVectorLinFE
Definition force.h:747
std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > d_cellsVselfBallsDofHandlerForce
Definition force.h:873
dealii::Tensor< 2, 3, double > d_stressKPoints
Definition force.h:792
void configForceLinFEFinalize()
void printAtomsForces()
prints the currently stored configurational forces on atoms and the Gaussian generator constant used ...
Definition MemoryStorage.h:33
Categorizes atoms into bins for efficient solution of nuclear electrostatic self-potential....
Definition vselfBinsManager.h:36
MemorySpace
Definition MemorySpaceType.h:33
Definition pseudoPotentialToDftfeConverter.cc:34
constexpr dftfe::uInt C_num1DQuad()
1d quadrature rule order
Definition constants.h:41
dealii::LinearAlgebra::distributed::Vector< elem_type, dealii::MemorySpace::Host > distributedCPUVec
Definition headers.h:92
std::uint32_t uInt
Definition TypeConfig.h:10
constexpr dftfe::uInt C_numCopies1DQuadNLPSP()
number of copies 1d quad rule non-local PSP
Definition constants.h:149
constexpr dftfe::uInt C_numCopies1DQuadLPSP()
number of copies 1d quad rule local PSP
Definition constants.h:164
constexpr dftfe::uInt C_num1DQuadLPSP()
1d quadrature rule order for local part of pseudopotential
Definition constants.h:157
constexpr dftfe::uInt C_rhoNodalPolyOrder()
rho nodal polynomial order
Definition constants.h:134
constexpr dftfe::uInt C_num1DQuadNLPSP()
1d quadrature rule order for non-local part of pseudopotential
Definition constants.h:142