Flow123d  DF_patch_fe_data_tables-f6c5b2e
element_values.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2019 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file fe_values.cc
15  * @brief Class ElementValues calculates data related to transformation
16  * of reference cell to actual cell (Jacobian, inverse Jacobian,
17  * determinant, point coordinates etc.).
18  * @author Jan Stebel
19  */
20 
21 #include "fem/mapping_p1.hh"
22 #include "quadrature/quadrature.hh"
23 #include "fem/element_values.hh"
24 #include "fem/fem_tools.hh"
25 
26 
27 
28 using namespace arma;
29 using namespace std;
30 
31 
32 
33 
34 
35 
36 
37 
38 template<unsigned int spacedim>
40  Quadrature &_quadrature,
41  unsigned int dim)
42 : dim_(dim),
43  n_points_(_quadrature.size()),
44  n_sides_(_quadrature.dim() == dim ? 0 : dim+1),
45  ref_data(nullptr),
46  side_ref_data(n_sides_)
47 {
48 
49 }
50 
51 
53  : n_points(np)
54 {
55  bar_coords.resize(np);
56  weights.resize(np);
57 }
58 
59 
60 template<unsigned int spacedim>
61 ElementData<spacedim>::ElementData(unsigned int size,
62  UpdateFlags flags,
63  unsigned int dim)
64 : dim_(dim),
65  JxW_values(flags & update_JxW_values ? size : 0),
66  side_JxW_values(flags & update_side_JxW_values ? size : 0),
67  jacobians(spacedim, dim, flags & update_jacobians ? size : 0),
68  determinants(flags & update_volume_elements ? size : 0),
69  inverse_jacobians(dim, spacedim, flags & update_inverse_jacobians ? size : 0),
70  points(spacedim, 1, flags & update_quadrature_points ? size : 0),
71  normal_vectors(spacedim, 1, flags & update_normal_vectors ? size : 0),
72  update_flags(flags)
73 {}
74 
75 
76 template<unsigned int spacedim>
78 {
79  if (cell.is_valid() || side.is_valid())
80  {
81  if (cell.is_valid())
82  printf("cell %d dim %d ", cell.elm_idx(), cell.dim());
83  else if (side.is_valid())
84  printf("cell %d dim %d side %d ", side.elem_idx(), side.dim(), side.side_idx());
85 
86  printf(" det[");
87  for (auto d : determinants) printf("%g ", d); printf("]");
88 
89  printf(" JxW[");
90  for (auto j : JxW_values) printf("%g ", j); printf("]");
91 
92  printf(" side_JxW[");
93  for (auto j : side_JxW_values) printf("%g ", j); printf("]");
94 
95  printf(" nv[");
96  for (unsigned int i=0; i<normal_vectors.size(); i++)
97  {
98  auto n = normal_vectors.vec<spacedim>(i);
99  printf(" [");
100  for (unsigned int c=0; c<spacedim; c++) printf("%g ", n[c]);
101  printf("]");
102  }
103  printf("]\n");
104  }
105 }
106 
107 
108 
109 
110 
111 template<unsigned int spacedim>
113 {
114  ASSERT( q.dim() == dim_ );
115  ASSERT( q.size() == n_points_ );
116  RefElementData *ref_data = new RefElementData(q.size());
117 
118  for (unsigned int i=0; i<q.size(); i++)
119  {
120  switch (q.dim())
121  {
122  case 0:
123  ref_data->bar_coords[i] = arma::vec( "0" );
124  //ref_data->bar_coords[i] = RefElement<0>::local_to_bary(q.point<0>(i));
125  break;
126  case 1:
127  ref_data->bar_coords[i] = RefElement<1>::local_to_bary(q.point<1>(i));
128  break;
129  case 2:
130  ref_data->bar_coords[i] = RefElement<2>::local_to_bary(q.point<2>(i));
131  break;
132  case 3:
133  ref_data->bar_coords[i] = RefElement<3>::local_to_bary(q.point<3>(i));
134  break;
135  default:
136  ASSERT_PERMANENT(false)(q.dim()).error("Unsupported dimension.\n");
137  break;
138  }
139  ref_data->weights[i] = q.weight(i);
140  }
141 
142  return ref_data;
143 }
144 
145 
146 
147 template<unsigned int spacedim>
149 {
150  switch (this->dim_)
151  {
152  case 0:
153  flags = MappingP1<0,spacedim>::update_each(flags);
154  break;
155  case 1:
156  flags = MappingP1<1,spacedim>::update_each(flags);
157  break;
158  case 2:
159  flags = MappingP1<2,spacedim>::update_each(flags);
160  break;
161  case 3:
162  flags = MappingP1<3,spacedim>::update_each(flags);
163  break;
164  default:
165  ASSERT_PERMANENT(false)(this->dim_).error("Unsupported dimension.\n");
166  break;
167  }
168  return flags;
169 }
170 
171 
172 
173 template<unsigned int spacedim>
175  Quadrature &_quadrature,
176  UpdateFlags _flags,
177  unsigned int dim)
178 : RefElementValues<spacedim>(_quadrature, dim),
179  data(this->n_points_, update_each(_flags), dim)
180 {
181  this->ref_initialize(_quadrature, dim);
182 }
183 
184 template<unsigned int spacedim>
185 void RefElementValues<spacedim>::ref_initialize(Quadrature &_quadrature, unsigned int dim)
186 {
187  //if (dim == 0) return; // avoid unnecessary allocation of dummy 0 dimensional objects
188  if ( _quadrature.dim() == dim )
189  {
190  // precompute element data
191  this->ref_data = this->init_ref_data(_quadrature);
192  }
193  else if ( _quadrature.dim() + 1 == dim )
194  {
195  // precompute side data
196  for (unsigned int sid = 0; sid < this->n_sides_; sid++)
197  {
198  Quadrature side_quad(dim);
199  // transform the side quadrature points to the cell quadrature points
200  switch (dim)
201  {
202  case 1:
203  side_quad = _quadrature.make_from_side<1>(sid);
204  break;
205  case 2:
206  side_quad = _quadrature.make_from_side<2>(sid);
207  break;
208  case 3:
209  side_quad = _quadrature.make_from_side<3>(sid);
210  break;
211  default:
212  ASSERT_PERMANENT(false)(dim).error("Unsupported dimension.\n");
213  break;
214  }
215  this->side_ref_data[sid] = this->init_ref_data(side_quad);
216  }
217  }
218 }
219 
220 
221 template<unsigned int spacedim>
223 {
224  if (ref_data) delete ref_data;
225 
226  for (unsigned int sid=0; sid<n_sides_; sid++)
227  delete side_ref_data[sid];
228 }
229 
230 
231 template<unsigned int spacedim>
233 {
234  ASSERT_EQ( this->dim_, cell.dim() );
235  data.cell = cell;
236 
237  // calculate Jacobian of mapping, JxW, inverse Jacobian
238  switch (this->dim_)
239  {
240  case 0:
241  if (cell.is_valid() && data.update_flags & update_quadrature_points)
242  data.points.set(0) = Armor::vec<spacedim>( MappingP1<0, spacedim>::element_map(cell) );
243  break;
244  case 1:
245  fill_data<1>();
246  break;
247  case 2:
248  fill_data<2>();
249  break;
250  case 3:
251  fill_data<3>();
252  break;
253  default:
254  ASSERT_PERMANENT(false)(this->dim_).error("Unsupported dimension.\n");
255  break;
256  }
257 }
258 
259 
260 template<unsigned int spacedim>
261 void ElementValues<spacedim>::reinit(const Side & cell_side)
262 {
263  ASSERT_EQ( this->dim_, cell_side.dim()+1 );
264  data.side = cell_side;
265 
266  // calculate Jacobian of mapping, JxW, inverse Jacobian, normal vector(s)
267  switch (this->dim_)
268  {
269  case 1:
270  fill_data<1>();
271  fill_side_data<1>();
272  break;
273  case 2:
274  fill_data<2>();
275  fill_side_data<2>();
276  break;
277  case 3:
278  fill_data<3>();
279  fill_side_data<3>();
280  break;
281  default:
282  ASSERT_PERMANENT(false)(this->dim_).error("Unsupported dimension.\n");
283  break;
284  }
285 }
286 
287 
288 template<unsigned int spacedim>
289 template<unsigned int dim>
291 {
292  typename MappingP1<dim,spacedim>::ElementMap coords;
293 
294  if ((data.update_flags & update_jacobians) |
295  (data.update_flags & update_volume_elements) |
296  (data.update_flags & update_JxW_values) |
297  (data.update_flags & update_inverse_jacobians) |
298  (data.update_flags & update_normal_vectors) |
299  (data.update_flags & update_quadrature_points))
300  {
301  if (cell().is_valid())
302  coords = MappingP1<dim,spacedim>::element_map(cell());
303  else
304  coords = MappingP1<dim,spacedim>::element_map(side().element());
305  /*
306  * Coords predpocitat do pomocne struktury.behem add_patch_xyz_points.
307  * ? Kde predpocitavat, metody bulk_points apod. jsou v sablonovanych assemblacich
308  * ? Uz nyni je problem se zavislostmi, kdy fem zavisi na kodu z coupling_lib
309  * Nasledne rozkopirovat do tabulek pro dim a bulk/side.
310  * K tomu nejspis potrebujeme zavest ElQ<Vector> coords pro jednotlive dimenze.
311  * Z coords uz nejspis pujde dopocitat vsechny ostatni veliciny pomoci vektorovych operaci.
312  * Z ElementValues potrebujeme sadu referencnich dat ref_data a side_ref_data pro kazdou dimenzi
313  */
314  }
315 
316  // calculation of Jacobian dependent data
317  if ((data.update_flags & update_jacobians) |
318  (data.update_flags & update_volume_elements) |
319  (data.update_flags & update_JxW_values) |
320  (data.update_flags & update_inverse_jacobians) |
321  (data.update_flags & update_normal_vectors))
322  {
323  arma::mat::fixed<spacedim,dim> jac = MappingP1<dim,spacedim>::jacobian(coords);
324 
325  // update Jacobians
326  if (data.update_flags & update_jacobians)
327  for (unsigned int i=0; i<this->n_points_; i++)
328  data.jacobians.set(i) = Armor::mat<spacedim,dim>( jac );
329 
330  // calculation of determinant dependent data
331  if ((data.update_flags & update_volume_elements) |
332  (data.update_flags & update_JxW_values))
333  {
334  double det = fabs(fe_tools::determinant(jac));
335 
336  // update determinants
337  if (data.update_flags & update_volume_elements)
338  for (unsigned int i=0; i<this->n_points_; i++)
339  data.determinants[i] = det;
340 
341  // update JxW values
342  if (data.update_flags & update_JxW_values)
343  for (unsigned int i=0; i<this->n_points_; i++)
344  data.JxW_values[i] = det*this->ref_data->weights[i];
345  }
346 
347  // update inverse Jacobians
348  if (data.update_flags & update_inverse_jacobians)
349  {
350  arma::mat::fixed<dim,spacedim> ijac;
351  ijac = fe_tools::inverse(jac);
352 // if (dim==spacedim)
353 // {
354 // ijac = inv(jac);
355 // }
356 // else
357 // {
358 // ijac = pinv(jac);
359 // }
360  for (unsigned int i=0; i<this->n_points_; i++)
361  data.inverse_jacobians.set(i) = Armor::mat<dim,spacedim>( ijac );
362  }
363  }
364 
365  // quadrature points in the actual cell coordinate system
366  if (cell().is_valid() && data.update_flags & update_quadrature_points)
367  {
368  for (unsigned int i=0; i<this->n_points_; i++)
369  data.points.set(i) = Armor::vec<spacedim>( coords*this->ref_data->bar_coords[i] );
370  }
371 }
372 
373 
374 template<unsigned int spacedim>
375 template<unsigned int dim>
377 {
378  const unsigned int side_idx = side().side_idx();
379 
380  // calculation of normal vectors to the side
381  if (data.update_flags & update_normal_vectors)
382  {
383  arma::vec::fixed<spacedim> n_cell;
384  n_cell = trans(data.inverse_jacobians.template mat<dim,spacedim>(0))*RefElement<dim>::normal_vector(side_idx);
385  n_cell = n_cell/norm(n_cell,2);
386  for (unsigned int i=0; i<this->n_points_; i++)
387  data.normal_vectors.set(i) = Armor::vec<spacedim>( n_cell );
388  }
389 
390  // Quadrature points in the actual cell coordinate system.
391  // The points location can vary from side to side.
392  if (data.update_flags & update_quadrature_points)
393  {
395  for (unsigned int i=0; i<this->n_points_; i++)
396  data.points.set(i) = Armor::vec<spacedim>( coords*this->side_ref_data[side_idx]->bar_coords[i] );
397  }
398 
399  if (data.update_flags & update_side_JxW_values)
400  {
401  double side_det;
402  if (dim <= 1)
403  {
404  side_det = 1;
405  }
406  else
407  {
408  arma::mat::fixed<spacedim,dim> side_coords;
409  arma::mat::fixed<spacedim, MatrixSizes<dim>::dim_minus_one > side_jac; // some compilers complain for case dim==0
410 
411  // calculation of side Jacobian
412  for (unsigned int n=0; n<dim; n++)
413  for (unsigned int c=0; c<spacedim; c++)
414  side_coords(c,n) = (*data.side.node(n))[c];
415  side_jac = MappingP1<MatrixSizes<dim>::dim_minus_one,spacedim>::jacobian(side_coords);
416 
417  // calculation of JxW
418  side_det = fabs(fe_tools::determinant(side_jac));
419  }
420  for (unsigned int i=0; i<this->n_points_; i++)
421  data.side_JxW_values[i] = side_det*this->side_ref_data[side_idx]->weights[i];
422  }
423 }
424 
425 
426 
427 
428 
429 template class RefElementValues<3>;
430 template class ElementValues<3>;
431 
#define ASSERT(expr)
Definition: asserts.hh:351
#define ASSERT_PERMANENT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
Definition: asserts.hh:348
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:333
ArmaVec< Type, nr > vec(uint mat_index) const
Definition: armor.hh:869
unsigned int size() const
Definition: armor.hh:776
bool is_valid() const
Definition: accessors.hh:184
unsigned int dim() const
Definition: accessors.hh:188
std::vector< double > determinants
Determinants of Jacobians at quadrature points.
ElementAccessor< spacedim > cell
Iterator to last updated cell.
ElementData(unsigned int size, UpdateFlags flags, unsigned int dim)
Resize the data arrays.
void print()
Print calculated data.
std::vector< double > JxW_values
Transformed quadrature weights.
Armor::array normal_vectors
Normal vectors (spacedim) to the element at the quadrature points lying on a side.
std::vector< double > side_JxW_values
JxW values for sides.
Side side
Iterator to last updated cell side.
Class for computation of data on cell and side.
ElementValues(Quadrature &_quadrature, UpdateFlags _flags, unsigned int dim)
Constructor.
void fill_side_data()
Calculates the mapping data on a side of a cell.
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
UpdateFlags update_each(UpdateFlags flags)
Determine quantities to be recomputed on each cell.
void fill_data()
Compute data from reference cell and using MappingP1.
Affine mapping between reference and actual cell.
Definition: mapping_p1.hh:59
static arma::mat::fixed< spacedim, dim > jacobian(const ElementMap &coords)
Definition: mapping_p1.cc:58
arma::mat::fixed< spacedim, dim+1 > ElementMap
Definition: mapping_p1.hh:64
static ElementMap element_map(ElementAccessor< 3 > elm)
Definition: mapping_p1.cc:48
static UpdateFlags update_each(UpdateFlags flags)
Determines which additional quantities have to be computed.
Definition: mapping_p1.cc:30
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
Quadrature make_from_side(unsigned int sid) const
Definition: quadrature.cc:47
unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:86
unsigned int dim() const
Definition: quadrature.hh:72
Armor::ArmaVec< double, point_dim > point(unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:151
double weight(unsigned int i) const
Returns the ith weight.
Definition: quadrature.hh:103
Structure for storing the precomputed element data.
RefElementData(unsigned int np)
Resize vectors to size np.
std::vector< arma::vec > bar_coords
Barycentric coordinates of quadrature points.
std::vector< double > weights
Quadrature weights.
void ref_initialize(Quadrature &_quadrature, unsigned int dim)
Initialize ref_data or side_ref_data.
RefElementData * init_ref_data(const Quadrature &q)
Precompute data on reference element.
~RefElementValues()
Correct deallocation of objects created by 'initialize' methods.
RefElementValues(Quadrature &_quadrature, unsigned int dim)
static BaryPoint local_to_bary(const LocalPoint &lp)
Converts from local to barycentric coordinates.
Definition: ref_element.cc:187
static LocalPoint normal_vector(unsigned int sid)
unsigned int elem_idx() const
Returns index of element in Mesh::element_vec_.
Definition: accessors.hh:448
unsigned int side_idx() const
Returns local index of the side on the element.
Definition: accessors.hh:444
bool is_valid() const
Returns true if the side has assigned element.
Definition: accessors.hh:452
unsigned int dim() const
Returns dimension of the side, that is dimension of the element minus one.
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell.
ArmaMat< double, N, M > mat
Definition: armor.hh:936
ArmaVec< double, N > vec
Definition: armor.hh:933
Eigen::Matrix< ArrayDbl, spacedim, dim > jacobian(const Eigen::Matrix< ArrayDbl, spacedim, dim+1 > &coords)
Definition: eigen_tools.hh:228
arma::mat::fixed< n, m > inverse(const arma::mat::fixed< m, n > &A)
Calculates inverse of rectangular matrix or pseudoinverse of non-rectangular matrix.
Definition: fem_tools.hh:149
double determinant(const T &M)
Calculates determinant of a rectangular matrix.
void printf(BasicWriter< Char > &w, BasicCStringRef< Char > format, ArgList args)
Definition: printf.h:444
Basic definitions of numerical quadrature rules.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
Definition: update_flags.hh:68
@ update_volume_elements
Determinant of the Jacobian.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_jacobians
Volume element.
@ update_side_JxW_values
Transformed quadrature weight for cell sides.
@ update_inverse_jacobians
Volume element.
@ update_quadrature_points
Transformed quadrature points.