DFT-FE 1.3.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
vectorUtilities.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 vectorUtilities_h
19#define vectorUtilities_h
20
21#include <headers.h>
22#include <operator.h>
23
24
25
26namespace dftfe
27{
28 /**
29 * @brief Contains generic utils functions related to custom partitioned flattened dealii vector
30 *
31 * @author Phani Motamarri, Sambit Das
32 */
33 namespace vectorTools
34 {
35 void
37 const dealii::DoFHandler<3> &dofHandlerPar,
38 dealii::AffineConstraints<double> &constraints);
39
40 /** @brief Create constraint matrix using serial mesh.
41 * Temporary fix for a bug (Issue #7053) in deal.ii until it is resolved.
42 *
43 * @param[in] serial Triangulation which must be exactly same as the
44 * parallel triangulation associated with dofHandlerPar
45 * @param[in] parallel DofHandler
46 * @param[out] periodic hanging constraints.
47 * @param[out] only hanging constraints
48 */
49 void
51 const dealii::Triangulation<3, 3> &serTria,
52 const dealii::DoFHandler<3> &dofHandlerPar,
53 const MPI_Comm &mpi_comm_parent,
54 const MPI_Comm &mpi_comm_domain,
55 const std::vector<std::vector<double>> &domainBoundingVectors,
56 dealii::AffineConstraints<double> &periodicHangingConstraints,
57 dealii::AffineConstraints<double> &onlyHangingConstraints,
58 const dftfe::Int verbosity,
59 const bool periodicX,
60 const bool periodicY,
61 const bool periodicZ);
62
63
64 /** @brief Creates a custom partitioned flattened dealii vector.
65 * stores multiple components asociated with a node sequentially.
66 *
67 * @param partitioner associated with single component vector
68 * @param blockSize number of components associated with each node
69 *
70 * @return flattenedArray custom partitioned dealii vector
71 */
72 template <typename T>
73 void
75 const std::shared_ptr<const dealii::Utilities::MPI::Partitioner>
76 &partitioner,
77 const dftfe::uInt blockSize,
78 distributedCPUVec<T> &flattenedArray);
79
80
81
82 /** @brief Creates a cell local index set map for flattened array
83 *
84 * @param partitioner associated with the flattened array
85 * @param matrix_free_data object pointer associated with the matrix free data structure
86 * @param blockSize number of components associated with each node
87 *
88 * @return flattenedArrayMacroCellLocalProcIndexId macrocell's subcell local proc index map
89 * @return flattenedArrayCellLocalProcIndexId cell local proc index map
90 */
91 void
93 const std::shared_ptr<
95 &partitioner,
96 const dealii::MatrixFree<3, double> &matrix_free_data,
97 const dftfe::uInt mfDofHandlerIndex,
98 const dftfe::uInt blockSize,
99 std::vector<std::vector<dftfe::uInt>>
100 &flattenedArrayMacroCellLocalProcIndexId,
101 std::vector<std::vector<dftfe::uInt>>
102 &flattenedArrayCellLocalProcIndexId);
103
104
105 /** @brief Creates a cell local index set map for flattened array
106 *
107 * @param partitioner associated with the flattened array
108 * @param matrix_free_data object pointer associated with the matrix free data structure
109 * @param blockSize number of components associated with each node
110 *
111 * @return flattenedArrayMacroCellLocalProcIndexId macrocell's subcell local proc index map
112 * @return flattenedArrayCellLocalProcIndexId cell local proc index map
113 */
114 void
116 const std::shared_ptr<
118 &partitioner,
119 const dealii::MatrixFree<3, double> &matrix_free_data,
120 const dftfe::uInt mfDofHandlerIndex,
121 const dftfe::uInt blockSize,
122 std::vector<dftfe::uInt> &flattenedArrayMacroCellLocalProcIndexId,
123 std::vector<dftfe::uInt> &normalCellIdToMacroCellIdMap,
124 std::vector<dftfe::uInt> &macroCellIdToNormalCellIdMap,
125 std::vector<dftfe::uInt> &flattenedArrayCellLocalProcIndexId);
126
127 /** @brief Creates a cell local index set map for flattened array
128 *
129 * @param partitioner associated with the flattened array
130 * @param matrix_free_data object pointer associated with the matrix free data structure
131 * @param blockSize number of components associated with each node
132 *
133 * @return flattenedArrayMacroCellLocalProcIndexId macrocell's subcell local proc index map
134 * @return flattenedArrayCellLocalProcIndexId cell local proc index map
135 */
136 void
138 const std::shared_ptr<
140 &partitioner,
141 const dealii::MatrixFree<3, double> &matrix_free_data,
142 const dftfe::uInt mfDofHandlerIndex,
143 const dftfe::uInt blockSize,
144 std::vector<dftfe::uInt> &flattenedArrayCellLocalProcIndexId);
145
146
147#ifdef USE_COMPLEX
148 /** @brief Copies a single field component from a flattenedArray STL
149 * vector containing multiple component fields to a 2-component field (real
150 * and complex) parallel distributed vector.
151 *
152 * @param[in] flattenedArray flattened STL vector with multiple component
153 * fields
154 * @param[in] totalNumberComponents total number of component fiels in
155 * flattenedArray
156 * @param[in] componentIndexRange desired range field components
157 * [componentIndexRange.first,componentIndexRange.second)
158 * @param[in] localProcDofIndicesReal local dof indices in the current
159 * processor which correspond to component-1 of 2-component parallel
160 * distributed array
161 * @param[in] localProcDofIndicesImag local dof indices in the current
162 * processor which correspond to component-2 of 2-component parallel
163 * distributed array
164 * @param[out] componentVectors vector of two component field parallel
165 * distributed vectors with the values corresponding to fields of
166 * componentIndexRange of flattenedArray. componentVectors is expected to be
167 * of the size componentIndexRange.second-componentIndexRange.first.
168 * Further, each entry of componentVectors is assumed to be already
169 * initialized with the 2-component version of the same single component
170 * partitioner used in the creation of the flattenedArray partitioner.
171 */
172 void
174 const std::complex<double> *flattenedArray,
175 const dftfe::uInt totalNumberComponents,
176 const dftfe::uInt localVectorSize,
177 const std::pair<dftfe::uInt, dftfe::uInt> componentIndexRange,
178 const std::vector<dealii::types::global_dof_index>
179 &localProcDofIndicesReal,
180 const std::vector<dealii::types::global_dof_index>
181 &localProcDofIndicesImag,
182 std::vector<distributedCPUVec<double>> &componentVectors);
183
184 void
186 const std::complex<double> *flattenedArray,
187 const dftfe::uInt totalNumberComponents,
188 const dftfe::uInt localVectorSize,
189 const std::pair<dftfe::uInt, dftfe::uInt> componentIndexRange,
190
191 std::vector<distributedCPUVec<double>> &componentVectors);
192
193#else
194 /** @brief Copies a single field component from a flattenedArray STL
195 * vector containing multiple component fields to a single field parallel
196 * distributed vector.
197 *
198 * @param[in] flattenedArray flattened STL vector with multiple component
199 * fields
200 * @param[in] totalNumberComponents total number of component fiels in
201 * flattenedArray
202 * @param[in] componentIndexRange desired range field components
203 * [componentIndexRange.first,componentIndexRange.second)
204 * @param[out] componentVectors vector of parallel distributed vectors with
205 * fields corresponding to componentIndexRange. componentVectors is expected
206 * to be of the size componentIndexRange.second-componentIndexRange.first.
207 * Further, each entry of componentVectors is assumed to be already
208 * initialized with the same single component partitioner used in the
209 * creation of the flattenedArray partitioner.
210 */
211 void
213 const double *flattenedArray,
214 const dftfe::uInt totalNumberComponents,
215 const dftfe::uInt localVectorSize,
216 const std::pair<dftfe::uInt, dftfe::uInt> componentIndexRange,
217 std::vector<distributedCPUVec<double>> &componentVectors);
218
219#endif
220
221#ifdef USE_COMPLEX
222 /** @brief Copies a single field component from a flattenedArray parallel distributed
223 * vector containing multiple component fields to a 2-component field (real
224 * and complex) parallel distributed vector.
225 *
226 * @param[in] flattenedArray flattened parallel distributed vector with
227 * multiple component fields
228 * @param[in] totalNumberComponents total number of component fiels in
229 * flattenedArray
230 * @param[in] componentIndexRange desired range field components
231 * [componentIndexRange.first,componentIndexRange.second)
232 * @param[in] localProcDofIndicesReal local dof indices in the current
233 * processor which correspond to component-1 of 2-component parallel
234 * distributed array
235 * @param[in] localProcDofIndicesImag local dof indices in the current
236 * processor which correspond to component-2 of 2-component parallel
237 * distributed array
238 * @param[out] componentVectors vector of two component field parallel
239 * distributed vectors with the values corresponding to fields of
240 * componentIndexRange of flattenedArray. componentVectors is expected to be
241 * of the size componentIndexRange.second-componentIndexRange.first.
242 * Further, each entry of componentVectors is assumed to be already
243 * initialized with the 2-component version of the same single component
244 * partitioner used in the creation of the flattenedArray partitioner.
245 * @param[in] isFlattenedDealiiGhostValuesUpdated default is false. Use
246 * true for optimization if update ghost values has already been called in
247 * the flattened dealii vec.
248 */
249 void
251 const distributedCPUVec<std::complex<double>> &flattenedArray,
252 const dftfe::uInt totalNumberComponents,
253 const std::pair<dftfe::uInt, dftfe::uInt> componentIndexRange,
254 const std::vector<dealii::types::global_dof_index>
255 &localProcDofIndicesReal,
256 const std::vector<dealii::types::global_dof_index>
257 &localProcDofIndicesImag,
258 std::vector<distributedCPUVec<double>> &componentVectors,
259 const bool isFlattenedDealiiGhostValuesUpdated = false);
260
261#else
262 /** @brief Copies a single field component from a flattenedArray parallel distributed
263 * vector containing multiple component fields to a single field parallel
264 * distributed vector.
265 *
266 * @param[in] flattenedArray flattened parallel distributed vector with
267 * multiple component fields
268 * @param[in] totalNumberComponents total number of component fiels in
269 * flattenedArray
270 * @param[in] componentIndexRange desired range field components
271 * [componentIndexRange.first,componentIndexRange.second)
272 * @param[out] componentVectors vector of parallel distributed vectors with
273 * fields corresponding to componentIndexRange. componentVectors is expected
274 * to be of the size componentIndexRange.second-componentIndexRange.first.
275 * Further, each entry of componentVectors is assumed to be already
276 * initialized with the same single component partitioner used in the
277 * creation of the flattenedArray partitioner.
278 * @param[in] isFlattenedDealiiGhostValuesUpdated default is false. Use
279 * true for optimization if update ghost values has already been called in
280 * the flattened dealii vec.
281 */
282 void
284 const distributedCPUVec<double> &flattenedArray,
285 const dftfe::uInt totalNumberComponents,
286 const std::pair<dftfe::uInt, dftfe::uInt> componentIndexRange,
287 std::vector<distributedCPUVec<double>> &componentVectors,
288 const bool isFlattenedDealiiGhostValuesUpdated = false);
289
290#endif
291
292#ifdef USE_COMPLEX
293 /** @brief Copies to a flattenedArray parallel distributed
294 * vector containing multiple component fields from a 2-component field
295 * (real and complex) parallel distributed vector.
296 *
297 * @param[out] flattenedArray flattened parallel distributed vector with
298 * multiple component fields
299 * @param[in] totalNumberComponents total number of component fiels in
300 * flattenedArray
301 * @param[in] componentIndexRange desired range field components
302 * [componentIndexRange.first,componentIndexRange.second)
303 * @param[in] localProcDofIndicesReal local dof indices in the current
304 * processor which correspond to component-1 of 2-component parallel
305 * distributed array
306 * @param[in] localProcDofIndicesImag local dof indices in the current
307 * processor which correspond to component-2 of 2-component parallel
308 * distributed array
309 * @param[in] componentVectors vector of two component field parallel
310 * distributed vectors with the values corresponding to fields of
311 * componentIndexRange of flattenedArray. componentVectors is expected to be
312 * of the size componentIndexRange.second-componentIndexRange.first.
313 * Further, each entry of componentVectors is assumed to be already
314 * initialized with the 2-component version of the same single component
315 * partitioner used in the creation of the flattenedArray partitioner.
316 */
317 void
319 distributedCPUVec<std::complex<double>> &flattenedArray,
320 const dftfe::uInt totalNumberComponents,
321 const std::pair<dftfe::uInt, dftfe::uInt> componentIndexRange,
322 const std::vector<dealii::types::global_dof_index>
323 &localProcDofIndicesReal,
324 const std::vector<dealii::types::global_dof_index>
325 &localProcDofIndicesImag,
326 const std::vector<distributedCPUVec<double>> &componentVectors);
327
328#else
329 /** @brief Copies to a flattenedArray parallel distributed
330 * vector containing multiple component fields from a single field parallel
331 * distributed vector.
332 *
333 * @param[out] flattenedArray flattened parallel distributed vector with
334 * multiple component fields
335 * @param[in] totalNumberComponents total number of component fiels in
336 * flattenedArray
337 * @param[in] componentIndexRange desired range field components
338 * [componentIndexRange.first,componentIndexRange.second)
339 * @param[in] componentVectors vector of parallel distributed vectors with
340 * fields corresponding to componentIndexRange. componentVectors is expected
341 * to be of the size componentIndexRange.second-componentIndexRange.first.
342 * Further, each entry of componentVectors is assumed to be already
343 * initialized with the same single component partitioner used in the
344 * creation of the flattenedArray partitioner.
345 */
346 void
348 distributedCPUVec<double> &flattenedArray,
349 const dftfe::uInt totalNumberComponents,
350 const std::pair<dftfe::uInt, dftfe::uInt> componentIndexRange,
351 const std::vector<distributedCPUVec<double>> &componentVectors);
352
353#endif
354
355#ifdef USE_COMPLEX
356 /** @brief Copies to a flattenedArray stl
357 * vector containing multiple component fields from a 2-component field
358 * (real and complex) parallel distributed vector.
359 *
360 * @param[out] flattenedArray flattened stl vector with multiple component
361 * fields
362 * @param[in] totalNumberComponents total number of component fiels in
363 * flattenedArray
364 * @param[in] componentIndexRange desired range field components
365 * [componentIndexRange.first,componentIndexRange.second)
366 * @param[in] localProcDofIndicesReal local dof indices in the current
367 * processor which correspond to component-1 of 2-component parallel
368 * distributed array
369 * @param[in] localProcDofIndicesImag local dof indices in the current
370 * processor which correspond to component-2 of 2-component parallel
371 * distributed array
372 * @param[in] componentVectors vector of two component field parallel
373 * distributed vectors with the values corresponding to fields of
374 * componentIndexRange of flattenedArray. componentVectors is expected to be
375 * of the size componentIndexRange.second-componentIndexRange.first.
376 * Further, each entry of componentVectors is assumed to be already
377 * initialized with the 2-component version of the same single component
378 * partitioner used in the creation of the flattenedArray partitioner.
379 */
380 void
382 std::vector<std::complex<double>> &flattenedArray,
383 const dftfe::uInt totalNumberComponents,
384 const std::pair<dftfe::uInt, dftfe::uInt> componentIndexRange,
385 const std::vector<dealii::types::global_dof_index>
386 &localProcDofIndicesReal,
387 const std::vector<dealii::types::global_dof_index>
388 &localProcDofIndicesImag,
389 const std::vector<distributedCPUVec<double>> &componentVectors);
390
391#else
392 /** @brief Copies to a flattenedArray stl
393 * vector containing multiple component fields from a single field parallel
394 * distributed vector.
395 *
396 * @param[out] flattenedArray flattened stl vector with multiple component
397 * fields
398 * @param[in] totalNumberComponents total number of component fiels in
399 * flattenedArray
400 * @param[in] componentIndexRange desired range field components
401 * [componentIndexRange.first,componentIndexRange.second)
402 * @param[in] componentVectors vector of parallel distributed vectors with
403 * fields corresponding to componentIndexRange. componentVectors is expected
404 * to be of the size componentIndexRange.second-componentIndexRange.first.
405 * Further, each entry of componentVectors is assumed to be already
406 * initialized with the same single component partitioner used in the
407 * creation of the flattenedArray partitioner.
408 */
409 void
411 std::vector<double> &flattenedArray,
412 const dftfe::uInt totalNumberComponents,
413 const std::pair<dftfe::uInt, dftfe::uInt> componentIndexRange,
414 const std::vector<distributedCPUVec<double>> &componentVectors);
415
416#endif
417
418 std::pair<dealii::Point<3>, dealii::Point<3>>
419 createBoundingBoxTriaLocallyOwned(const dealii::DoFHandler<3> &dofHandler);
420
421
422 void
424 const dealii::MatrixFree<3, double> &matrix_free_data,
425 const dftfe::uInt mfDofHandlerIndex,
426 std::vector<dftfe::uInt> &nodesPerCellClassificationMap);
427
428
429 void
431 const dealii::MatrixFree<3, double> &matrix_free_data,
432 const dftfe::uInt mfDofHandlerIndex,
433 const dealii::AffineConstraints<double> &constraintMatrix,
434 std::vector<dftfe::uInt> &nodesPerCellClassificationMap,
435 std::vector<dftfe::uInt> &globalArrayClassificationMap);
436
437
438
439 } // namespace vectorTools
440} // namespace dftfe
441#endif
A class template to store the communication pattern (i.e., which entries/nodes to receive from which ...
Definition MPIPatternP2P.h:57
Contains generic utils functions related to custom partitioned flattened dealii vector.
Definition vectorUtilities.h:34
void classifyInteriorSurfaceNodesInCell(const dealii::MatrixFree< 3, double > &matrix_free_data, const dftfe::uInt mfDofHandlerIndex, std::vector< dftfe::uInt > &nodesPerCellClassificationMap)
void copyFlattenedSTLVecToSingleCompVec(const double *flattenedArray, const dftfe::uInt totalNumberComponents, const dftfe::uInt localVectorSize, const std::pair< dftfe::uInt, dftfe::uInt > componentIndexRange, std::vector< distributedCPUVec< double > > &componentVectors)
Copies a single field component from a flattenedArray STL vector containing multiple component fields...
void computeCellLocalIndexSetMap(const std::shared_ptr< const utils::mpi::MPIPatternP2P< dftfe::utils::MemorySpace::HOST > > &partitioner, const dealii::MatrixFree< 3, double > &matrix_free_data, const dftfe::uInt mfDofHandlerIndex, const dftfe::uInt blockSize, std::vector< std::vector< dftfe::uInt > > &flattenedArrayMacroCellLocalProcIndexId, std::vector< std::vector< dftfe::uInt > > &flattenedArrayCellLocalProcIndexId)
Creates a cell local index set map for flattened array.
void createDealiiVector(const std::shared_ptr< const dealii::Utilities::MPI::Partitioner > &partitioner, const dftfe::uInt blockSize, distributedCPUVec< T > &flattenedArray)
Creates a custom partitioned flattened dealii vector. stores multiple components asociated with a nod...
void makeAffineConstraintsConsistentInParallel(const dealii::DoFHandler< 3 > &dofHandlerPar, dealii::AffineConstraints< double > &constraints)
void createParallelConstraintMatrixFromSerial(const dealii::Triangulation< 3, 3 > &serTria, const dealii::DoFHandler< 3 > &dofHandlerPar, const MPI_Comm &mpi_comm_parent, const MPI_Comm &mpi_comm_domain, const std::vector< std::vector< double > > &domainBoundingVectors, dealii::AffineConstraints< double > &periodicHangingConstraints, dealii::AffineConstraints< double > &onlyHangingConstraints, const dftfe::Int verbosity, const bool periodicX, const bool periodicY, const bool periodicZ)
Create constraint matrix using serial mesh. Temporary fix for a bug (Issue #7053) in deal....
void copySingleCompVecToFlattenedDealiiVec(distributedCPUVec< double > &flattenedArray, const dftfe::uInt totalNumberComponents, const std::pair< dftfe::uInt, dftfe::uInt > componentIndexRange, const std::vector< distributedCPUVec< double > > &componentVectors)
Copies to a flattenedArray parallel distributed vector containing multiple component fields from a si...
std::pair< dealii::Point< 3 >, dealii::Point< 3 > > createBoundingBoxTriaLocallyOwned(const dealii::DoFHandler< 3 > &dofHandler)
void classifyInteriorSurfaceNodesInGlobalArray(const dealii::MatrixFree< 3, double > &matrix_free_data, const dftfe::uInt mfDofHandlerIndex, const dealii::AffineConstraints< double > &constraintMatrix, std::vector< dftfe::uInt > &nodesPerCellClassificationMap, std::vector< dftfe::uInt > &globalArrayClassificationMap)
void copySingleCompVecToFlattenedSTLVec(std::vector< double > &flattenedArray, const dftfe::uInt totalNumberComponents, const std::pair< dftfe::uInt, dftfe::uInt > componentIndexRange, const std::vector< distributedCPUVec< double > > &componentVectors)
Copies to a flattenedArray stl vector containing multiple component fields from a single field parall...
void copyFlattenedDealiiVecToSingleCompVec(const distributedCPUVec< double > &flattenedArray, const dftfe::uInt totalNumberComponents, const std::pair< dftfe::uInt, dftfe::uInt > componentIndexRange, std::vector< distributedCPUVec< double > > &componentVectors, const bool isFlattenedDealiiGhostValuesUpdated=false)
Copies a single field component from a flattenedArray parallel distributed vector containing multiple...
Definition pseudoPotentialToDftfeConverter.cc:34
dealii::LinearAlgebra::distributed::Vector< elem_type, dealii::MemorySpace::Host > distributedCPUVec
Definition headers.h:92
std::uint32_t uInt
Definition TypeConfig.h:10
std::int32_t Int
Definition TypeConfig.h:11