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