Flow123d  jenkins-Flow123d-linux-release-multijob-282
mh_fe_values.cc
Go to the documentation of this file.
1 /*
2  * mh_fe_values.cc
3  *
4  * Created on: Jul 2, 2012
5  * Author: jb
6  */
7 
8 
9 #include "system/system.hh"
10 #include "system/math_fce.h"
11 #include "mesh/mesh.h"
12 
13 #include "mh_fe_values.hh"
14 
15 
17  // allocate temp. arrays
18  loc_matrix_ = new double [4*4 + 4*4 + 4*4];
20  bas_alfa = inv_loc_matrix_ + 4*4;
21  bas_beta = bas_alfa + 4;
22  bas_gama = bas_beta +4;
23  bas_delta = bas_gama +4;
24 }
25 
26 
27 
29  delete [] loc_matrix_;
30 }
31 
32 
33 
34 void MHFEValues::update(ElementFullIter ele, FieldType &anisotropy, FieldType_Scalar &cross_section, FieldType_Scalar &conductivity) {
35 
36  ASSERT(!( ele == NULL ),"NULL as argument of function local_matrix()\n");
37 
38  double scale = 1/ conductivity.value( ele->centre(), ele->element_accessor() ) / cross_section.value( ele->centre(), ele->element_accessor() );
39  switch( ele->dim() ) {
40  case 1:
41  local_matrix_line( ele, anisotropy , scale);
42  break;
43  case 2:
44  local_matrix_triangle( ele, anisotropy, scale);
45  break;
46  case 3:
47  local_matrix_tetrahedron( ele, anisotropy, scale );
48  break;
49  }
50 
51  // matrix inverse
52 
53  double det = MatrixInverse(loc_matrix_, inv_loc_matrix_, ele->n_sides() );
54  if (fabs(det) < NUM_ZERO) {
55  xprintf(Warn,"Singular local matrix of the element %d\n",ele.id());
56  PrintSmallMatrix(loc_matrix_, ele->n_sides());
57  xprintf(Err,"det: %30.18e \n",det);
58  }
59 }
60 
61 
62 
64  return loc_matrix_;
65 }
66 
67 
68 
70  return inv_loc_matrix_;
71 }
72 
73 
74 
76  switch( ele->dim() ) {
77  case 1:
78  {
79 
80  arma::vec3 line_vec = ele->node[1]->point() - ele->node[0]->point();
81  if (face == 0) {
82  return - arma::norm( point - ele->node[1]->point(), 2) * line_vec / arma::dot( line_vec, line_vec) ;
83  } else {
84  return arma::norm( point - ele->node[0]->point(), 2) * line_vec / arma::dot( line_vec, line_vec) ;
85  }
86  }
87  case 2:
88  {
89  // make rotated coordinate system with triangle in plane XY, origin in A and axes X == AB
90  arma::vec3 ex(ele->node[1]->point() - ele->node[0]->point());
91  ex /= norm(ex,2);
92 
93  arma::vec3 ac(ele->node[2]->point() - ele->node[0]->point());
94  arma::vec3 ez = cross(ex, ac);
95  ez /= norm(ez,2);
96 
97  arma::vec3 ey = cross(ez,ex);
98  ey /= norm(ey, 2);
99 
100  // compute point in new coordinate system
101  arma::vec3 u = point - ele->node[0]->point();
102 
103  // compute vector value from the base function
104  return bas_gama[ face ] * (
105  (dot(u, ex) - bas_alfa[ face ] ) * ex
106  +(dot(u, ey) - bas_beta[ face ] ) * ey
107  );
108  }
109  case 3:
110  {
111  arma::vec3 RT0_Y;
112  RT0_Y[0] = bas_alfa[ face ];
113  RT0_Y[1] = bas_beta[ face ];
114  RT0_Y[2] = bas_gama[ face ];
115 
116  return bas_delta[ face ]*(point - RT0_Y);
117  }
118  }
119 
120  return arma::vec3("0 0 0");
121 }
122 
123 
124 
125 
126 
127 //=============================================================================
128 // CALCULATE LOCAL MARIX FOR LINEAR ELEMENT
129 //=============================================================================
130 void MHFEValues::local_matrix_line(ElementFullIter ele, FieldType &anisotropy, double scale )
131 {
132  double val;
134 
135  //getting vector on the line and normalizing it
136  //it is the transformation matrix from 3D to 1D line of the 1D element
137  arma::vec3 line_vec = ele->node[1]->point() - ele->node[0]->point();
138  line_vec /= arma::norm(line_vec,2);
139 
140  //transforming the conductivity in 3D to resistivity in 1D
141  //computing v_transpose * K_inverse * v
142  val = scale * arma::dot(line_vec,
143  (anisotropy.value(ele->centre(), ele->element_accessor() )).i() * line_vec
144  ) * ele->measure() / 3.0;
145 
146  loc[0][0] = val;
147  loc[1][1] = val;
148  loc[0][1] = - val / 2.0;
149  loc[1][0] = - val / 2.0;
150 
151  bas_alfa[0] = bas_alfa[1] = -1.0 / ele->measure();
152  bas_beta[0] = 1.0;
153  bas_beta[1] = 0.0;
154 }
155 
156 
157 //=============================================================================
158 // CALCULATE LOCAL MATRIX FOR TRIANGULAR ELEMENT
159 //=============================================================================
160 void MHFEValues::local_matrix_triangle( ElementFullIter ele, FieldType &anisotropy, double scale )
161 {
162  double midpoint[ 3 ][ 2 ]; // Midpoints of element's sides
163  double alfa[ 3 ]; //
164  double beta[ 3 ]; // |- Parametrs of basis functions
165  double gama[ 3 ]; // /
166  double poly[ 6 ]; // Koeficients of polynomial
167  double p[ 3 ]; // Polynomial's values in midpoints
168  int i, j; // Loops' counters
169  double nod_coor[ 3 ][ 2 ]; // Coordinates of nodes;
171 
172  // make rotated coordinate system with triangle in plane XY, origin in A and axes X == AB
173  arma::vec3 ex(ele->node[1]->point() - ele->node[0]->point());
174  ex /= arma::norm(ex,2);
175 
176  arma::vec3 ac(ele->node[2]->point() - ele->node[0]->point());
177  arma::vec3 ez = arma::cross(ex, ac);
178  ez /= norm(ez,2);
179 
180  arma::vec3 ey = arma::cross(ez,ex);
181  ey /= arma::norm(ey, 2);
182 
183  //transformation matrix from 3D to 2D plane of the 2D element
184  arma::mat r(3,2);
185  r.col(0) = ex;
186  r.col(1) = ey;
187 
188  //transforming 3D conductivity tensor to 2D resistance tensor
189  arma::mat resistance_tensor = r.t() * ((anisotropy.value(ele->centre(), ele->element_accessor() )).i() * r);
190 
191  node_coordinates_triangle( ele, nod_coor );
192  side_midpoint_triangle( nod_coor, midpoint );
193  basis_functions_triangle( nod_coor, alfa, beta, gama );
194  for( i = 0; i < 3; i++ )
195  for( j = i; j < 3; j++ ) {
196  calc_polynom_triangle( alfa[ i ], beta[ i ],
197  alfa[ j ], beta[ j ],
198  resistance_tensor, poly );
199  p[ 0 ] = polynom_value_triangle( poly, midpoint[ 0 ] );
200  p[ 1 ] = polynom_value_triangle( poly, midpoint[ 1 ] );
201  p[ 2 ] = polynom_value_triangle( poly, midpoint[ 2 ] );
202  loc[i][j] = scale * gama[i] * gama[j] * ele->measure() *
203  ( p[0] + p[1] + p[2] ) / 3.0;
204  loc[j][i] = loc[i][j];
205  }
206  for(i = 0; i < 3; i++) {
207  bas_alfa[ i ] = alfa[ i ];
208  bas_beta[ i ] = beta[ i ];
209  bas_gama[ i ] = gama[ i ];
210  }
211 }
212 
213 
214 /*
215  * Computes coordinates of vertices of the triangle in the local orthogonal system
216  * that has normalized vector V0 to V1 as the first vector of the basis.
217  */
219 {
220  double u[ 3 ], v[ 3 ]; // vectors of sides 0 and 2
221  double t[ 3 ]; // u normalized (tangenta)
222  double n[ 3 ]; // normal
223  double b[ 3 ]; // binormal
224 
225  u[ 0 ] = ele->node[ 1 ]->getX() - ele->node[ 0 ]->getX();
226  u[ 1 ] = ele->node[ 1 ]->getY() - ele->node[ 0 ]->getY();
227  u[ 2 ] = ele->node[ 1 ]->getZ() - ele->node[ 0 ]->getZ();
228  v[ 0 ] = ele->node[ 2 ]->getX() - ele->node[ 0 ]->getX();
229  v[ 1 ] = ele->node[ 2 ]->getY() - ele->node[ 0 ]->getY();
230  v[ 2 ] = ele->node[ 2 ]->getZ() - ele->node[ 0 ]->getZ();
231 
232  vector_product( u, v, n );
233  t[ 0 ] = u[ 0 ];
234  t[ 1 ] = u[ 1 ];
235  t[ 2 ] = u[ 2 ];
236  normalize_vector( t );
237  normalize_vector( n );
238  vector_product( n, u, b );
239  normalize_vector( b );
240 
241  nod[ 0 ][ 0 ] = 0.0;
242  nod[ 0 ][ 1 ] = 0.0;
243  nod[ 1 ][ 0 ] = scalar_product( t, u );
244  nod[ 1 ][ 1 ] = scalar_product( b, u );
245  nod[ 2 ][ 0 ] = scalar_product( t, v );
246  nod[ 2 ][ 1 ] = scalar_product( b, v );
247 }
248 
249 
250 /**
251  * Computes coordinates of the side midpoints in the local orthogonal coordinate system.
252  */
253 void MHFEValues::side_midpoint_triangle( double nod[ 3 ][ 2 ], double midpoint[ 3 ][ 2 ] )
254 {
255  midpoint[ 0 ][ 0 ] = ( nod[ 0 ][ 0 ] + nod[ 1 ][ 0 ] ) / 2.0;
256  midpoint[ 0 ][ 1 ] = ( nod[ 0 ][ 1 ] + nod[ 1 ][ 1 ] ) / 2.0;
257  midpoint[ 1 ][ 0 ] = ( nod[ 1 ][ 0 ] + nod[ 2 ][ 0 ] ) / 2.0;
258  midpoint[ 1 ][ 1 ] = ( nod[ 1 ][ 1 ] + nod[ 2 ][ 1 ] ) / 2.0;
259  midpoint[ 2 ][ 0 ] = ( nod[ 2 ][ 0 ] + nod[ 0 ][ 0 ] ) / 2.0;
260  midpoint[ 2 ][ 1 ] = ( nod[ 2 ][ 1 ] + nod[ 0 ][ 1 ] ) / 2.0;
261 }
262 
263 
264 
265 /*
266  * Computes coefficients of the RT basis functions.
267  */
268 void MHFEValues::basis_functions_triangle( double nod[ 3 ][ 2 ], double alfa[],
269  double beta[], double gama[] )
270 {
271  bas_func_0_triangle( nod[ 0 ][ 0 ], nod[ 0 ][ 1 ],
272  nod[ 1 ][ 0 ], nod[ 1 ][ 1 ],
273  nod[ 2 ][ 0 ], nod[ 2 ][ 1 ],
274  alfa, beta, gama );
275  bas_func_0_triangle( nod[ 1 ][ 0 ], nod[ 1 ][ 1 ],
276  nod[ 2 ][ 0 ], nod[ 2 ][ 1 ],
277  nod[ 0 ][ 0 ], nod[ 0 ][ 1 ],
278  alfa + 1, beta + 1, gama + 1 );
279  bas_func_0_triangle( nod[ 2 ][ 0 ], nod[ 2 ][ 1 ],
280  nod[ 0 ][ 0 ], nod[ 0 ][ 1 ],
281  nod[ 1 ][ 0 ], nod[ 1 ][ 1 ],
282  alfa + 2, beta + 2, gama + 2 );
283 }
284 //=============================================================================
285 // GENERATE BASIS FUNCTION FOR SIDE 0
286 //=============================================================================
287 void MHFEValues::bas_func_0_triangle( double x0, double y0,
288  double x1, double y1,
289  double x2, double y2,
290  double *alfa, double *beta, double *gama)
291 {
292  *alfa = x2;
293  *beta = y2;
294  *gama = 1.0 / fabs( ( x0 - x2 ) * ( y1 - y0 ) + ( y0 - y2 ) * ( x0 - x1 ) );
295 }
296 //=============================================================================
297 // CALCULATE POLYNOM OF SCALAR PRODUCT
298 //=============================================================================
299 void MHFEValues::calc_polynom_triangle( double al_i, double be_i, double al_j, double be_j,
300  arma::mat::fixed<2,2> a, double poly[] )
301 {
302  poly[ 0 ] = a( 0,0 ) * al_i * al_j +
303  a( 0,1 ) * be_i * al_j +
304  a( 1,0 ) * al_i * be_j +
305  a( 1,1 ) * be_i * be_j;
306 
307  poly[ 1 ] = ( a( 0,0 ) * al_i +
308  a( 0,0 ) * al_j +
309  a( 0,1 ) * be_i +
310  a( 1,0 ) * be_j ) * -1.0;
311 
312  poly[ 2 ] = ( a( 1,1 ) * be_i +
313  a( 1,1 ) * be_j +
314  a( 1,0 ) * al_i +
315  a( 0,1) * al_j ) * -1.0;
316  poly[ 3 ] = a( 0,0 );
317  poly[ 4 ] = a( 0,1 ) + a( 1,0 );
318  poly[ 5 ] = a( 1,1 );
319 }
320 //=============================================================================
321 // CALCULATE VALUE OF POLYNOM IN GIVEN POINT
322 //=============================================================================
323 double MHFEValues::polynom_value_triangle( double poly[], double point[] )
324 {
325  double rc;
326 
327  rc = poly[ 0 ] +
328  poly[ 1 ] * point[ 0 ] +
329  poly[ 2 ] * point[ 1 ] +
330  poly[ 3 ] * point[ 0 ] * point[ 0 ] +
331  poly[ 4 ] * point[ 0 ] * point[ 1 ] +
332  poly[ 5 ] * point[ 1 ] * point[ 1 ];
333  return rc;
334 }
335 
336 //=============================================================================
337 // CALCULATE LOCAL MARIX FOR SIMPLEX ELEMENT
338 //=============================================================================
339 void MHFEValues::local_matrix_tetrahedron( ElementFullIter ele, FieldType &anisotropy, double scale )
340 {
341  double alfa[ 4 ];
342  double beta[ 4 ]; // | Parametrs of basis functions
343  double gama[ 4 ]; // |
344  double delta[ 4 ]; // /
345  double poly[ 10 ]; // Koeficients of polynomial
346  int i, j; // Loops' counters
348 
349  //transforming 3D conductivity tensor to 3D resistance tensor
350  arma::mat resistance_tensor = (anisotropy.value(ele->centre(), ele->element_accessor() )).i();
351 
352  //OBSOLETE
353  //SmallMtx3 resistance_tensor = (SmallMtx3)(ele->material->hydrodynamic_resistence);
354  //compares old and new resistance tensor
355  //DBGMSG("my: %f %f %f %f\t orig: %f %f %f %f\n",
356  // my_resistance_tensor(0,0), my_resistance_tensor(0,1), my_resistance_tensor(1,0), my_resistance_tensor(1,1),
357  // resistance_tensor[0][0], resistance_tensor[0][1], resistance_tensor[1][0], resistance_tensor[1][1] );
358 
359  basis_functions_tetrahedron( ele, alfa, beta, gama, delta);
360  for( i = 0; i < 4; i++ )
361  for( j = i; j < 4; j++ ) {
363  alfa[ i ], beta[ i ], gama[ i ],
364  alfa[ j ], beta[ j ], gama[ j ],
365  resistance_tensor, poly );
366  loc[i][j] = scale * delta[i] * delta[j] * polynom_integral_tetrahedron( ele, poly );
367  loc[j][i] = loc[i][j];
368  }
369  for(i = 0; i < 4; i++) {
370  bas_alfa[ i ] = alfa[ i ];
371  bas_beta[ i ] = beta[ i ];
372  bas_gama[ i ] = gama[ i ];
373  bas_delta[ i ] = delta[ i ];
374  }
375 }
376 //=============================================================================
377 //
378 //=============================================================================
379 void MHFEValues::basis_functions_tetrahedron( ElementFullIter ele, double alfa[], double beta[],double gama[], double delta[])
380 {
381  int li;
382  SideIter pSid;
383 
384  for( li = 0; li < 4; li++ ) {
385  alfa[ li ] = ele->node[ li ]->getX();
386  beta[ li ] = ele->node[ li ]->getY();
387  gama[ li ] = ele->node[ li ]->getZ();
388  pSid = ele->side( li );
389  delta[ li ] = 1.0 / ( pSid->measure() *
390  ( pSid->normal()[ 0 ] * pSid->centre()[ 0 ] +
391  pSid->normal()[ 1 ] * pSid->centre()[ 1 ] +
392  pSid->normal()[ 2 ] * pSid->centre()[ 2 ] -
393  pSid->normal()[ 0 ] * alfa[ li ] -
394  pSid->normal()[ 1 ] * beta[ li ] -
395  pSid->normal()[ 2 ] * gama[ li ] ) );
396  }
397 }
398 //=============================================================================
399 // CALCULATE POLYNOM OF SCALAR PRODUCT
400 //=============================================================================
401 void MHFEValues::calc_polynom_tetrahedron( double al_i, double be_i, double ga_i,
402  double al_j, double be_j, double ga_j,
403  arma::mat::fixed<3,3> a, double poly[] )
404 {
405  // Constant term
406  poly[ 0 ] = a( 0,0 ) * al_i * al_j +
407  a( 0,1 ) * be_i * al_j +
408  a( 0,2 ) * ga_i * al_j +
409  a( 1,0 ) * al_i * be_j +
410  a( 1,1 ) * be_i * be_j +
411  a( 1,2 ) * ga_i * be_j +
412  a( 2,0 ) * al_i * ga_j +
413  a( 2,1 ) * be_i * ga_j +
414  a( 2,2 ) * ga_i * ga_j;
415  // Term with x
416  poly[ 1 ] = ( a( 0,0 ) * al_i +
417  a( 0,1 ) * be_i +
418  a( 0,2 ) * ga_i +
419  a( 0,0 ) * al_j +
420  a( 1,0 ) * be_j +
421  a( 2,0 ) * ga_j ) * -1.0;
422  // Term with y
423  poly[ 2 ] = ( a( 1,0 ) * al_i +
424  a( 1,1 ) * be_i +
425  a( 1,2 ) * ga_i +
426  a( 0,1 ) * al_j +
427  a( 1,1 ) * be_j +
428  a( 2,1 ) * ga_j ) * -1.0;
429  // Term with z
430  poly[ 3 ] = ( a( 2,0 ) * al_i +
431  a( 2,1 ) * be_i +
432  a( 2,2 ) * ga_i +
433  a( 0,2 ) * al_j +
434  a( 1,2 ) * be_j +
435  a( 2,2 ) * ga_j ) * -1.0;
436  // Term with xy
437  poly[ 4 ] = a( 0,1 ) + a( 1,0 );
438  // Term with xz
439  poly[ 5 ] = a( 0,2 ) + a( 2,0 );
440  // Term with yz
441  poly[ 6 ] = a( 1,2 ) + a( 2,1 );
442  // Term with x^2
443  poly[ 7 ] = a( 0,0 );
444  // Term with y^2
445  poly[ 8 ] = a( 1,1 );
446  // Term with z^2
447  poly[ 9 ] = a( 2,2 );
448 }
449 //=============================================================================
450 // CALCULATE INTEGRAL OF QUADRATICAL POLYNOM OVER TETRAHEDRON
451 //=============================================================================
453 {
454  double rc;
455  double v;
456  double t[ 3 ];
457 
458  rc = 0.0;
459  v = ele->measure();
460  t[ 0 ] = ele->centre()[ 0 ];
461  t[ 1 ] = ele->centre()[ 1 ];
462  t[ 2 ] = ele->centre()[ 2 ];
463  // Constant
464  rc += poly[ 0 ] * v;
465  // Term with x
466  rc += poly[ 1 ] * v * t[ 0 ];
467  // Term with y
468  rc += poly[ 2 ] * v * t[ 1 ];
469  // Term with z
470  rc += poly[ 3 ] * v * t[ 2 ];
471  // Term with xy
472  rc += poly[ 4 ] * v * ( t[ 0 ] * t[ 1 ] + 1.0 / 20.0 * (
473  ( ele->node[ 0 ]->getX() - t[ 0 ] ) * ( ele->node[ 0 ]->getY() - t[ 1 ] ) +
474  ( ele->node[ 1 ]->getX() - t[ 0 ] ) * ( ele->node[ 1 ]->getY() - t[ 1 ] ) +
475  ( ele->node[ 2 ]->getX() - t[ 0 ] ) * ( ele->node[ 2 ]->getY() - t[ 1 ] ) +
476  ( ele->node[ 3 ]->getX() - t[ 0 ] ) * ( ele->node[ 3 ]->getY() - t[ 1 ] )
477  ) );
478  // Term with xz
479  rc += poly[ 5 ] * v * ( t[ 0 ] * t[ 2 ] + 1.0 / 20.0 * (
480  ( ele->node[ 0 ]->getX() - t[ 0 ] ) * ( ele->node[ 0 ]->getZ() - t[ 2 ] ) +
481  ( ele->node[ 1 ]->getX() - t[ 0 ] ) * ( ele->node[ 1 ]->getZ() - t[ 2 ] ) +
482  ( ele->node[ 2 ]->getX() - t[ 0 ] ) * ( ele->node[ 2 ]->getZ() - t[ 2 ] ) +
483  ( ele->node[ 3 ]->getX() - t[ 0 ] ) * ( ele->node[ 3 ]->getZ() - t[ 2 ] )
484  ) );
485  // Term with yz
486  rc += poly[ 6 ] * v * ( t[ 1 ] * t[ 2 ] + 1.0 / 20.0 * (
487  ( ele->node[ 0 ]->getY() - t[ 1 ] ) * ( ele->node[ 0 ]->getZ() - t[ 2 ] ) +
488  ( ele->node[ 1 ]->getY() - t[ 1 ] ) * ( ele->node[ 1 ]->getZ() - t[ 2 ] ) +
489  ( ele->node[ 2 ]->getY() - t[ 1 ] ) * ( ele->node[ 2 ]->getZ() - t[ 2 ] ) +
490  ( ele->node[ 3 ]->getY() - t[ 1 ] ) * ( ele->node[ 3 ]->getZ() - t[ 2 ] )
491  ) );
492  // Term with x^2
493  rc += poly[ 7 ] * v * ( t[ 0 ] * t[ 0 ] + 1.0 / 20.0 * (
494  ( ele->node[ 0 ]->getX() - t[ 0 ] ) * ( ele->node[ 0 ]->getX() - t[ 0 ] ) +
495  ( ele->node[ 1 ]->getX() - t[ 0 ] ) * ( ele->node[ 1 ]->getX() - t[ 0 ] ) +
496  ( ele->node[ 2 ]->getX() - t[ 0 ] ) * ( ele->node[ 2 ]->getX() - t[ 0 ] ) +
497  ( ele->node[ 3 ]->getX() - t[ 0 ] ) * ( ele->node[ 3 ]->getX() - t[ 0 ] )
498  ) );
499  // Term with y^2
500  rc += poly[ 8 ] * v * ( t[ 1 ] * t[ 1 ] + 1.0 / 20.0 * (
501  ( ele->node[ 0 ]->getY() - t[ 1 ] ) * ( ele->node[ 0 ]->getY() - t[ 1 ] ) +
502  ( ele->node[ 1 ]->getY() - t[ 1 ] ) * ( ele->node[ 1 ]->getY() - t[ 1 ] ) +
503  ( ele->node[ 2 ]->getY() - t[ 1 ] ) * ( ele->node[ 2 ]->getY() - t[ 1 ] ) +
504  ( ele->node[ 3 ]->getY() - t[ 1 ] ) * ( ele->node[ 3 ]->getY() - t[ 1 ] )
505  ) );
506  // Term with z^2
507  rc += poly[ 9 ] * v * ( t[ 2 ] * t[ 2 ] + 1.0 / 20.0 * (
508  ( ele->node[ 0 ]->getZ() - t[ 2 ] ) * ( ele->node[ 0 ]->getZ() - t[ 2 ] ) +
509  ( ele->node[ 1 ]->getZ() - t[ 2 ] ) * ( ele->node[ 1 ]->getZ() - t[ 2 ] ) +
510  ( ele->node[ 2 ]->getZ() - t[ 2 ] ) * ( ele->node[ 2 ]->getZ() - t[ 2 ] ) +
511  ( ele->node[ 3 ]->getZ() - t[ 2 ] ) * ( ele->node[ 3 ]->getZ() - t[ 2 ] )
512  ) );
513  return rc;
514 }
void calc_polynom_tetrahedron(double al_i, double be_i, double ga_i, double al_j, double be_j, double ga_j, arma::mat::fixed< 3, 3 > a, double poly[])
void vector_product(double u[], double v[], double x[])
Definition: math_fce.cc:80
void node_coordinates_triangle(ElementFullIter ele, double nod[3][2])
double polynom_value_triangle(double poly[], double point[])
void local_matrix_tetrahedron(ElementFullIter ele, FieldType &cond_anisothropy, double scale)
double * bas_gama
Definition: mh_fe_values.hh:68
double polynom_integral_tetrahedron(ElementFullIter ele, double poly[])
void bas_func_0_triangle(double x0, double y0, double x1, double y1, double x2, double y2, double *alfa, double *beta, double *gama)
void basis_functions_triangle(double nod[3][2], double alfa[], double beta[], double gama[])
SmallVec3_t * SmallMtx3
Definition: math_fce.h:43
arma::vec3 centre() const
Definition: sides.cc:135
???
void local_matrix_line(ElementFullIter ele, FieldType &cond_anisothropy, double scale)
SmallVec4_t * SmallMtx4
Definition: math_fce.h:44
void basis_functions_tetrahedron(ElementFullIter ele, double alfa[], double beta[], double gama[], double delta[])
double * bas_delta
Definition: mh_fe_values.hh:69
void normalize_vector(double u[])
Definition: math_fce.cc:60
int id()
Returns id of the iterator. See VectorId documentation.
Definition: sys_vector.hh:179
#define ASSERT(...)
Definition: global_defs.h:121
Definition: system.hh:72
double * inv_loc_matrix_
Definition: mh_fe_values.hh:62
void update(ElementFullIter ele, FieldType &cond_anisothropy, FieldType_Scalar &cross_section, FieldType_Scalar &conductivity)
Definition: mh_fe_values.cc:34
double * inv_local_matrix()
Definition: mh_fe_values.cc:69
double scalar_product(double u[], double v[])
Definition: math_fce.cc:53
void PrintSmallMatrix(double *mtx, int size)
Definition: math_fce.cc:226
#define xprintf(...)
Definition: system.hh:100
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:312
void side_midpoint_triangle(double nod[3][2], double midpoint[3][2])
double * bas_beta
Definition: mh_fe_values.hh:67
double measure() const
Definition: sides.cc:42
#define NUM_ZERO
Numerical helpers.
Definition: math_fce.h:53
double * local_matrix()
Definition: mh_fe_values.cc:63
void calc_polynom_triangle(double al_i, double be_i, double al_j, double be_j, arma::mat::fixed< 2, 2 > a, double poly[])
double MatrixInverse(double *a, double *b, int size)
Definition: math_fce.cc:116
arma::vec3 RT0_value(ElementFullIter ele, arma::vec3 point, unsigned int face)
Definition: mh_fe_values.cc:75
void local_matrix_triangle(ElementFullIter ele, FieldType &cond_anisothropy, double scale)
arma::vec3 normal() const
Definition: sides.cc:65
SmallVec2_t * SmallMtx2
Definition: mh_fe_values.hh:22
double * bas_alfa
Definition: mh_fe_values.hh:66
double * loc_matrix_
Definition: mh_fe_values.hh:61
Definition: system.hh:72