forked from dftfeDevelopers/dftfe
-
Notifications
You must be signed in to change notification settings - Fork 0
/
constraintMatrixInfo.cc
324 lines (260 loc) · 9.18 KB
/
constraintMatrixInfo.cc
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
// ---------------------------------------------------------------------
//
// Copyright (c) 2017-2018 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.
//
// ---------------------------------------------------------------------
//
// @author Phani Motamarri
//
#include <constraintMatrixInfo.h>
#include <linearAlgebraOperations.h>
namespace dftfe {
//
//Declare dftUtils functions
//
namespace dftUtils
{
//
//wrapper function to call blas function daxpy or zapxy depending
//on the data type (complex or double)
//
void callaxpy(const unsigned int *n,
const double *alpha,
double *x,
const unsigned int *incx,
double *y,
const unsigned int *incy)
{
daxpy_(n,
alpha,
x,
incx,
y,
incy);
}
void callaxpy(const unsigned int *n,
const std::complex<double> *alpha,
std::complex<double> *x,
const unsigned int *incx,
std::complex<double> *y,
const unsigned int *incy)
{
zaxpy_(n,
alpha,
x,
incx,
y,
incy);
}
//
//constructor
//
constraintMatrixInfo::constraintMatrixInfo()
{
}
//
//destructor
//
constraintMatrixInfo::~constraintMatrixInfo()
{
}
//
//store constraintMatrix row data in STL vector
//
void constraintMatrixInfo::initialize(const std::shared_ptr< const dealii::Utilities::MPI::Partitioner > & partitioner,
const dealii::ConstraintMatrix & constraintMatrixData)
{
clear();
const dealii::IndexSet & locally_owned_dofs = partitioner->locally_owned_range();
const dealii::IndexSet & ghost_dofs = partitioner->ghost_indices();
for(dealii::IndexSet::ElementIterator it = locally_owned_dofs.begin(); it != locally_owned_dofs.end();++it)
{
if(constraintMatrixData.is_constrained(*it))
{
const dealii::types::global_dof_index lineDof = *it;
d_rowIdsLocal.push_back(partitioner->global_to_local(lineDof));
d_rowIdsGlobal.push_back(lineDof);
d_inhomogenities.push_back(constraintMatrixData.get_inhomogeneity(lineDof));
const std::vector<std::pair<dealii::types::global_dof_index, double > > * rowData=constraintMatrixData.get_constraint_entries(lineDof);
d_rowSizes.push_back(rowData->size());
for(unsigned int j = 0; j < rowData->size();++j)
{
Assert((*rowData)[j].first<partitioner->size(),
dealii::ExcMessage("Index out of bounds"));
d_columnIdsGlobal.push_back((*rowData)[j].first);
d_columnIdsLocal.push_back(partitioner->global_to_local((*rowData)[j].first));
d_columnValues.push_back((*rowData)[j].second);
}
}
}
for(dealii::IndexSet::ElementIterator it = ghost_dofs.begin(); it != ghost_dofs.end();++it)
{
if(constraintMatrixData.is_constrained(*it))
{
const dealii::types::global_dof_index lineDof = *it;
d_rowIdsLocal.push_back(partitioner->global_to_local(lineDof));
d_rowIdsGlobal.push_back(lineDof);
d_inhomogenities.push_back(constraintMatrixData.get_inhomogeneity(lineDof));
const std::vector<std::pair<dealii::types::global_dof_index, double > > * rowData=constraintMatrixData.get_constraint_entries(lineDof);
d_rowSizes.push_back(rowData->size());
for(unsigned int j = 0; j < rowData->size();++j)
{
Assert((*rowData)[j].first<partitioner->size(),
dealii::ExcMessage("Index out of bounds"));
d_columnIdsGlobal.push_back((*rowData)[j].first);
d_columnIdsLocal.push_back(partitioner->global_to_local((*rowData)[j].first));
d_columnValues.push_back((*rowData)[j].second);
}
}
}
}
void constraintMatrixInfo::precomputeMaps(const std::shared_ptr< const dealii::Utilities::MPI::Partitioner> & unFlattenedPartitioner,
const std::shared_ptr< const dealii::Utilities::MPI::Partitioner> & flattenedPartitioner,
const unsigned int blockSize)
{
//
//Get required sizes
//
const unsigned int n_ghosts = unFlattenedPartitioner->n_ghost_indices();
const unsigned int localSize = unFlattenedPartitioner->local_size();
const unsigned int totalSize = n_ghosts + localSize;
d_localIndexMapUnflattenedToFlattened.clear();
d_localIndexMapUnflattenedToFlattened.resize(totalSize);
//
//fill the data array
//
for(unsigned int ilocalDof = 0; ilocalDof < totalSize; ++ilocalDof)
{
const dealii::types::global_dof_index globalIndex = unFlattenedPartitioner->local_to_global(ilocalDof);
d_localIndexMapUnflattenedToFlattened[ilocalDof] = flattenedPartitioner->global_to_local(globalIndex*blockSize);
}
}
//
//set the constrained degrees of freedom to values so that constraints
//are satisfied
//
void constraintMatrixInfo::distribute(dealii::parallel::distributed::Vector<double> &fieldVector) const
{
fieldVector.update_ghost_values();
unsigned int count = 0;
for(unsigned int i = 0; i < d_rowIdsLocal.size();++i)
{
double new_value = d_inhomogenities[i];
for(unsigned int j = 0; j < d_rowSizes[i]; ++j)
{
new_value += fieldVector.local_element(d_columnIdsLocal[count])*d_columnValues[count];
count++;
}
fieldVector.local_element(d_rowIdsLocal[i]) = new_value;
}
}
template<typename T>
void constraintMatrixInfo::distribute(dealii::parallel::distributed::Vector<T> &fieldVector,
const unsigned int blockSize) const
{
fieldVector.update_ghost_values();
unsigned int count = 0;
const unsigned int inc = 1;
std::vector<T> newValuesBlock(blockSize,0.0);
for(unsigned int i = 0; i < d_rowIdsLocal.size(); ++i)
{
std::fill(newValuesBlock.begin(),
newValuesBlock.end(),
d_inhomogenities[i]);
const dealii::types::global_dof_index startingLocalDofIndexRow = d_localIndexMapUnflattenedToFlattened[d_rowIdsLocal[i]];
for(unsigned int j = 0; j < d_rowSizes[i]; ++j)
{
Assert(count<d_columnIdsGlobal.size(),
dealii::ExcMessage("Overloaded distribute for flattened array has indices out of bounds"));
const dealii::types::global_dof_index startingLocalDofIndexColumn = d_localIndexMapUnflattenedToFlattened[d_columnIdsLocal[count]];
T alpha = d_columnValues[count];
callaxpy(&blockSize,
&alpha,
fieldVector.begin()+startingLocalDofIndexColumn,
&inc,
&newValuesBlock[0],
&inc);
count++;
}
std::copy(&newValuesBlock[0],
&newValuesBlock[0]+blockSize,
fieldVector.begin()+startingLocalDofIndexRow);
}
}
//
//set the constrained degrees of freedom to values so that constraints
//are satisfied for flattened array
//
template<typename T>
void constraintMatrixInfo::distribute_slave_to_master(dealii::parallel::distributed::Vector<T> & fieldVector,
const unsigned int blockSize) const
{
unsigned int count = 0;
const unsigned int inc = 1;
for(unsigned int i = 0; i < d_rowIdsLocal.size(); ++i)
{
const dealii::types::global_dof_index startingLocalDofIndexRow = d_localIndexMapUnflattenedToFlattened[d_rowIdsLocal[i]];
for(unsigned int j = 0; j < d_rowSizes[i]; ++j)
{
const dealii::types::global_dof_index startingLocalDofIndexColumn=d_localIndexMapUnflattenedToFlattened[d_columnIdsLocal[count]];
T alpha = d_columnValues[count];
callaxpy(&blockSize,
&alpha,
fieldVector.begin()+startingLocalDofIndexRow,
&inc,
fieldVector.begin()+startingLocalDofIndexColumn,
&inc);
count++;
}
//
//set slave contribution to zero
//
std::fill(fieldVector.begin()+startingLocalDofIndexRow,
fieldVector.begin()+startingLocalDofIndexRow+blockSize,
0.0);
}
}
template<typename T>
void constraintMatrixInfo::set_zero(dealii::parallel::distributed::Vector<T> & fieldVector,
const unsigned int blockSize) const
{
for(unsigned int i = 0; i < d_rowIdsLocal.size(); ++i)
{
const dealii::types::global_dof_index startingLocalDofIndexRow = d_localIndexMapUnflattenedToFlattened[d_rowIdsLocal[i]];
//set constrained nodes to zero
std::fill(fieldVector.begin()+startingLocalDofIndexRow,
fieldVector.begin()+startingLocalDofIndexRow+blockSize,
0.0);
}
}
//
//
//clear the data variables
//
void constraintMatrixInfo::clear()
{
d_rowIdsGlobal.clear();
d_rowIdsLocal.clear();
d_columnIdsLocal.clear();
d_columnIdsGlobal.clear();
d_columnValues.clear();
d_inhomogenities.clear();
d_rowSizes.clear();
}
template void constraintMatrixInfo::distribute(dealii::parallel::distributed::Vector<dataTypes::number> & fieldVector,
const unsigned int blockSize) const;
template void constraintMatrixInfo::distribute_slave_to_master(dealii::parallel::distributed::Vector<dataTypes::number> & fieldVector,
const unsigned int blockSize) const;
template void constraintMatrixInfo::set_zero(dealii::parallel::distributed::Vector<dataTypes::number> & fieldVector,
const unsigned int blockSize) const;
}
}