45 const MPI_Comm &mpi_comm_domain,
46 const MPI_Comm &mpi_intercomm_kpts,
68 std::vector<
const dealii::AffineConstraints<double> *> &constraintsVector,
69 const dealii::AffineConstraints<double> &onlyHangingNodeConstraints,
70 const dealii::DoFHandler<3> &dofHandler,
71 const dealii::AffineConstraints<double> &constraintMatrix,
72 const std::vector<std::vector<double>> &atomLocations,
73 const std::vector<std::vector<double>> &imagePositions,
74 const std::vector<dftfe::Int> &imageIds,
75 const std::vector<double> &imageCharges,
76 const double radiusAtomBall);
95 std::vector<
const dealii::AffineConstraints<double> *> &constraintsVector,
96 const dealii::AffineConstraints<double> &onlyHangingNodeConstraints,
97 const dealii::DoFHandler<3> &dofHandler,
98 const dealii::AffineConstraints<double> &constraintMatrix,
99 const std::vector<std::vector<double>> &atomLocations,
100 const std::vector<std::vector<double>> &imagePositions,
101 const std::vector<dftfe::Int> &imageIds,
102 const std::vector<double> &imageCharges,
103 const bool vselfPerturbationUpdateForStress =
false);
123 const std::shared_ptr<
125 FEBasisOperations<double, double, dftfe::utils::MemorySpace::HOST>>
129 const dealii::AffineConstraints<double> &hangingPeriodicConstraintMatrix,
130 const std::vector<std::vector<double>> &imagePositions,
131 const std::vector<dftfe::Int> &imageIds,
132 const std::vector<double> &imageCharges,
133 std::vector<std::vector<double>> &localVselfs,
134 std::map<dealii::CellId, std::vector<double>> &bQuadValuesAllAtoms,
135 std::map<dealii::CellId, std::vector<dftfe::Int>> &bQuadAtomIdsAllAtoms,
136 std::map<dealii::CellId, std::vector<dftfe::Int>>
137 &bQuadAtomIdsAllAtomsImages,
138 std::map<dealii::CellId, std::vector<dftfe::uInt>>
139 &bCellNonTrivialAtomIds,
140 std::vector<std::map<dealii::CellId, std::vector<dftfe::uInt>>>
141 &bCellNonTrivialAtomIdsBins,
142 std::map<dealii::CellId, std::vector<dftfe::uInt>>
143 &bCellNonTrivialAtomImageIds,
144 std::vector<std::map<dealii::CellId, std::vector<dftfe::uInt>>>
145 &bCellNonTrivialAtomImageIdsBins,
146 const std::vector<double> &smearingWidths,
147 std::vector<double> &smearedChargeScaling,
149 const bool useSmearedCharges =
false,
150 const bool isVselfPerturbationSolve =
false);
152# ifdef DFTFE_WITH_DEVICE
169 solveVselfInBinsDevice(
170 const std::shared_ptr<
172 FEBasisOperations<double, double, dftfe::utils::MemorySpace::HOST>>
179 &cellGradNIGradNJIntergralDevice,
180 const std::shared_ptr<
183 const dealii::AffineConstraints<double> &hangingPeriodicConstraintMatrix,
184 const std::vector<std::vector<double>> &imagePositions,
185 const std::vector<dftfe::Int> &imageIds,
186 const std::vector<double> &imageCharges,
187 std::vector<std::vector<double>> &localVselfs,
188 std::map<dealii::CellId, std::vector<double>> &bQuadValuesAllAtoms,
189 std::map<dealii::CellId, std::vector<dftfe::Int>> &bQuadAtomIdsAllAtoms,
190 std::map<dealii::CellId, std::vector<dftfe::Int>>
191 &bQuadAtomIdsAllAtomsImages,
192 std::map<dealii::CellId, std::vector<dftfe::uInt>>
193 &bCellNonTrivialAtomIds,
194 std::vector<std::map<dealii::CellId, std::vector<dftfe::uInt>>>
195 &bCellNonTrivialAtomIdsBins,
196 std::map<dealii::CellId, std::vector<dftfe::uInt>>
197 &bCellNonTrivialAtomImageIds,
198 std::vector<std::map<dealii::CellId, std::vector<dftfe::uInt>>>
199 &bCellNonTrivialAtomImageIdsBins,
200 const std::vector<double> &smearingWidths,
201 std::vector<double> &smearedChargeScaling,
203 const bool useSmearedCharges =
false,
204 const bool isVselfPerturbationSolve =
false);
215 const std::shared_ptr<
217 FEBasisOperations<double, double, dftfe::utils::MemorySpace::HOST>>
222# ifdef DFTFE_WITH_DEVICE
225 &cellGradNIGradNJIntergralDevice,
226 const std::shared_ptr<
230 const dealii::AffineConstraints<double> &hangingPeriodicConstraintMatrix,
231 const std::vector<std::vector<double>> &imagePositions,
232 const std::vector<dftfe::Int> &imageIds,
233 const std::vector<double> &imageCharges,
234 const std::vector<double> &smearingWidths,
236 const bool useSmearedCharges =
false);
239 const std::map<dftfe::Int, std::set<dftfe::Int>> &
243 const std::map<dftfe::Int, std::set<dftfe::Int>> &
248 const std::vector<std::map<dealii::types::global_dof_index, dftfe::Int>> &
253 const std::vector<std::map<dealii::types::global_dof_index, dftfe::Int>> &
258 const std::vector<std::map<dealii::types::global_dof_index, dftfe::Int>> &
264 std::map<dealii::types::global_dof_index, dealii::Point<3>>> &
268 const std::vector<distributedCPUVec<double>> &
272 const std::vector<distributedCPUVec<double>> &
278 const std::vector<distributedCPUVec<double>> &
282 const std::map<dftfe::uInt, dftfe::uInt> &
304 const dealii::DoFHandler<3> &dofHandler,
305 const dealii::AffineConstraints<double> &onlyHangingNodeConstraints);
318 std::map<dftfe::Int, std::set<dftfe::Int>>
d_bins;
325 std::vector<std::map<dealii::types::global_dof_index, dftfe::Int>>
330 std::vector<std::map<dealii::types::global_dof_index, dftfe::Int>>
334 std::vector<std::map<dealii::types::global_dof_index, dealii::Point<3>>>
338 std::vector<std::map<dealii::types::global_dof_index, double>>
342 std::vector<std::map<dealii::types::global_dof_index, dftfe::Int>>
364 std::vector<std::map<dealii::types::global_dof_index, double>>
d_atomsInBin;
void solveVselfInBins(const std::shared_ptr< dftfe::basis::FEBasisOperations< double, double, dftfe::utils::MemorySpace::HOST > > &basisOperationsPtr, const dftfe::uInt offset, const dftfe::uInt matrixFreeQuadratureIdAX, const dealii::AffineConstraints< double > &hangingPeriodicConstraintMatrix, const std::vector< std::vector< double > > &imagePositions, const std::vector< dftfe::Int > &imageIds, const std::vector< double > &imageCharges, std::vector< std::vector< double > > &localVselfs, std::map< dealii::CellId, std::vector< double > > &bQuadValuesAllAtoms, std::map< dealii::CellId, std::vector< dftfe::Int > > &bQuadAtomIdsAllAtoms, std::map< dealii::CellId, std::vector< dftfe::Int > > &bQuadAtomIdsAllAtomsImages, std::map< dealii::CellId, std::vector< dftfe::uInt > > &bCellNonTrivialAtomIds, std::vector< std::map< dealii::CellId, std::vector< dftfe::uInt > > > &bCellNonTrivialAtomIdsBins, std::map< dealii::CellId, std::vector< dftfe::uInt > > &bCellNonTrivialAtomImageIds, std::vector< std::map< dealii::CellId, std::vector< dftfe::uInt > > > &bCellNonTrivialAtomImageIdsBins, const std::vector< double > &smearingWidths, std::vector< double > &smearedChargeScaling, const dftfe::uInt smearedChargeQuadratureId, const bool useSmearedCharges=false, const bool isVselfPerturbationSolve=false)
Solve nuclear electrostatic self-potential of atoms in each bin one-by-one.