jacobi_oo_STL
getmatrix.cpp
Go to the documentation of this file.
1 #include "getmatrix.h"
2 #include "userset.h"
3 
4 #include <algorithm>
5 #include <cassert>
6 #include <cmath>
7 #include <iomanip>
8 #include <iostream>
9 #include <list>
10 #include <vector>
11 using namespace std;
12 
13 
14 // general routine for lin. triangular elements
15 
16 void CalcElem(int const ial[3], double const xc[], double ske[3][3], double fe[3])
17 //void CalcElem(const int* __restrict__ ial, const double* __restrict__ xc, double* __restrict__ ske[3], double* __restrict__ fe)
18 {
19  const int i1 = 2 * ial[0], i2 = 2 * ial[1], i3 = 2 * ial[2];
20  const double x13 = xc[i3 + 0] - xc[i1 + 0], y13 = xc[i3 + 1] - xc[i1 + 1],
21  x21 = xc[i1 + 0] - xc[i2 + 0], y21 = xc[i1 + 1] - xc[i2 + 1],
22  x32 = xc[i2 + 0] - xc[i3 + 0], y32 = xc[i2 + 1] - xc[i3 + 1];
23  const double jac = fabs(x21 * y13 - x13 * y21);
24 
25  ske[0][0] = 0.5 / jac * (y32 * y32 + x32 * x32);
26  ske[0][1] = 0.5 / jac * (y13 * y32 + x13 * x32);
27  ske[0][2] = 0.5 / jac * (y21 * y32 + x21 * x32);
28  ske[1][0] = ske[0][1];
29  ske[1][1] = 0.5 / jac * (y13 * y13 + x13 * x13);
30  ske[1][2] = 0.5 / jac * (y21 * y13 + x21 * x13);
31  ske[2][0] = ske[0][2];
32  ske[2][1] = ske[1][2];
33  ske[2][2] = 0.5 / jac * (y21 * y21 + x21 * x21);
34 
35  const double xm = (xc[i1 + 0] + xc[i2 + 0] + xc[i3 + 0]) / 3.0,
36  ym = (xc[i1 + 1] + xc[i2 + 1] + xc[i3 + 1]) / 3.0;
37  //fe[0] = fe[1] = fe[2] = 0.5 * jac * FunctF(xm, ym) / 3.0;
38  fe[0] = fe[1] = fe[2] = 0.5 * jac * fNice(xm, ym) / 3.0;
39 }
40 
41 
42 // general routine for lin. triangular elements,
43 // non-symm. matrix
44 // node numbering in element: a s c e n d i n g indices !!
45 // GH: deprecated
46 void AddElem(int const ial[3], double const ske[3][3], double const fe[3],
47  int const id[], int const ik[], double sk[], double f[])
48 {
49  for (int i = 0; i < 3; ++i)
50  {
51  const int ii = ial[i], // row ii (global index)
52  id1 = id[ii], // start and
53  id2 = id[ii + 1]; // end of row ii in matrix
54  int ip = id1;
55  for (int j = 0; j < 3; ++j) // no symmetry assumed
56  {
57  const int jj = ial[j];
58  bool not_found = true;
59  do // find entry jj (global index) in row ii
60  {
61  not_found = (ik[ip] != jj);
62  ++ip;
63  }
64  while (not_found && ip < id2);
65 
66 #ifndef NDEBUG // compiler option -DNDEBUG switches off the check
67  if (not_found) // no entry found !!
68  {
69  cout << "Error in AddElem: (" << ii << "," << jj << ") ["
70  << ial[0] << "," << ial[1] << "," << ial[2] << "]\n";
71  assert(!not_found);
72  }
73 #endif
74  sk[ip - 1] += ske[i][j];
75  }
76  f[ii] += fe[i];
77  }
78 }
79 
80 
81 // ----------------------------------------------------------------------------
82 
83 
84 
85 
86 // ####################################################################
87 
89  : _mesh(mesh), _nrows(0), _nnz(0), _id(0), _ik(0), _sk(0)
90 {
92  return;
93 }
94 
96 {
97  int const nelem(_mesh.Nelems());
98  int const ndof_e(_mesh.NdofsElement());
99  auto const &ia(_mesh.GetConnectivity());
100 // Determine the number of matrix rows
101  _nrows = *max_element(ia.cbegin(), ia.cbegin() + ndof_e * nelem);
102  ++_nrows; // node numberng: 0 ... nnode-1
103  assert(*min_element(ia.cbegin(), ia.cbegin() + ndof_e * nelem) == 0); // numbering starts with 0 ?
104 
105 // Collect for each node those nodes it is connected to (multiple entries)
106 // Detect the neighboring nodes
107  vector< list<int> > cc(_nrows); // cc[i] is the list of nodes a node i is connected to
108  for (int i = 0; i < nelem; ++i)
109  {
110  int const idx = ndof_e * i;
111  for (int k = 0; k < ndof_e; ++k)
112  {
113  list<int> &cck = cc.at(ia[idx + k]);
114  cck.insert( cck.end(), ia.cbegin() + idx, ia.cbegin() + idx + ndof_e );
115  }
116  }
117 // Delete the multiple entries
118  _nnz = 0;
119  for (auto &it : cc)
120  {
121  it.sort();
122  it.unique();
123  _nnz += static_cast<int>(it.size());
124  // cout << it.size() << " :: "; copy(it->begin(),it->end(), ostream_iterator<int,char>(cout," ")); cout << endl;
125  }
126 
127 // CSR data allocation
128  _id.resize(_nrows + 1); // Allocate memory for CSR row pointer
129  _ik.resize(_nnz); // Allocate memory for CSR column index vector
130 
131 // copy CSR data
132  _id[0] = 0; // begin of first row
133  for (size_t i = 0; i < cc.size(); ++i)
134  {
135  //cout << i << " " << nid.at(i) << endl;;
136  const list<int> &ci = cc.at(i);
137  const auto nci = static_cast<int>(ci.size());
138  _id[i + 1] = _id[i] + nci; // begin of next line
139  copy(ci.begin(), ci.end(), _ik.begin() + _id[i] );
140  }
141 
142  assert(_nnz == _id[_nrows]);
143  _sk.resize(_nnz); // Allocate memory for CSR column index vector
144  return;
145 }
146 
147 
148 void CRS_Matrix::Debug() const
149 {
150 // ID points to first entry of row
151 // no symmetry assumed
152  cout << "\nMatrix (nnz = " << _id[_nrows] << ")\n";
153 
154  for (int row = 0; row < _nrows; ++row)
155  {
156  cout << "Row " << row << " : ";
157  int const id1 = _id[row];
158  int const id2 = _id[row + 1];
159  for (int j = id1; j < id2; ++j)
160  {
161  cout.setf(ios::right, ios::adjustfield);
162  cout << "[" << setw(2) << _ik[j] << "] " << setw(4) << _sk[j] << " ";
163  }
164  cout << endl;
165  }
166  return;
167 }
168 
169 void CRS_Matrix::CalculateLaplace(vector<double> &f)
170 {
171  assert(_mesh.NdofsElement() == 3); // only for triangular, linear elements
172  //cout << _nnz << " vs. " << _id[_nrows] << " " << _nrows<< endl;
173  assert(_nnz == _id[_nrows]);
174 
175  for (int k = 0; k < _nrows; ++k)
176  {
177  _sk[k] = 0.0;
178  }
179  for (int k = 0; k < _nrows; ++k)
180  {
181  f[k] = 0.0;
182  }
183 
184  double ske[3][3], fe[3];
185  // Loop over all elements
186  auto const nelem = _mesh.Nelems();
187  auto const &ia = _mesh.GetConnectivity();
188  auto const &xc = _mesh.GetCoords();
189 
190  for (int i = 0; i < nelem; ++i)
191  {
192  CalcElem(ia.data() + 3 * i, xc.data(), ske, fe);
193  //AddElem(ia.data()+3 * i, ske, fe, _id.data(), _ik.data(), _sk.data(), f.data()); // GH: deprecated
194  AddElem_3(ia.data() + 3 * i, ske, fe, f);
195  }
196 
197  //Debug();
198 
199  return;
200 }
201 
202 void CRS_Matrix::ApplyDirichletBC(std::vector<double> const &u, std::vector<double> &f)
203 {
204  double const PENALTY = 1e6;
205  auto const idx = _mesh.Index_DirichletNodes();
206  int const nidx = static_cast<int>(idx.size());
207 
208  for (int row = 0; row < nidx; ++row)
209  {
210  int const k = idx[row];
211  int const id1 = fetch(k, k); // Find diagonal entry of row
212  assert(id1 >= 0);
213  _sk[id1] += PENALTY; // matrix weighted scaling feasible
214  f[k] += PENALTY * u[k];
215  }
216 
217  return;
218 }
219 
220 void CRS_Matrix::GetDiag(vector<double> &d) const
221 {
222  assert( _nrows == static_cast<int>(d.size()) );
223 
224  for (int row = 0; row < _nrows; ++row)
225  {
226  const int ia = fetch(row, row); // Find diagonal entry of row
227  assert(ia >= 0);
228  d[row] = _sk[ia];
229  }
230  return;
231 }
232 
233 bool CRS_Matrix::Compare2Old(int nnode, int const id[], int const ik[], double const sk[]) const
234 {
235  bool bn = (nnode == _nrows); // number of rows
236  if (!bn)
237  {
238  cout << "######### Error: " << "number of rows" << endl;
239  }
240 
241  bool bz = (id[nnode] == _nnz); // number of non zero elements
242  if (!bz)
243  {
244  cout << "######### Error: " << "number of non zero elements" << endl;
245  }
246 
247  bool bd = equal(id, id + nnode + 1, _id.cbegin()); // row starts
248  if (!bd)
249  {
250  cout << "######### Error: " << "row starts" << endl;
251  }
252 
253  bool bk = equal(ik, ik + id[nnode], _ik.cbegin()); // column indices
254  if (!bk)
255  {
256  cout << "######### Error: " << "column indices" << endl;
257  }
258 
259  bool bv = equal(sk, sk + id[nnode], _sk.cbegin()); // values
260  if (!bv)
261  {
262  cout << "######### Error: " << "values" << endl;
263  }
264 
265  return bn && bz && bd && bk && bv;
266 }
267 
268 
269 void CRS_Matrix::Mult(vector<double> &w, vector<double> const &u) const
270 {
271  assert( _nrows == static_cast<int>(w.size()) );
272  assert( w.size() == u.size() );
273 
274  for (int row = 0; row < _nrows; ++row)
275  {
276  double wi = 0.0;
277  for (int ij = _id[row]; ij < _id[row + 1]; ++ij)
278  {
279  wi += _sk[ij] * u[ _ik[ij] ];
280  }
281  w[row] = wi;
282  }
283  return;
284 }
285 
286 void CRS_Matrix::Defect(vector<double> &w,
287  vector<double> const &f, vector<double> const &u) const
288 {
289  assert( _nrows == static_cast<int>(w.size()) );
290  assert( w.size() == u.size() && u.size() == f.size() );
291 
292  for (int row = 0; row < _nrows; ++row)
293  {
294  double wi = f[row];
295  for (int ij = _id[row]; ij < _id[row + 1]; ++ij)
296  {
297  wi -= _sk[ij] * u[ _ik[ij] ];
298  }
299  w[row] = wi;
300  }
301  return;
302 }
303 
304 int CRS_Matrix::fetch(int const row, int const col) const
305 {
306  int const id2 = _id[row + 1]; // end and
307  int ip = _id[row]; // start of recent row (global index)
308 
309  while (ip < id2 && _ik[ip] != col) // find index col (global index)
310  {
311  ++ip;
312  }
313  if (ip >= id2)
314  {
315  ip = -1;
316 #ifndef NDEBUG // compiler option -DNDEBUG switches off the check
317  cout << "No column " << col << " in row " << row << endl;
318  assert(ip >= id2);
319 #endif
320  }
321  return ip;
322 }
323 
324 
325 void CRS_Matrix::AddElem_3(int const ial[3], double const ske[3][3], double const fe[3], vector<double> &f)
326 {
327  for (int i = 0; i < 3; ++i)
328  {
329  const int ii = ial[i]; // row ii (global index)
330  for (int j = 0; j < 3; ++j) // no symmetry assumed
331  {
332  const int jj = ial[j]; // column jj (global index)
333  int ip = fetch(ii, jj); // find column entry jj in row ii
334 #ifndef NDEBUG // compiler option -DNDEBUG switches off the check
335  if (ip < 0) // no entry found !!
336  {
337  cout << "Error in AddElem: (" << ii << "," << jj << ") ["
338  << ial[0] << "," << ial[1] << "," << ial[2] << "]\n";
339  assert(ip >= 0);
340  }
341 #endif
342  _sk[ip] += ske[i][j];
343  }
344  f[ii] += fe[i];
345  }
346 }
347 
xc
xc
Definition: ascii_read_meshvector.m:30
fNice
double fNice(double const x, double const y)
Definition: userset.h:27
AddElem
void AddElem(int const ial[3], double const ske[3][3], double const fe[3], int const id[], int const ik[], double sk[], double f[])
Definition: getmatrix.cpp:46
Mesh::Index_DirichletNodes
virtual std::vector< int > Index_DirichletNodes() const =0
CRS_Matrix::Compare2Old
bool Compare2Old(int nnode, int const id[], int const ik[], double const sk[]) const
Definition: getmatrix.cpp:233
Mesh::GetConnectivity
const std::vector< int > & GetConnectivity() const
Definition: geom.h:93
CRS_Matrix::CRS_Matrix
CRS_Matrix(Mesh const &mesh)
Definition: getmatrix.cpp:88
Mesh
Definition: geom.h:11
userset.h
CRS_Matrix::ApplyDirichletBC
void ApplyDirichletBC(std::vector< double > const &u, std::vector< double > &f)
Definition: getmatrix.cpp:202
CRS_Matrix::Mult
void Mult(std::vector< double > &w, std::vector< double > const &u) const
Definition: getmatrix.cpp:269
CRS_Matrix::Debug
void Debug() const
Definition: getmatrix.cpp:148
CRS_Matrix::CalculateLaplace
void CalculateLaplace(std::vector< double > &f)
Definition: getmatrix.cpp:169
CRS_Matrix::GetDiag
void GetDiag(std::vector< double > &d) const
Definition: getmatrix.cpp:220
CalcElem
void CalcElem(int const ial[3], double const xc[], double ske[3][3], double fe[3])
Definition: getmatrix.cpp:16
CRS_Matrix::AddElem_3
void AddElem_3(int const ial[3], double const ske[3][3], double const fe[3], std::vector< double > &f)
Definition: getmatrix.cpp:325
ia
ia
Definition: ascii_read_meshvector.m:35
nnode
function vertex minimal boundary edge info in an ASCII file Matlab indexing is stored(starts with 1). % % The output file format is compatible with Mesh_2d_3_matlab nnode
Definition: ascii_write_mesh.m:23
getmatrix.h
CRS_Matrix::Defect
void Defect(std::vector< double > &w, std::vector< double > const &f, std::vector< double > const &u) const
Definition: getmatrix.cpp:286
CRS_Matrix::Derive_Matrix_Pattern
void Derive_Matrix_Pattern()
Definition: getmatrix.cpp:95
nelem
nelem
Definition: ascii_read_meshvector.m:24
Mesh::GetCoords
const std::vector< double > & GetCoords() const
Definition: geom.h:124
Mesh::Nelems
int Nelems() const
Definition: geom.h:35
CRS_Matrix::fetch
int fetch(int row, int col) const
Definition: getmatrix.cpp:304
Mesh::NdofsElement
int NdofsElement() const
Definition: geom.h:53