openCARP
Doxygen code documentation for the open cardiac electrophysiology simulator openCARP
SF_fem_utils.h
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 // openCARP is an open cardiac electrophysiology simulator.
3 //
4 // Copyright (C) 2020 openCARP project
5 //
6 // This program is licensed under the openCARP Academic Public License (APL)
7 // v1.0: You can use and redistribute it and/or modify it in non-commercial
8 // academic environments under the terms of APL as published by the openCARP
9 // project v1.0, or (at your option) any later version. Commercial use requires
10 // a commercial license (info@opencarp.org).
11 //
12 // This program is distributed without any warranty; see the openCARP APL for
13 // more details.
14 //
15 // You should have received a copy of the openCARP APL along with this program
16 // and can find it online: http://www.opencarp.org/license
17 // ----------------------------------------------------------------------------
18 
32 #ifndef _SF_FEM_H
33 #define _SF_FEM_H
34 
35 #include <cassert>
36 #include <cmath>
37 #include <cstring>
38 #include <mpi.h>
39 
40 #define SF_MAX_ELEM_NODES 10
41 
42 #include "SF_container.h"
43 
44 #include "SF_abstract_vector.h"
45 #include "SF_abstract_matrix.h"
46 #include "SF_abstract_lin_solver.h"
47 
48 namespace SF {
49 
58 inline short num_dof(elem_t type, short order)
59 {
60  short ndof = -1;
61  switch(type) {
62  case Line:
63  if (order == 1) ndof = 2;
64  else if (order == 2) ndof = 3;
65  else {
66  fprintf(stderr, "Line element order %d not implemented yet.\n", order);
67  exit(1);
68  }
69  break;
70  case Tri:
71  if (order == 1) ndof = 3;
72  else if (order == 2) ndof = 6;
73  else {
74  fprintf(stderr, "Tri element order %d not implemented yet.\n", order);
75  exit(1);
76  }
77  break;
78  case Quad:
79  if (order == 1) ndof = 4;
80  else if (order == 2) ndof = 8;
81  else {
82  fprintf(stderr, "Quad element order %d not implemented yet.\n", order);
83  exit(1);
84  }
85  break;
86  case Tetra:
87  if (order == 1) ndof = 4;
88  else if (order == 2) ndof = 10;
89  else {
90  fprintf(stderr, "Tetra element order %d not implemented yet.\n", order);
91  exit(1);
92  }
93  break;
94  case Hexa:
95  if (order == 1) ndof = 8;
96  else if (order == 2) ndof = 20; // serendipity element
97  else {
98  fprintf(stderr, "Hexa element order %d not implemented yet.\n", order);
99  exit(1);
100  }
101  break;
102 
103  default:
104  fprintf(stderr, "%s error: Unsupported element type.\n", __func__);
105  exit(1);
106  }
107 
108  return ndof;
109 }
110 
111 
121 inline void general_integration_points(const elem_t type,
122  const short order,
123  Point* ip,
124  double* w,
125  int & nint)
126 {
127  switch (type) {
128  case Line:
129  if (order == 1) {
130  ip[0].x = 0.; ip[0].y = 0.; ip[0].z = 0.;
131  w[0] = 1.; nint = 1;
132  } else if (order == 2) {
133  const double sqrt3 = 0.577350269189626; //=sqrt(1./3.)
134  ip[0].x = -sqrt3; ip[0].y = 0.; ip[0].z = 0.;
135  ip[1].x = +sqrt3; ip[1].y = 0.; ip[1].z = 0.;
136  w[0] = 1.; w[1] = 1.; nint = 2;
137  }
138  break;
139 
140  case Tri:
141  if (order == 1) {
142  ip[0].x = 1. / 3.; ip[0].y = 1. / 3.; ip[0].z = 0.;
143  w[0] = 1. / 2.;
144  nint = 1;
145  } else if (order == 2) {
146  ip[0].x = 1. / 6.; ip[0].y = 1. / 6.; ip[0].z = 0.;
147  ip[1].x = 4. / 6.; ip[1].y = 1. / 6.; ip[1].z = 0.;
148  ip[2].x = 1. / 6.; ip[2].y = 4. / 6.; ip[2].z = 0.;
149  w[0] = 1. / 6.; w[1] = 1. / 6.; w[2] = 1. / 6.;
150  nint = 3;
151  } else if (order == 3 || order == 4) {
152  // Gauss computed
153  ip[0].x = 0.188409405952072339650, ip[0].y = 0.787659461760847001700, ip[0].z = 0;
154  ip[1].x = 0.523979067720100721850, ip[1].y = 0.409466864440734712450, ip[1].z = 0;
155  ip[2].x = 0.808694385677669824730, ip[2].y = 0.088587959512703928766, ip[2].z = 0;
156  ip[3].x = 0.106170269119576471390, ip[3].y = 0.787659461760847001700, ip[3].z = 0;
157  ip[4].x = 0.295266567779632616020, ip[4].y = 0.409466864440734712450, ip[4].z = 0;
158  ip[5].x = 0.455706020243648035620, ip[5].y = 0.088587959512703928766, ip[5].z = 0;
159  ip[6].x = 0.023931132287080617016, ip[6].y = 0.787659461760847001700, ip[6].z = 0;
160  ip[7].x = 0.066554067839164496312, ip[7].y = 0.409466864440734712450, ip[7].z = 0;
161  ip[8].x = 0.102717654809626260380, ip[8].y = 0.088587959512703928766, ip[8].z = 0;
162 
163  w[0] = 0.019396383305959434551;
164  w[1] = 0.063678085099884929043;
165  w[2] = 0.055814420483044288601;
166  w[3] = 0.03103421328953510222;
167  w[4] = 0.1018849361598159059;
168  w[5] = 0.089303072772870889517;
169  w[6] = 0.019396383305959434551;
170  w[7] = 0.063678085099884929043;
171  w[8] = 0.055814420483044288601;
172  nint = 9;
173  }
174  break;
175 
176  case Quad:
177  if (order == 1) {
178  ip[0].x = 0.; ip[0].y = 0.; ip[0].z = 0.;
179  w[0] = 4.;
180  nint = 1;
181  } else if (order == 2 || order == 3) {
182  const double sqrt3 = 0.577350269189626; //=sqrt(1./3.)
183  ip[0].x = -sqrt3; ip[0].y = -sqrt3; ip[0].z = 0.;
184  ip[1].x = +sqrt3; ip[1].y = -sqrt3; ip[1].z = 0.;
185  ip[2].x = +sqrt3; ip[2].y = +sqrt3; ip[2].z = 0.;
186  ip[3].x = -sqrt3; ip[3].y = +sqrt3; ip[3].z = 0.;
187  w[0] = 1.; w[1] = 1.; w[2] = 1.; w[3] = 1.;
188  nint = 4;
189  } else if (order == 4) {
190  // Gauss computed
191  ip[0].x = 0.88729833462074170214, ip[0].y = 0.88729833462074170214, ip[0].z = 0;
192  ip[1].x = 0.88729833462074170214, ip[1].y = 0.5, ip[1].z = 0;
193  ip[2].x = 0.88729833462074170214, ip[2].y = 0.11270166537925829786, ip[2].z = 0;
194  ip[3].x = 0.5, ip[3].y = 0.88729833462074170214, ip[3].z = 0;
195  ip[4].x = 0.5, ip[4].y = 0.5, ip[4].z = 0;
196  ip[5].x = 0.5, ip[5].y = 0.11270166537925829786, ip[5].z = 0;
197  ip[6].x = 0.11270166537925829786, ip[6].y = 0.88729833462074170214, ip[6].z = 0;
198  ip[7].x = 0.11270166537925829786, ip[7].y = 0.5, ip[7].z = 0;
199  ip[8].x = 0.11270166537925829786, ip[8].y = 0.11270166537925829786, ip[8].z = 0;
200 
201  w[0] = 0.077160493827160225866;
202  w[1] = 0.12345679012345638081;
203  w[2] = 0.077160493827160225866;
204  w[3] = 0.12345679012345638081;
205  w[4] = 0.19753086419753024261;
206  w[5] = 0.12345679012345638081;
207  w[6] = 0.077160493827160225866;
208  w[7] = 0.12345679012345638081;
209  w[8] = 0.077160493827160225866;
210 
211  nint = 9;
212  }
213  break;
214 
215  case Tetra:
216  if (order == 0 || order == 1) {
217  // exact for linears
218  ip[0].x = 0.25; ip[0].y = 0.25; ip[0].z = 0.25;
219  w[0] = .16666666666666666666;
220  nint = 1;
221  } else if (order == 2) {
222  // exact for quadratics
223  double gauss1 = 0.13819660112501051518; // (5-sqrt(5))/20
224  double gauss2 = 0.58541019662496845446; // (5+3sqrt(5))/20
225  double weight = 0.0416666666666666666666666666666666; // 1/24
226  ip[0].x = gauss1; ip[0].y = gauss1; ip[0].z = gauss1;
227  ip[1].x = gauss2; ip[1].y = gauss1; ip[1].z = gauss1;
228  ip[2].x = gauss1; ip[2].y = gauss2; ip[2].z = gauss1;
229  ip[3].x = gauss1; ip[3].y = gauss1; ip[3].z = gauss2;
230  w[0] = weight; w[1] = weight; w[2] = weight; w[3] = weight;
231  nint = 4;
232  } else if (order == 3 || order == 4) {
233  // Gauss computed
234  ip[0].x = 0.12764656212038541505, ip[0].y = 0.29399880063162286969, ip[0].z = 0.54415184401122529412;
235  ip[1].x = 0.2457133252117133515, ip[1].y = 0.56593316507280100325, ip[1].z = 0.12251482265544139105;
236  ip[2].x = 0.30377276481470755209, ip[2].y = 0.070679724159396897787, ip[2].z = 0.54415184401122529412;
237  ip[3].x = 0.58474756320489440498, ip[3].y = 0.13605497680284600603, ip[3].z = 0.12251482265544139105;
238  ip[4].x = 0.034202793236766414198, ip[4].y = 0.29399880063162286969, ip[4].z = 0.54415184401122529412;
239  ip[5].x = 0.065838687060044420729, ip[5].y = 0.56593316507280100325, ip[5].z = 0.12251482265544139105;
240  ip[6].x = 0.081395667014670256001, ip[6].y = 0.070679724159396897787, ip[6].z = 0.54415184401122529412;
241  ip[7].x = 0.15668263733681833672, ip[7].y = 0.13605497680284600603, ip[7].z = 0.12251482265544139105;
242 
243  w[0] = 0.0091694299214797256314;
244  w[1] = 0.0211570064545240181520;
245  w[2] = 0.0160270405984766287080;
246  w[3] = 0.0369798563588529458080;
247  w[4] = 0.0091694299214797291009;
248  w[5] = 0.0211570064545240285600;
249  w[6] = 0.0160270405984766356470;
250  w[7] = 0.0369798563588529666250;
251 
252  nint = 8;
253  }
254  break;
255 
256  case Hexa:
257  if (order == 0) {
258  ip[0].x = 0.; ip[0].y = 0.; ip[0].z = 0.;
259  w[0] = 8.;
260  nint = 1;
261  } else if (order == 1 || order == 2) {
262  const double sqrt3 = 0.577350269189626; //=sqrt(1./3.)
263  ip[0].x = -sqrt3; ip[0].y = -sqrt3; ip[0].z = -sqrt3;
264  ip[1].x = +sqrt3; ip[1].y = -sqrt3; ip[1].z = -sqrt3;
265  ip[2].x = +sqrt3; ip[2].y = +sqrt3; ip[2].z = -sqrt3;
266  ip[3].x = -sqrt3; ip[3].y = +sqrt3; ip[3].z = -sqrt3;
267  ip[4].x = -sqrt3; ip[4].y = -sqrt3; ip[4].z = +sqrt3;
268  ip[5].x = +sqrt3; ip[5].y = -sqrt3; ip[5].z = +sqrt3;
269  ip[6].x = +sqrt3; ip[6].y = +sqrt3; ip[6].z = +sqrt3;
270  ip[7].x = -sqrt3; ip[7].y = +sqrt3; ip[7].z = +sqrt3;
271  w[0] = 1.; w[1] = 1.; w[2] = 1.; w[3] = 1.;
272  w[4] = 1.; w[5] = 1.; w[6] = 1.; w[7] = 1.;
273  nint = 8;
274  }
275  break;
276 
277  case Prism:
278  if (order == 0) {
279  ip[0].x = 0.33333333333333331, ip[0].y = 0.33333333333333337, ip[0].z = 0.5;
280  w[0] = 0.5;
281  nint = 1;
282  } else if(order == 1 || order == 2) {
283 #if 1
284  ip[0].x = 0.8168475628; ip[0].y = 0.0915762112; ip[0].z = 0.2113248706;
285  ip[1].x = 0.8168475628; ip[1].y = 0.0915762112; ip[1].z = 0.7886751294;
286  ip[2].x = 0.0915762112; ip[2].y = 0.8168475628; ip[2].z = 0.2113248706;
287  ip[3].x = 0.0915762112; ip[3].y = 0.8168475628; ip[3].z = 0.7886751294;
288  ip[4].x = 0.0915762112; ip[4].y = 0.0915762112; ip[4].z = 0.2113248706;
289  ip[5].x = 0.0915762112; ip[5].y = 0.0915762112; ip[5].z = 0.7886751294;
290  w[0] = 0.0833333358;
291  w[1] = 0.0833333358;
292  w[2] = 0.0833333358;
293  w[3] = 0.0833333358;
294  w[4] = 0.0833333358;
295  w[5] = 0.0833333358;
296  nint = 6;
297 #else
298  ip[0].x = 0.28001991549907407, ip[0].y = 0.64494897427831788, ip[0].z = 0.78867513459481287;
299  ip[1].x = 0.28001991549907407, ip[1].y = 0.64494897427831788, ip[1].z = 0.21132486540518713;
300  ip[2].x = 0.66639024601470143, ip[2].y = 0.15505102572168217, ip[2].z = 0.78867513459481287;
301  ip[3].x = 0.66639024601470143, ip[3].y = 0.15505102572168217, ip[3].z = 0.21132486540518713;
302  ip[4].x = 0.075031110222608124, ip[4].y = 0.64494897427831788, ip[4].z = 0.78867513459481287;
303  ip[5].x = 0.075031110222608124, ip[5].y = 0.64494897427831788, ip[5].z = 0.21132486540518713;
304  ip[6].x = 0.17855872826361643, ip[6].y = 0.15505102572168217, ip[6].z = 0.78867513459481287;
305  ip[7].x = 0.17855872826361643, ip[7].y = 0.15505102572168217, ip[7].z = 0.21132486540518713;
306  w[0] = 0.045489654564005534;
307  w[1] = 0.045489654564005548;
308  w[2] = 0.079510345435994237;
309  w[3] = 0.079510345435994265;
310  w[4] = 0.045489654564005548;
311  w[5] = 0.045489654564005562;
312  w[6] = 0.079510345435994265;
313  w[7] = 0.079510345435994292;
314  nint = 8;
315 #endif
316  }
317  break;
318 
319  case Pyramid:
320  if (order == 1) {
321  ip[0].x = -0.433013; ip[0].y = -0.433013; ip[0].z = 0.25;
322  ip[1].x = 0.433013; ip[1].y = -0.433013; ip[1].z = 0.25;
323  ip[2].x = 0.433013; ip[2].y = 0.433013; ip[2].z = 0.25;
324  ip[3].x = -0.433013; ip[3].y = 0.433013; ip[3].z = 0.25;
325  w[0] = 1. / 3;
326  w[1] = 1. / 3;
327  w[2] = 1. / 3;
328  w[3] = 1. / 3;
329  nint = 4;
330  } else if (order == 2) {
331  ip[0 ].x = 0.040086493940919059, ip[0 ].y = 0.040086493940919059, ip[0 ].z = 0.93056815579702623;
332  ip[1 ].x = 0.19053106107826956, ip[1 ].y = 0.19053106107826956, ip[1 ].z = 0.66999052179242813;
333  ip[2 ].x = 0.38681920811135617, ip[2 ].y = 0.38681920811135617, ip[2 ].z = 0.33000947820757187;
334  ip[3 ].x = 0.53726377524870672, ip[3 ].y = 0.53726377524870672, ip[3 ].z = 0.069431844202973714;
335  ip[4 ].x = 0.040086493940919059, ip[4 ].y = -0.040086493940919059, ip[4 ].z = 0.93056815579702623;
336  ip[5 ].x = 0.19053106107826956, ip[5 ].y = -0.19053106107826956, ip[5 ].z = 0.66999052179242813;
337  ip[6 ].x = 0.38681920811135617, ip[6 ].y = -0.38681920811135617, ip[6 ].z = 0.33000947820757187;
338  ip[7 ].x = 0.53726377524870672, ip[7 ].y = -0.53726377524870672, ip[7 ].z = 0.069431844202973714;
339  ip[8 ].x = -0.040086493940919059, ip[8 ].y = 0.040086493940919059, ip[8 ].z = 0.93056815579702623;
340  ip[9 ].x = -0.19053106107826956, ip[9 ].y = 0.19053106107826956, ip[9 ].z = 0.66999052179242813;
341  ip[10].x = -0.38681920811135617, ip[10].y = 0.38681920811135617, ip[10].z = 0.33000947820757187;
342  ip[11].x = -0.53726377524870672, ip[11].y = 0.53726377524870672, ip[11].z = 0.069431844202973714;
343  ip[12].x = -0.040086493940919059, ip[12].y = -0.040086493940919059, ip[12].z = 0.93056815579702623;
344  ip[13].x = -0.19053106107826956, ip[13].y = -0.19053106107826956, ip[13].z = 0.66999052179242813;
345  ip[14].x = -0.38681920811135617, ip[14].y = -0.38681920811135617, ip[14].z = 0.33000947820757187;
346  ip[15].x = -0.53726377524870672, ip[15].y = -0.53726377524870672, ip[15].z = 0.069431844202973714;
347 
348  w[0 ] = 0.000838466012258970;
349  w[1 ] = 0.035511343496716564;
350  w[2 ] = 0.14636983865620462;
351  w[3 ] = 0.15061368516815288;
352  w[4 ] = 0.000838466012258971;
353  w[5 ] = 0.035511343496716585;
354  w[6 ] = 0.1463698386562047;
355  w[7 ] = 0.15061368516815296;
356  w[8 ] = 0.000838466012258971;
357  w[9 ] = 0.035511343496716585;
358  w[10] = 0.1463698386562047;
359  w[11] = 0.15061368516815296;
360  w[12] = 0.000838466012258971;
361  w[13] = 0.035511343496716599;
362  w[14] = 0.14636983865620476;
363  w[15] = 0.15061368516815302;
364  nint = 16;
365  }
366  break;
367 
368  default:
369  fprintf(stderr, "%s error: Unsupported element type.\n", __func__);
370  exit(1);
371  }
372 }
373 
374 
383 inline void reference_shape(const elem_t type,
384  const Point ip,
385  dmat<double> & rshape)
386 {
387  switch (type) {
388  case Line:
389  {
390  rshape[0][0] = 1.0 - ip.x;
391  rshape[0][1] = ip.x;
392 
393  rshape[1][0] = -1.0;
394  rshape[1][1] = 1.0;
395 
396  rshape[2][0] = 0.;
397  rshape[2][1] = 0.;
398 
399  rshape[3][0] = 0.;
400  rshape[3][1] = 0.;
401  break;
402  }
403 
404  case Tri:
405  {
406  double lam0 = (1.0 - ip.x - ip.y);
407  double lam1 = ip.x;
408  double lam2 = ip.y;
409 
410  rshape[0][0] = lam0;
411  rshape[0][1] = lam1;
412  rshape[0][2] = lam2;
413 
414  // derivative w.r.t. xi
415  rshape[1][0] = -1.0;
416  rshape[1][1] = 1.0;
417  rshape[1][2] = 0.0;
418 
419  // derivative w.r.t. eta
420  rshape[2][0] = -1.0;
421  rshape[2][1] = 0.0;
422  rshape[2][2] = 1.0;
423 
424  // derivative w.r.t. zeta (is zero)
425  rshape[3][0] = 0.;
426  rshape[3][1] = 0.;
427  rshape[3][2] = 0.;
428  break;
429  }
430 
431  case Quad: {
432  double qrtr = 0.25;
433  const double node[4][2] =
434  {
435  { -1.0, -1.0},
436  { 1.0, -1.0},
437  { 1.0, 1.0},
438  { -1.0, 1.0}
439  };
440 
441  // adjust to ordering given in carp manual
442  int v[4] = {0, 1, 2, 3};
443 
444  for (int i = 0; i < 4; i++) {
445  // shape function
446  rshape[0][i] = qrtr * (1. + ip.x * node[v[i]][0]) * (1. + ip.y * node[v[i]][1]);
447  // derivative w.r.t. xi
448  rshape[1][i] = node[v[i]][0] * qrtr * (1. + ip.y * node[v[i]][1]);
449  // derivative w.r.t. eta
450  rshape[2][i] = node[v[i]][1] * qrtr * (1. + ip.x * node[v[i]][0]);
451  // derivative w.r.t. zeta (is zero)
452  rshape[3][i] = 0.;
453  }
454  break;
455  }
456 
457  case Tetra:
458  {
459  double lam0 = 1.0 - ip.x - ip.y - ip.z;
460  double lam1 = ip.x;
461  double lam2 = ip.y;
462  double lam3 = ip.z;
463  //static int v[10] = {0,3,1,2,7,8,4,6,9,5};
464 
465  // shape function
466  rshape[0][0] = lam0;
467  rshape[0][1] = lam1;
468  rshape[0][2] = lam2;
469  rshape[0][3] = lam3;
470 
471  // derivative w.r.t. xi
472  rshape[1][0] = -1.0;
473  rshape[1][1] = 1.0;
474  rshape[1][2] = 0.0;
475  rshape[1][3] = 0.0;
476 
477  // derivative w.r.t. eta
478  rshape[2][0] = -1.0;
479  rshape[2][1] = 0.0;
480  rshape[2][2] = 1.0;
481  rshape[2][3] = 0.0;
482 
483  // derivative w.r.t. zeta
484  rshape[3][0] = -1.0;
485  rshape[3][1] = 0.0;
486  rshape[3][2] = 0.0;
487  rshape[3][3] = 1.0;
488  break;
489  }
490 
491  case Hexa:
492  {
493  const double oito = 1.0 / 8.0;
494  static const double node[8][3] =
495  {
496  { -1.0, -1.0, -1.0},
497  { 1.0, -1.0, -1.0},
498  { 1.0, 1.0, -1.0},
499  { -1.0, 1.0, -1.0},
500 
501  { -1.0, -1.0, 1.0},
502  { 1.0, -1.0, 1.0},
503  { 1.0, 1.0, 1.0},
504  { -1.0, 1.0, 1.0}
505  };
506 
507  // adjust to ordering given in carp manual
508  static int v[8] = {4, 7, 6, 5, 0, 1, 2, 3};
509 
510  for (int i = 0; i < 8; i++) {
511  // shape function
512  rshape[0][i] = oito * (1. + ip.x * node[v[i]][0]) * (1. + ip.y * node[v[i]][1]) * (1. + ip.z * node[v[i]][2]);
513  // derivative w.r.t. xi
514  rshape[1][i] = node[v[i]][0] * oito * (1. + ip.y * node[v[i]][1]) * (1. + ip.z * node[v[i]][2]);
515  // derivative w.r.t. eta
516  rshape[2][i] = node[v[i]][1] * oito * (1. + ip.x * node[v[i]][0]) * (1. + ip.z * node[v[i]][2]);
517  // derivative w.r.t zeta
518  rshape[3][i] = node[v[i]][2] * oito * (1. + ip.x * node[v[i]][0]) * (1. + ip.y * node[v[i]][1]);
519  }
520 
521  break;
522  }
523 
524  case Prism:
525  {
526  rshape[0][0] = (1.0 - ip.x - ip.y) * ip.z;
527  rshape[0][1] = ip.y * ip.z;
528  rshape[0][2] = ip.x * ip.z;
529  rshape[0][3] = (1.0 - ip.x - ip.y) * (1. - ip.z);
530  rshape[0][4] = ip.x * (1. - ip.z);
531  rshape[0][5] = ip.y * (1. - ip.z);
532 
533  // derivative w.r.t. x
534  rshape[1][0] = -ip.z;
535  rshape[1][1] = 0.0;
536  rshape[1][2] = ip.z;
537  rshape[1][3] = ip.z - 1.0;
538  rshape[1][4] = 1. - ip.z;
539  rshape[1][5] = 0.0;
540 
541  // derivative w.r.t. y
542  rshape[2][0] = -ip.z;
543  rshape[2][1] = ip.z;
544  rshape[2][2] = 0.0;
545  rshape[2][3] = ip.z - 1.0;
546  rshape[2][4] = 0.0;
547  rshape[2][5] = 1. - ip.z;
548 
549  // derivative w.r.t. z
550  rshape[3][0] = 1.0 - ip.x - ip.y;
551  rshape[3][1] = ip.y;
552  rshape[3][2] = ip.x;
553  rshape[3][3] = ip.x + ip.y - 1.0;
554  rshape[3][4] = -ip.x;
555  rshape[3][5] = -ip.y;
556  break;
557  }
558 
559  case Pyramid:
560  {
561  const double qrtr = 0.25;
562  // shape function
563  const double lterm0 = ip.x * ip.y * ip.z / (1.0 - ip.z);
564  rshape[0][0] = qrtr * ( (1.0 + ip.x) * (1.0 + ip.y) - ip.z + lterm0 );
565  rshape[0][1] = qrtr * ( (1.0 - ip.x) * (1.0 + ip.y) - ip.z - lterm0 );
566  rshape[0][2] = qrtr * ( (1.0 - ip.x) * (1.0 - ip.y) - ip.z + lterm0 );
567  rshape[0][3] = qrtr * ( (1.0 + ip.x) * (1.0 - ip.y) - ip.z - lterm0 );
568  rshape[0][4] = ip.z;
569 
570  // derivative w.r.t. xi
571  const double lterm1 = (ip.y * ip.z) / (1.0 - ip.z);
572  rshape[1][0] = qrtr * ( (1.0 + ip.y) + lterm1 );
573  rshape[1][1] = qrtr * ( -(1.0 + ip.y) - lterm1 );
574  rshape[1][2] = qrtr * ( -(1.0 - ip.y) + lterm1 );
575  rshape[1][3] = qrtr * ( (1.0 - ip.y) - lterm1 );
576  rshape[1][4] = 0.0;
577 
578  // derivative w.r.t. eta
579  const double lterm2 = (ip.x * ip.z) / (1.0 - ip.z);
580  rshape[2][0] = qrtr * ( (1.0 + ip.x) + lterm2 );
581  rshape[2][1] = qrtr * ( (1.0 - ip.x) - lterm2 );
582  rshape[2][2] = qrtr * ( -(1.0 - ip.x) + lterm2 );
583  rshape[2][3] = qrtr * ( -(1.0 + ip.x) - lterm2 );
584  rshape[2][4] = 0.0;
585 
586  // derivative w.r.t zeta
587  const double lterm3 = ((ip.x * ip.y * ip.z) / (1.0 - ip.z) * (1.0 - ip.z)) + (ip.x * ip.y / (1.0 - ip.z));
588  rshape[3][0] = qrtr * ( -1.0 + lterm3 );
589  rshape[3][1] = qrtr * ( -1.0 - lterm3 );
590  rshape[3][2] = qrtr * ( -1.0 + lterm3 );
591  rshape[3][3] = qrtr * ( -1.0 - lterm3 );
592  rshape[3][4] = 1.0;
593  break;
594  }
595 
596  default:
597  fprintf(stderr, "%s: Unimplemented element type! Aborting!\n", __func__);
598  exit(1);
599  }
600 }
601 
610 inline void jacobian_matrix(const dmat<double> & rshape,
611  const int npts,
612  const Point* pts,
613  double *J)
614 {
615  // zero the 3x3 jacobian matrix
616  memset (J, 0, 9 * sizeof(double) );
617  // note that ref_shape hold shape functions in [0][i] and derivatives in [k][i]
618  // if element is 2 dimensional then [3][i] is set to 0!
619  for (int i = 0; i < npts; i++)
620  {
621  J[0] += rshape[1][i] * pts[i].x;
622  J[1] += rshape[1][i] * pts[i].y;
623  J[2] += rshape[1][i] * pts[i].z;
624  J[3] += rshape[2][i] * pts[i].x;
625  J[4] += rshape[2][i] * pts[i].y;
626  J[5] += rshape[2][i] * pts[i].z;
627  J[6] += rshape[3][i] * pts[i].x;
628  J[7] += rshape[3][i] * pts[i].y;
629  J[8] += rshape[3][i] * pts[i].z;
630  }
631 }
632 
633 inline void invert_jacobian_matrix(const elem_t type, double* J, double & detJ)
634 {
635  switch(type)
636  {
637  // all 3D elems
638  default:
639  case Tetra:
640  invert_3x3(J, detJ);
641  break;
642 
643  // 2D elems
644  case Quad:
645  case Tri: {
646  double J2[4] = {J[0], J[1], J[3], J[4]};
647 
648  invert_2x2(J2, detJ);
649  detJ = fabs(detJ);
650 
651  J[0] = J2[0]; J[1] = J2[1]; J[2] = 0.0;
652  J[3] = J2[2]; J[4] = J2[3]; J[5] = 0.0;
653  J[6] = 0.0; J[7] = 0.0; J[8] = 0.0;
654  break;
655  }
656 
657  // 1D elem
658  case Line:
659  detJ = J[0];
660  J[0] = 1.0 / detJ;
661  break;
662  }
663 }
664 
674 inline void shape_deriv(const double *iJ,
675  const dmat<double> & rshape,
676  const int ndof,
677  dmat<double> & shape)
678 {
679  for (int in = 0; in < ndof; in++)
680  {
681  shape[1][in] = iJ[0] * rshape[1][in] +
682  iJ[1] * rshape[2][in] +
683  iJ[2] * rshape[3][in];
684 
685  shape[2][in] = iJ[3] * rshape[1][in] +
686  iJ[4] * rshape[2][in] +
687  iJ[5] * rshape[3][in];
688 
689  shape[3][in] = iJ[6] * rshape[1][in] +
690  iJ[7] * rshape[2][in] +
691  iJ[8] * rshape[3][in];
692  }
693 }
694 
695 
703 template<class T, class S>
705 {
706  private:
707  const meshdata<T, S> & _mesh;
708 // meshdata<T, S> & _mesh; ///< mesh reference, ACHTUNG: HAB HIER DEN CNST entferent, um die coords updaten zu können
709  const vector<T> & _glob_numbr;
710  T _esize;
711  const T* _offset_con;
712  size_t _eidx;
713  int _rank;
714 
715  public:
721  element_view(const meshdata<T, S> & mesh, const SF_nbr nbr) :
722  _mesh(mesh), _glob_numbr(_mesh.get_numbering(nbr))
723  {
724  if(mesh.l_numelem) this->set_elem(0);
725  MPI_Comm_rank(_mesh.comm, &_rank);
726  }
727 
733  inline void set_elem(size_t eidx)
734  {
735  _eidx = eidx;
736 
737  T offset = _mesh.dsp[_eidx];
738  _esize = _mesh.dsp[_eidx+1] - offset;
739  _offset_con = _mesh.con.data() + offset;
740  }
741 
747  inline bool next()
748  {
749  if(_eidx < (_mesh.l_numelem - 1) )
750  {
751  this->set_elem(_eidx + 1);
752  return true;
753  }
754  else
755  return false;
756  }
757 
763  inline T num_nodes() const
764  {
765  return _esize;
766  }
767 
773  inline elem_t type() const
774  {
775  return _mesh.type[_eidx];
776  }
777 
783  inline T tag() const
784  {
785  return _mesh.tag[_eidx];
786  }
787 
795  inline const T & node(short nidx) const
796  {
797  return _offset_con[nidx];
798  }
799 
807  inline const T & global_node(short nidx) const
808  {
809  return _glob_numbr[_offset_con[nidx]];
810  }
811 
820  inline const T & global_node(short nidx, SF_nbr nbr) const
821  {
822  const vector<T> & n = _mesh.get_numbering(nbr);
823  return n[_offset_con[nidx]];
824  }
825 
831  inline const T* nodes() const
832  {
833  return _offset_con;
834  }
842  inline Point coord(short nidx) const
843  {
844  T idx = _offset_con[nidx];
845  return {_mesh.xyz[idx*3+0], _mesh.xyz[idx*3+1], _mesh.xyz[idx*3+2]};
846  }
854  inline Point coord_upd(short nidx) const
855  {
856  T idx = _offset_con[nidx];
857  return {_mesh.xyz_Euler[idx*3+0], _mesh.xyz_Euler[idx*3+1], _mesh.xyz_Euler[idx*3+2]};
858  }
859 
861  inline void upd_coords(short nidx, double * pos_upd)
862  {
863  T idx = _offset_con[nidx];
864  _mesh.xyz_upd[3*idx + 0] = pos_upd[0];
865  _mesh.xyz_upd[3*idx + 1] = pos_upd[1];
866  _mesh.xyz_upd[3*idx + 2] = pos_upd[2];
867  }
868 
874  inline Point fiber() const
875  {
876  return {_mesh.fib[_eidx*3+0], _mesh.fib[_eidx*3+1], _mesh.fib[_eidx*3+2]};
877  }
883  inline Point sheet() const
884  {
885  if(_mesh.she.size())
886  return {_mesh.she[_eidx*3+0], _mesh.she[_eidx*3+1], _mesh.she[_eidx*3+2]};
887  else
888  return {0,0,0};
889  }
890 
896  bool has_sheet() const
897  {
898  return _mesh.she.size() > 0;
899  }
900 
906  inline size_t element_index() const
907  {
908  return _eidx;
909  }
910 
916  inline size_t global_element_index() const
917  {
918  return _mesh.epl.algebraic_layout()[_rank] + _eidx;
919  }
926  inline size_t global_element_index(SF_nbr nbr) const
927  {
928  const vector<T> & en = _mesh.get_numbering(nbr);
929  return en[_eidx];
930  }
931 
932  inline short num_dof(short order) const
933  {
934  return SF::num_dof(_mesh.type[_eidx], order);
935  }
936 
937  inline void integration_points(const short order, Point* ip, double* w, int & nint) const
938  {
939  general_integration_points(_mesh.type[_eidx], order, ip, w, nint);
940  }
941 
942  inline short dimension() const
943  {
944  switch(_mesh.type[_eidx])
945  {
946  default:
947  case Tetra:
948  return 3;
949 
950  case Tri:
951  case Quad:
952  return 2;
953 
954  case Line:
955  return 1;
956  }
957  }
958 };
959 
966 template<class T, class S>
968 {
969  public:
971  virtual void operator() (const element_view<T,S> & elem, dmat<double> & buff) = 0;
973 // virtual void operator() (const element_view<T,S> & elem, const element_view<T,S> & elem_upd, dmat<double> & buff) = 0; //needed to calculate deformation gradient
975  virtual void dpn(T & row_dpn, T & col_dpn) = 0;
976 };
977 
984 template<class T, class S>
986 {
987  protected:
988  void zero_buff(double* buff, T nrows)
989  {
990  for(T i=0; i<nrows; i++) buff[i] = 0.0;
991  }
992 
993  public:
995  virtual void operator() (const element_view<T,S> & elem, double* buff) = 0;
997  virtual void dpn(T & dpn) = 0;
998 };
999 
1000 
1010 template<class T, class V>
1011 inline void canonic_indices(const T* nidx,
1012  const T* nbr,
1013  const T esize,
1014  const short dpn,
1015  V* cidx)
1016 {
1017  for(T i=0; i<esize; i++)
1018  for(short j=0; j<dpn; j++)
1019  cidx[i*dpn + j] = nbr[nidx[i]]*dpn + j;
1020 }
1021 
1032 template<class T, class S>
1036 {
1037  int rank;
1038  MPI_Comm_rank(domain.comm, &rank);
1039 
1040  // we want to make sure that the element integrator fits the matrix
1041  T row_dpn, col_dpn;
1042  integrator.dpn(row_dpn, col_dpn);
1043  assert(mat.dpn_row() == row_dpn && mat.dpn_col() == col_dpn);
1044 
1045  // allocate row / col index buffers
1046  vector<SF_int> row_idx(SF_MAX_ELEM_NODES * row_dpn), col_idx(SF_MAX_ELEM_NODES * col_dpn);
1047 
1048  // allocate elem index buffer
1049  dmat<SF_real> ebuff(SF_MAX_ELEM_NODES * row_dpn, SF_MAX_ELEM_NODES * col_dpn);
1050 
1051  const vector<SF_int> & petsc_nbr = domain.get_numbering(NBR_PETSC);
1052 
1053  // start with assembly
1055  for(size_t eidx=0; eidx < domain.l_numelem; eidx++)
1056  {
1057  // set element view to current element
1058  view.set_elem(eidx);
1059 
1060  // calculate row/col indices of entries
1061  SF_int nnodes = view.num_nodes();
1062 
1063  row_idx.resize(nnodes*row_dpn);
1064  col_idx.resize(nnodes*col_dpn);
1065  canonic_indices(view.nodes(), petsc_nbr.data(), nnodes, row_dpn, row_idx.data());
1066  canonic_indices(view.nodes(), petsc_nbr.data(), nnodes, col_dpn, col_idx.data());
1067 
1068  // call integrator
1069  integrator(view, ebuff);
1070 
1071  // add values into system matrix
1072  const bool add = true;
1073 
1074 
1075  mat.set_values(row_idx, col_idx, ebuff.data(), add);
1076 
1077  }
1078  // finish assembly and progress output
1079  mat.finish_assembly();
1080 }
1081 
1082 template<class T, class S>
1086 {
1087  int rank;
1088  MPI_Comm_rank(domain.comm, &rank);
1089 
1090  // we want to make sure that the element integrator fits the matrix
1091  T row_dpn, col_dpn;
1092  integrator.dpn(row_dpn, col_dpn);
1093 
1094  assert(row_dpn == 1 && row_dpn == mat.dpn_row());
1095 
1096  // allocate row / col index buffers
1097  vector<SF_int> row_idx(SF_MAX_ELEM_NODES * row_dpn);
1098 
1099  // allocate elem index buffer
1100  dmat<S> ebuff(SF_MAX_ELEM_NODES * row_dpn, SF_MAX_ELEM_NODES * col_dpn);
1101  const vector<SF_int> & petsc_nbr = domain.get_numbering(NBR_PETSC);
1102 
1103  // start with assembly
1105  for(size_t eidx=0; eidx < domain.l_numelem; eidx++)
1106  {
1107  // set element view to current element
1108  view.set_elem(eidx);
1109 
1110  // calculate row/col indices of entries
1111  SF_int nnodes = view.num_nodes();
1112  row_idx.resize(nnodes*row_dpn);
1113  canonic_indices(view.nodes(), petsc_nbr.data(), nnodes, row_dpn, row_idx.data());
1114 
1115  // call integrator
1116  integrator(view, ebuff);
1117 
1118  // do lumping
1119  for(int i=0; i<nnodes; i++) {
1120  S rowsum = 0.0;
1121 
1122  for(int j=0; j<nnodes; j++)
1123  rowsum += ebuff[i][j];
1124 
1125  // add values into system matrix
1126  mat.set_value(row_idx[i], row_idx[i], rowsum, true);
1127  }
1128  }
1129 
1130  // finish assembly and progress output
1131  mat.finish_assembly();
1132 }
1133 
1145 template<class T, class S>
1149 {
1150  int rank;
1151  MPI_Comm_rank(vec.mesh->comm, &rank);
1152 
1153  // we want to make sure that the element integrator fits the matrix
1154  T dpn;
1155  integrator.dpn(dpn);
1156  assert(vec.dpn == dpn);
1157 
1158  // allocate row / col index buffers
1159  vector<SF_int> idx(SF_MAX_ELEM_NODES * dpn);
1160 
1161  // allocate elem index buffer
1162  vector<S> ebuff(SF_MAX_ELEM_NODES * dpn);
1163  const vector<SF_int> & petsc_nbr = domain.get_numbering(NBR_PETSC);
1164 
1165  // start with assembly
1167  for(size_t eidx=0; eidx < domain.l_numelem; eidx++)
1168  {
1169  // set element view to current element
1170  view.set_elem(eidx);
1171 
1172  // calculate row/col indices of entries
1173  SF_int nnodes = view.num_nodes();
1174  canonic_indices(view.nodes(), petsc_nbr.data(), nnodes, dpn, idx.data());
1175 
1176  // call integrator
1177  integrator(view, ebuff.data());
1178 
1179  // add values into system matrix
1180  vec.set(idx, ebuff, ADD_VALUES);
1181  }
1182  // finish assembly and progress output
1183  vec.finish_assembly();
1184 }
1185 
1186 template<class T, class S> inline
1188 {
1189  // the dpn has to be set for this to work
1190  int dpn = vec.dpn;
1192  assert(dpn > 0);
1193  assert(vec.layout == nodaltype);
1194 
1195  for(int i=0; i<view.num_nodes(); i++)
1196  for(int j=0; j<dpn; j++) {
1197  int idx = view.node(i)*dpn+j;
1198  buffer[i*dpn+j] = vec.get(idx);
1199  }
1200 }
1201 
1202 template<class T, class S> inline
1204 {
1205  // the dpn has to be set for this to work
1206  int dpn = vec.dpn;
1208  assert(dpn > 0);
1209  assert(vec.layout == nodaltype);
1210 
1211  SF_real* pvec = vec.ptr();
1212 
1213  for(int i=0; i<view.num_nodes(); i++)
1214  for(int j=0; j<dpn; j++) {
1215  int idx = view.node(i)*dpn+j;
1216  pvec[idx] = buffer[i*dpn+j];
1217  }
1218 
1219  vec.release_ptr(pvec);
1220 }
1221 
1222 template<class T, class S> inline
1223 void get_transformed_pts(const element_view<T,S> & view, Point* loc_pts, Point & trsf_fibre)
1224 {
1225  const elem_t type = view.type();
1226 
1227  switch(type) {
1228  case Line: {
1229  Point p0 = view.coord(0), p1 = view.coord(1);
1230 
1231  loc_pts[0] = {0, 0, 0};
1232  loc_pts[1] = {mag(p1 - p0), 0, 0};
1233  trsf_fibre = {1, 0, 0};
1234  break;
1235  }
1236 
1237  case Tri: {
1238  Point p0 = view.coord(0), p1 = view.coord(1), p2 = view.coord(2);
1239  Point f = trsf_fibre;
1240 
1241  Point p01 = p1 - p0, p02 = p2 - p0;
1242 
1243  Point x = normalize(p01);
1244  Point z = normalize(cross(p01,p02));
1245  Point y = cross(z, x);
1246 
1247  loc_pts[0] = {0, 0, 0};
1248  loc_pts[1] = {mag(p01), 0, 0};
1249  loc_pts[2] = {inner_prod(p02, x), inner_prod(p02, y), 0};
1250  trsf_fibre = {inner_prod(f, x), inner_prod(f, y), 0};
1251 
1252  if((fabs(trsf_fibre.x) + fabs(trsf_fibre.y)) < 1e-8) {
1253 // fprintf(stderr, "Fibre direction is orthogonal to triangle. Assigning (1,0,0) fiber direction.\n");
1254  trsf_fibre = {1, 0, 0};
1255  }
1256  else trsf_fibre = normalize(trsf_fibre);
1257  break;
1258  }
1259 
1260  case Quad: {
1261  Point p0 = view.coord(0), p1 = view.coord(1), p2 = view.coord(2), p3 = view.coord(3);
1262  Point f = trsf_fibre;
1263 
1264  Point p01 = p1 - p0, p02 = p2 - p0, p03 = p3 - p0;
1265 
1266  Point x = normalize(p01);
1267  Point z = normalize(cross(p01, p02));
1268  Point y = cross(z, x);
1269 
1270  loc_pts[0] = {0, 0, 0};
1271  loc_pts[1] = {mag(p01), 0, 0};
1272  loc_pts[2] = {inner_prod(p02, x), inner_prod(p02, y), 0};
1273  loc_pts[3] = {inner_prod(p03, x), inner_prod(p03, y), 0};
1274 
1275  trsf_fibre = {inner_prod(f, x), inner_prod(f, y), 0};
1276 
1277  if((fabs(trsf_fibre.x) + fabs(trsf_fibre.y)) < 1e-8) {
1278  fprintf(stderr, "Fibre direction is orthogonal to quad. Assigning (1,0,0) fiber direction.\n");
1279  trsf_fibre = {1, 0, 0};
1280  }
1281  else trsf_fibre = normalize(trsf_fibre);
1282  break;
1283  }
1284 
1285  // for 3D elements we just copy over the coords
1286  default:
1287  for(int i=0; i<view.num_nodes(); i++)
1288  loc_pts[i] = view.coord(i);
1289 
1290  break;
1291  }
1292 }
1293 
1294 
1296 template<class T, class S> inline
1297 void get_transformed_upd_pts(const element_view<T,S> & view, Point* loc_pts, Point & trsf_fibre)
1298 {
1299  const elem_t type = view.type();
1300 
1301  switch(type) {
1302  case Line: {
1303  Point p0 = view.coord_upd(0), p1 = view.coord_upd(1);
1304 
1305  loc_pts[0] = {0, 0, 0};
1306  loc_pts[1] = {mag(p1 - p0), 0, 0};
1307  trsf_fibre = {1, 0, 0};
1308  break;
1309  }
1310 
1311  case Tri: {
1312  Point p0 = view.coord_upd(0), p1 = view.coord_upd(1), p2 = view.coord_upd(2);
1313  Point f = trsf_fibre;
1314 
1315  Point p01 = p1 - p0, p02 = p2 - p0;
1316 
1317  Point x = normalize(p01);
1318  Point z = normalize(cross(p01,p02));
1319  Point y = cross(z, x);
1320 
1321  loc_pts[0] = {0, 0, 0};
1322  loc_pts[1] = {mag(p01), 0, 0};
1323  loc_pts[2] = {inner_prod(p02, x), inner_prod(p02, y), 0};
1324  trsf_fibre = {inner_prod(f, x), inner_prod(f, y), 0};
1325 
1326  if((fabs(trsf_fibre.x) + fabs(trsf_fibre.y)) < 1e-8) {
1327 // fprintf(stderr, "Fibre direction is orthogonal to triangle. Assigning (1,0,0) fiber direction.\n");
1328  trsf_fibre = {1, 0, 0};
1329  }
1330  else trsf_fibre = normalize(trsf_fibre);
1331  break;
1332  }
1333 
1334  case Quad: {
1335  Point p0 = view.coord_upd(0), p1 = view.coord_upd(1), p2 = view.coord_upd(2), p3 = view.coord_upd(3);
1336  Point f = trsf_fibre;
1337 
1338  Point p01 = p1 - p0, p02 = p2 - p0, p03 = p3 - p0;
1339 
1340  Point x = normalize(p01);
1341  Point z = normalize(cross(p01, p02));
1342  Point y = cross(z, x);
1343 
1344  loc_pts[0] = {0, 0, 0};
1345  loc_pts[1] = {mag(p01), 0, 0};
1346  loc_pts[2] = {inner_prod(p02, x), inner_prod(p02, y), 0};
1347  loc_pts[3] = {inner_prod(p03, x), inner_prod(p03, y), 0};
1348 
1349  trsf_fibre = {inner_prod(f, x), inner_prod(f, y), 0};
1350 
1351  if((fabs(trsf_fibre.x) + fabs(trsf_fibre.y)) < 1e-8) {
1352  fprintf(stderr, "Fibre direction is orthogonal to quad. Assigning (1,0,0) fiber direction.\n");
1353  trsf_fibre = {1, 0, 0};
1354  }
1355  else trsf_fibre = normalize(trsf_fibre);
1356  break;
1357  }
1358 
1359  // for 3D elements we just copy over the coords
1360  default:
1361  for(int i=0; i<view.num_nodes(); i++)
1362  loc_pts[i] = view.coord_upd(i);
1363 
1364  break;
1365  }
1366 }
1367 
1368 template<class T, class S> inline
1370 {
1371  SF::Point a, b;
1372  elem_t type = elem.type();
1373 
1374  switch (type) {
1375  case Tri:
1376  area = 0.5 * mag(cross(elem.coord(1) - elem.coord(0), elem.coord(2) - elem.coord(0)));
1377  break;
1378 
1379  case Quad:
1380  area = 0.5 * mag(cross(elem.coord(3) - elem.coord(1), elem.coord(2) - elem.coord(0)));
1381  break;
1382 
1383  default:
1384  fprintf(stderr, "You are trying to calculate the area of a surface element, but the input is not a 2-dimensional element");
1385  break;
1386  }
1387 }
1388 
1389 }
1390 
1391 #endif
Basic containers.
#define SF_MAX_ELEM_NODES
max #nodes defining an element
Definition: SF_fem_utils.h:40
double SF_real
Use the general double as real type.
Definition: SF_globals.h:38
std::int32_t SF_int
Use the general std::int32_t as int type.
Definition: SF_globals.h:37
virtual void finish_assembly()=0
virtual void set_values(const vector< T > &row_idx, const vector< T > &col_idx, const vector< S > &vals, bool add)=0
virtual void set_value(T row_idx, T col_idx, S val, bool add)=0
virtual void get(const vector< T > &idx, S *out)=0
virtual S * ptr()=0
virtual void release_ptr(S *&p)=0
ltype layout
used vector layout (nodal, algebraic, unset)
virtual void set(const vector< T > &idx, const vector< S > &vals, const bool additive=false)=0
int dpn
d.o.f. per mesh vertex.
virtual void finish_assembly()=0
const meshdata< mesh_int_t, mesh_real_t > * mesh
the connected mesh
S * data()
Definition: dense_mat.hpp:264
Comfort class. Provides getter functions to access the mesh member variables more comfortably.
Definition: SF_fem_utils.h:705
const T & node(short nidx) const
Access the connectivity information.
Definition: SF_fem_utils.h:795
Point fiber() const
Get element fiber direction.
Definition: SF_fem_utils.h:874
void integration_points(const short order, Point *ip, double *w, int &nint) const
Definition: SF_fem_utils.h:937
short num_dof(short order) const
Definition: SF_fem_utils.h:932
const T * nodes() const
Access the connectivity information.
Definition: SF_fem_utils.h:831
bool next()
Select next element if possible.
Definition: SF_fem_utils.h:747
void set_elem(size_t eidx)
Set the view to a new element.
Definition: SF_fem_utils.h:733
element_view(const meshdata< T, S > &mesh, const SF_nbr nbr)
Constructor. Initializes to element index 0.
Definition: SF_fem_utils.h:721
elem_t type() const
Getter function for the element type.
Definition: SF_fem_utils.h:773
const T & global_node(short nidx) const
Access the connectivity information.
Definition: SF_fem_utils.h:807
Point coord_upd(short nidx) const
Access updated vertex coordinates.
Definition: SF_fem_utils.h:854
const T & global_node(short nidx, SF_nbr nbr) const
Access the connectivity information.
Definition: SF_fem_utils.h:820
T tag() const
Getter function for the element tag.
Definition: SF_fem_utils.h:783
size_t global_element_index(SF_nbr nbr) const
Get currently selected element index.
Definition: SF_fem_utils.h:926
Point coord(short nidx) const
Access vertex coordinates.
Definition: SF_fem_utils.h:842
size_t global_element_index() const
Get currently selected element index.
Definition: SF_fem_utils.h:916
size_t element_index() const
Get currently selected element index.
Definition: SF_fem_utils.h:906
short dimension() const
Definition: SF_fem_utils.h:942
bool has_sheet() const
Check if a sheet direction is present.
Definition: SF_fem_utils.h:896
T num_nodes() const
Getter function for the number of nodes.
Definition: SF_fem_utils.h:763
void upd_coords(short nidx, double *pos_upd)
ACHTUNG: Update der Koordinaten durch const Qualifier des meshes hier nicht möglich.
Definition: SF_fem_utils.h:861
Point sheet() const
Get element sheet direction.
Definition: SF_fem_utils.h:883
Abstract matrix integration base class.
Definition: SF_fem_utils.h:968
virtual void operator()(const element_view< T, S > &elem, dmat< double > &buff)=0
compute the element matrix for a given element.
virtual void dpn(T &row_dpn, T &col_dpn)=0
ACHTUNG: HIER AENDERUNG VON MIR.
The mesh storage class. It contains both element and vertex data.
Definition: SF_container.h:385
vector< T > dsp
connectivity starting index of each element
Definition: SF_container.h:405
vector< S > she
sheet direction
Definition: SF_container.h:410
vector< S > fib
fiber direction
Definition: SF_container.h:409
size_t l_numelem
local number of elements
Definition: SF_container.h:388
vector< elem_t > type
element type
Definition: SF_container.h:407
vector< T > con
Definition: SF_container.h:401
vector< S > xyz_Euler
node coordinates in Euler formulation
Definition: SF_container.h:417
vector< S > xyz
node coordinates in Lagrange formulation
Definition: SF_container.h:416
MPI_Comm comm
the parallel mesh is defined on a MPI world
Definition: SF_container.h:393
vector< T > & get_numbering(SF_nbr nbr_type)
Get the vector defining a certain numbering.
Definition: SF_container.h:454
vector< T > tag
element tag
Definition: SF_container.h:406
non_overlapping_layout< T > epl
element parallel layout
Definition: SF_container.h:420
Abstract vector integration base class.
Definition: SF_fem_utils.h:986
void zero_buff(double *buff, T nrows)
Definition: SF_fem_utils.h:988
virtual void operator()(const element_view< T, S > &elem, double *buff)=0
compute the element matrix for a given element.
virtual void dpn(T &dpn)=0
return (by reference) the row and column dimensions
A vector storing arbitrary data.
Definition: SF_vector.h:43
size_t size() const
The current size of the vector.
Definition: SF_vector.h:104
void resize(size_t n)
Resize a vector.
Definition: SF_vector.h:209
T * data()
Pointer to the vector's start.
Definition: SF_vector.h:91
Definition: dense_mat.hpp:34
double mag(const Point &vect)
vector magnitude
Definition: SF_container.h:112
void get_transformed_pts(const element_view< T, S > &view, Point *loc_pts, Point &trsf_fibre)
void invert_3x3(S *ele, S &det)
Definition: dense_mat.hpp:452
void shape_deriv(const double *iJ, const dmat< double > &rshape, const int ndof, dmat< double > &shape)
Compute shape derivatives for an element, based on the shape derivatives of the associated reference ...
Definition: SF_fem_utils.h:674
double inner_prod(const Point &a, const Point &b)
Definition: SF_container.h:91
void assemble_vector(abstract_vector< T, S > &vec, meshdata< mesh_int_t, mesh_real_t > &domain, vector_integrator< mesh_int_t, mesh_real_t > &integrator)
Generalized vector assembly.
void get_area_of_surface_elem(SF::element_view< T, S > elem, S &area)
void extract_element_data(const element_view< mesh_int_t, mesh_real_t > &view, abstract_vector< T, S > &vec, SF_real *buffer)
void canonic_indices(const T *nidx, const T *nbr, const T esize, const short dpn, V *cidx)
Compute canonical indices from nodal indices and dpn.
dmat< S > invert_2x2(const dmat< S > &m)
Definition: dense_mat.hpp:508
void assemble_matrix(abstract_matrix< T, S > &mat, meshdata< mesh_int_t, mesh_real_t > &domain, matrix_integrator< mesh_int_t, mesh_real_t > &integrator)
Generalized matrix assembly.
void jacobian_matrix(const dmat< double > &rshape, const int npts, const Point *pts, double *J)
Compute Jacobian matrix from the real element to the reference element.
Definition: SF_fem_utils.h:610
short num_dof(elem_t type, short order)
Get number of d.o.f. for an element type and an Ansatz function order.
Definition: SF_fem_utils.h:58
Point normalize(const Point &vect)
Definition: SF_container.h:117
void assemble_lumped_matrix(abstract_matrix< T, S > &mat, meshdata< mesh_int_t, mesh_real_t > &domain, matrix_integrator< mesh_int_t, mesh_real_t > &integrator)
void set_element_data(const element_view< mesh_int_t, mesh_real_t > &view, SF_real *buffer, abstract_vector< T, S > &vec)
Point cross(const Point &a, const Point &b)
cross product
Definition: SF_container.h:85
void get_transformed_upd_pts(const element_view< T, S > &view, Point *loc_pts, Point &trsf_fibre)
ACHTUNG: DAS WURDE VON MIR HINZUGEFÜGT.
elem_t
element type enum
Definition: SF_container.h:54
@ Line
Definition: SF_container.h:62
@ Tri
Definition: SF_container.h:61
@ Prism
Definition: SF_container.h:59
@ Pyramid
Definition: SF_container.h:58
@ Tetra
Definition: SF_container.h:55
@ Quad
Definition: SF_container.h:60
@ Hexa
Definition: SF_container.h:56
void general_integration_points(const elem_t type, const short order, Point *ip, double *w, int &nint)
Compute the integration point locations and weights.
Definition: SF_fem_utils.h:121
void reference_shape(const elem_t type, const Point ip, dmat< double > &rshape)
Compute shape function and its derivatives on a reference element.
Definition: SF_fem_utils.h:383
SF_nbr
Enumeration encoding the different supported numberings.
Definition: SF_container.h:189
@ NBR_PETSC
PETSc numbering of nodes.
Definition: SF_container.h:192
void invert_jacobian_matrix(const elem_t type, double *J, double &detJ)
Definition: SF_fem_utils.h:633
Point and vector struct.
Definition: SF_container.h:66
double y
Definition: SF_container.h:68
double z
Definition: SF_container.h:69
double x
Definition: SF_container.h:67