DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
constraintMatrixInfo.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (c) 2017-2025 The Regents of the University of Michigan and DFT-FE
4// authors.
5//
6// This file is part of the DFT-FE code.
7//
8// The DFT-FE code is free software; you can use it, redistribute
9// it, and/or modify it under the terms of the GNU Lesser General
10// Public License as published by the Free Software Foundation; either
11// version 2.1 of the License, or (at your option) any later version.
12// The full text of the license can be found in the file LICENSE at
13// the top level of the DFT-FE distribution.
14//
15// ---------------------------------------------------------------------
16//
17
18
19#ifndef constraintMatrixInfo_H_
20#define constraintMatrixInfo_H_
21
22#include <vector>
23#include <MemoryStorage.h>
24#include "headers.h"
25
26namespace dftfe
27{
28 //
29 // Declare dftUtils functions
30 //
31 namespace dftUtils
32 {
33 /**
34 * @brief Overloads dealii's distribute and distribute_local_to_global functions associated with constraints class.
35 * Stores the dealii's constraint matrix data into STL vectors for faster
36 * memory access costs
37 *
38 * @author Phani Motamarri
39 *
40 */
41 template <dftfe::utils::MemorySpace memorySpace>
43 {
44 public:
45 /**
46 * class constructor
47 */
49
50 /**
51 * class destructor
52 */
54
55 /**
56 * @brief convert a given constraintMatrix to simple arrays (STL) for fast access
57 *
58 * @param partitioner associated with the dealii vector
59 * @param constraintMatrixData dealii constraint matrix from which the data is extracted
60 */
61 void
63 const std::shared_ptr<const dealii::Utilities::MPI::Partitioner>
64 & partitioner,
65 const dealii::AffineConstraints<double> &constraintMatrixData);
66
67 /**
68 * @brief overloaded dealii internal function "distribute" which sets the slave node
69 * field values from master nodes
70 *
71 * @param fieldVector parallel dealii vector
72 */
73 void
75
76 /**
77 * @brief overloaded dealii internal function distribute for flattened dealii array which sets
78 * the slave node field values from master nodes
79 *
80 * @param blockSize number of components for a given node
81 */
82 template <typename T>
83 void
85 const unsigned int blockSize) const;
86
87 template <typename T>
88 void
91
92 /**
93 * @brief transfers the contributions of slave nodes to master nodes using the constraint equation
94 * slave nodes are the nodes which are to the right of the constraint
95 * equation and master nodes are the nodes which are left of the
96 * constraint equation.
97 *
98 * @param fieldVector parallel dealii vector which is the result of matrix-vector product(vmult) withot taking
99 * care of constraints
100 * @param blockSize number of components for a given node
101 */
102 template <typename T>
103 void
105 const unsigned int blockSize) const;
106
107 template <typename T>
108 void
111
112 /**
113 * @brief Scales the constraints with the inverse diagonal mass matrix so that the scaling of the vector can be done at the cell level
114 *
115 * @param invSqrtMassVec the inverse diagonal mass matrix
116 */
117 void
119 const distributedCPUVec<double> &invSqrtMassVec);
120
121 void
124
125
126 /**
127 * @brief sets field values at constrained nodes to be zero
128 *
129 * @param fieldVector parallel dealii vector with fields stored in a flattened format
130 * @param blockSize number of field components for a given node
131 */
132 template <typename T>
133 void
135 const unsigned int blockSize) const;
136 template <typename T>
137 void
140
141 /**
142 * clear data members
143 */
144 void
146
147
148 private:
149 std::vector<dealii::types::global_dof_index> d_rowIdsGlobal;
150 std::vector<dealii::types::global_dof_index> d_rowIdsLocal;
151 std::vector<dealii::types::global_dof_index> d_columnIdsLocal;
152 std::vector<dealii::types::global_dof_index> d_columnIdsGlobal;
153 std::vector<double> d_columnValues;
154 std::vector<double> d_inhomogenities;
155 std::vector<dealii::types::global_dof_index> d_rowSizes;
156 std::vector<dealii::types::global_dof_index>
158 };
159
160#if defined(DFTFE_WITH_DEVICE)
161 /**
162 * @brief Overloads dealii's distribute and distribute_local_to_global functions associated with constraints class.
163 * Stores the dealii's constraint matrix data into STL vectors for faster
164 * memory access costs
165 *
166 * @author Sambit Das, Phani Motamarri
167 *
168 */
169 template <>
171 {
172 public:
173 /**
174 * class constructor
175 */
177
178 /**
179 * class destructor
180 */
182
183 /**
184 * @brief convert a given constraintMatrix to simple arrays (STL) for fast access
185 *
186 * @param partitioner associated with the dealii vector
187 * @param constraintMatrixData dealii constraint matrix from which the data is extracted
188 */
189 void
191 const std::shared_ptr<const dealii::Utilities::MPI::Partitioner>
192 & partitioner,
193 const dealii::AffineConstraints<double> &constraintMatrixData,
194 const bool useInhomogeneties = true);
195
196 void
197 distribute(distributedCPUVec<double> &fieldVector) const;
198
199 template <typename T>
200 void
201 distribute(distributedCPUVec<T> &fieldVector,
202 const unsigned int blockSize) const;
203
204 /**
205 * @brief overloaded dealii internal function distribute for flattened dealii array which sets
206 * the slave node field values from master nodes
207 *
208 * @param blockSize number of components for a given node
209 */
210 template <typename NumberType>
211 void
215 &fieldVector) const;
216
217
218 /**
219 * @brief Scales the constraints with the inverse diagonal mass matrix so that the scaling of the vector can be done at the cell level
220 *
221 * @param invSqrtMassVec the inverse diagonal mass matrix
222 */
223 void
225 const dftfe::utils::MemoryStorage<double,
227 &invSqrtMassVec);
228
229
230 /**
231 * @brief transfers the contributions of slave nodes to master nodes using the constraint equation
232 * slave nodes are the nodes which are to the right of the constraint
233 * equation and master nodes are the nodes which are left of the
234 * constraint equation.
235 *
236 * @param fieldVector parallel dealii vector which is the result of matrix-vector product(vmult) withot taking
237 * care of constraints
238 * @param blockSize number of components for a given node
239 */
240 template <typename NumberType>
241 void
243 distributedDeviceVec<NumberType> &fieldVector) const;
244
245 template <typename T>
246 void
248 const unsigned int blockSize) const;
249
250 void
252 const distributedCPUVec<double> &invSqrtMassVec);
253
254 /**
255 * @brief sets field values at constrained nodes to be zero
256 *
257 * @param fieldVector parallel dealii vector with fields stored in a flattened format
258 * @param blockSize number of field components for a given node
259 */
260 template <typename NumberType>
261 void
262 set_zero(distributedDeviceVec<NumberType> &fieldVector) const;
263
264 template <typename T>
265 void
266 set_zero(distributedCPUVec<T> &fieldVector,
267 const unsigned int blockSize) const;
268
269 /**
270 * clear data members
271 */
272 void
273 clear();
274
275
276 private:
277 std::vector<unsigned int> d_rowIdsLocal;
278 std::vector<unsigned int> d_columnIdsLocal;
279 std::vector<double> d_columnValues;
280 std::vector<double> d_inhomogenities;
281 std::vector<unsigned int> d_rowSizes;
282 std::vector<unsigned int> d_rowSizesAccumulated;
283 std::vector<dealii::types::global_dof_index>
285
286 dftfe::utils::MemoryStorage<unsigned int,
288 d_rowIdsLocalDevice;
289 dftfe::utils::MemoryStorage<unsigned int,
291 d_columnIdsLocalDevice;
293 d_columnValuesDevice;
295 d_inhomogenitiesDevice;
296 dftfe::utils::MemoryStorage<unsigned int,
298 d_rowSizesDevice;
299 dftfe::utils::MemoryStorage<unsigned int,
301 d_rowSizesAccumulatedDevice;
302 dftfe::utils::MemoryStorage<dealii::types::global_dof_index,
304 d_localIndexMapUnflattenedToFlattenedDevice;
305
306 std::vector<unsigned int> d_rowIdsLocalBins;
307 std::vector<unsigned int> d_columnIdsLocalBins;
308 std::vector<unsigned int> d_columnIdToRowIdMapBins;
309 std::vector<double> d_columnValuesBins;
310 std::vector<unsigned int> d_binColumnSizes;
311 std::vector<unsigned int> d_binColumnSizesAccumulated;
312
313 dftfe::utils::MemoryStorage<unsigned int,
315 d_rowIdsLocalBinsDevice;
316 dftfe::utils::MemoryStorage<unsigned int,
318 d_columnIdsLocalBinsDevice;
319 dftfe::utils::MemoryStorage<unsigned int,
321 d_columnIdToRowIdMapBinsDevice;
323 d_columnValuesBinsDevice;
324
325 unsigned int d_numConstrainedDofs;
326 };
327#endif
328
329 } // namespace dftUtils
330
331} // namespace dftfe
332#endif
Overloads dealii's distribute and distribute_local_to_global functions associated with constraints cl...
Definition constraintMatrixInfo.h:43
std::vector< double > d_inhomogenities
Definition constraintMatrixInfo.h:154
std::vector< dealii::types::global_dof_index > d_rowIdsLocal
Definition constraintMatrixInfo.h:150
void distribute(dftfe::linearAlgebra::MultiVector< T, memorySpace > &fieldVector) const
void set_zero(distributedCPUVec< T > &fieldVector, const unsigned int blockSize) const
sets field values at constrained nodes to be zero
std::vector< double > d_columnValues
Definition constraintMatrixInfo.h:153
std::vector< dealii::types::global_dof_index > d_rowSizes
Definition constraintMatrixInfo.h:155
void distribute(distributedCPUVec< T > &fieldVector, const unsigned int blockSize) const
overloaded dealii internal function distribute for flattened dealii array which sets the slave node f...
std::vector< dealii::types::global_dof_index > d_rowIdsGlobal
Definition constraintMatrixInfo.h:149
void initialize(const std::shared_ptr< const dealii::Utilities::MPI::Partitioner > &partitioner, const dealii::AffineConstraints< double > &constraintMatrixData)
convert a given constraintMatrix to simple arrays (STL) for fast access
void set_zero(dftfe::linearAlgebra::MultiVector< T, memorySpace > &fieldVector) const
void initializeScaledConstraints(const distributedCPUVec< double > &invSqrtMassVec)
Scales the constraints with the inverse diagonal mass matrix so that the scaling of the vector can be...
std::vector< dealii::types::global_dof_index > d_columnIdsGlobal
Definition constraintMatrixInfo.h:152
std::vector< dealii::types::global_dof_index > d_localIndexMapUnflattenedToFlattened
Definition constraintMatrixInfo.h:157
void distribute_slave_to_master(dftfe::linearAlgebra::MultiVector< T, memorySpace > &fieldVector) const
void distribute(distributedCPUVec< double > &fieldVector) const
overloaded dealii internal function "distribute" which sets the slave node field values from master n...
void distribute_slave_to_master(distributedCPUVec< T > &fieldVector, const unsigned int blockSize) const
transfers the contributions of slave nodes to master nodes using the constraint equation slave nodes ...
std::vector< dealii::types::global_dof_index > d_columnIdsLocal
Definition constraintMatrixInfo.h:151
void initializeScaledConstraints(const dftfe::utils::MemoryStorage< double, memorySpace > &invSqrtMassVec)
An class template to encapsulate a MultiVector. A MultiVector is a collection of vectors belonging t...
Definition MultiVector.h:127
Definition MemoryStorage.h:33
Contains repeatedly used functions in the KSDFT calculations.
Definition CompositeData.h:29
@ DEVICE
Definition MemorySpaceType.h:36
Definition pseudoPotentialToDftfeConverter.cc:34
dealii::LinearAlgebra::distributed::Vector< elem_type, dealii::MemorySpace::Host > distributedCPUVec
Definition headers.h:92