Flow123d  3.9.1-8f37ba8
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 
25 
26 
27 using namespace arma;
28 using namespace std;
29 
30 
31 
32 
33 
34 
35 
36 
38  : n_points(np)
39 {
40  bar_coords.resize(np);
41  weights.resize(np);
42 }
43 
44 
45 template<unsigned int spacedim>
47  UpdateFlags flags,
48  unsigned int dim)
49 : dim_(dim),
50  JxW_values(flags & update_JxW_values ? size : 0),
51  side_JxW_values(flags & update_side_JxW_values ? size : 0),
52  jacobians(spacedim, dim, flags & update_jacobians ? size : 0),
53  determinants(flags & update_volume_elements ? size : 0),
54  inverse_jacobians(dim, spacedim, flags & update_inverse_jacobians ? size : 0),
55  points(spacedim, 1, flags & update_quadrature_points ? size : 0),
56  normal_vectors(spacedim, 1, flags & update_normal_vectors ? size : 0),
57  update_flags(flags)
58 {}
59 
60 
61 template<unsigned int spacedim>
63 {
64  if (cell.is_valid() || side.is_valid())
65  {
66  if (cell.is_valid())
67  printf("cell %d dim %d ", cell.elm_idx(), cell.dim());
68  else if (side.is_valid())
69  printf("cell %d dim %d side %d ", side.elem_idx(), side.dim(), side.side_idx());
70 
71  printf(" det[");
72  for (auto d : determinants) printf("%g ", d); printf("]");
73 
74  printf(" JxW[");
75  for (auto j : JxW_values) printf("%g ", j); printf("]");
76 
77  printf(" side_JxW[");
78  for (auto j : side_JxW_values) printf("%g ", j); printf("]");
79 
80  printf(" nv[");
81  for (unsigned int i=0; i<normal_vectors.size(); i++)
82  {
83  auto n = normal_vectors.vec<spacedim>(i);
84  printf(" [");
85  for (unsigned int c=0; c<spacedim; c++) printf("%g ", n[c]);
86  printf("]");
87  }
88  printf("]\n");
89  }
90 }
91 
92 
93 
94 
95 
96 template<unsigned int spacedim>
98 {
99  ASSERT( q.dim() == dim_ );
100  ASSERT( q.size() == n_points_ );
101  RefElementData *ref_data = new RefElementData(q.size());
102 
103  for (unsigned int i=0; i<q.size(); i++)
104  {
105  switch (q.dim())
106  {
107  case 1:
108  ref_data->bar_coords[i] = RefElement<1>::local_to_bary(q.point<1>(i));
109  break;
110  case 2:
111  ref_data->bar_coords[i] = RefElement<2>::local_to_bary(q.point<2>(i));
112  break;
113  case 3:
114  ref_data->bar_coords[i] = RefElement<3>::local_to_bary(q.point<3>(i));
115  break;
116  default:
117  ASSERT_PERMANENT(false)(q.dim()).error("Unsupported dimension.\n");
118  break;
119  }
120  ref_data->weights[i] = q.weight(i);
121  }
122 
123  return ref_data;
124 }
125 
126 
127 
128 template<unsigned int spacedim>
130 {
131  switch (dim_)
132  {
133  case 0:
134  flags = MappingP1<0,spacedim>::update_each(flags);
135  break;
136  case 1:
137  flags = MappingP1<1,spacedim>::update_each(flags);
138  break;
139  case 2:
140  flags = MappingP1<2,spacedim>::update_each(flags);
141  break;
142  case 3:
143  flags = MappingP1<3,spacedim>::update_each(flags);
144  break;
145  default:
146  ASSERT_PERMANENT(false)(dim_).error("Unsupported dimension.\n");
147  break;
148  }
149  return flags;
150 }
151 
152 
153 
154 template<unsigned int spacedim>
156  Quadrature &_quadrature,
157  UpdateFlags _flags,
158  unsigned int dim)
159 : dim_(dim),
160  n_points_(_quadrature.size()),
161  n_sides_(_quadrature.dim() == dim ? 0 : dim+1),
162  ref_data(nullptr),
163  side_ref_data(n_sides_),
164  data(n_points_, update_each(_flags), dim)
165 {
166  if (dim == 0) return; // avoid unnecessary allocation of dummy 0 dimensional objects
167  if ( _quadrature.dim() == dim )
168  {
169  // precompute element data
170  ref_data = init_ref_data(_quadrature);
171  }
172  else if ( _quadrature.dim() + 1 == dim )
173  {
174  // precompute side data
175  for (unsigned int sid = 0; sid < n_sides_; sid++)
176  {
177  Quadrature side_quad(dim);
178  // transform the side quadrature points to the cell quadrature points
179  switch (dim)
180  {
181  case 1:
182  side_quad = _quadrature.make_from_side<1>(sid);
183  break;
184  case 2:
185  side_quad = _quadrature.make_from_side<2>(sid);
186  break;
187  case 3:
188  side_quad = _quadrature.make_from_side<3>(sid);
189  break;
190  default:
191  ASSERT_PERMANENT(false)(dim).error("Unsupported dimension.\n");
192  break;
193  }
194  side_ref_data[sid] = init_ref_data(side_quad);
195  }
196  }
197 }
198 
199 
200 template<unsigned int spacedim>
202 {
203  if (ref_data) delete ref_data;
204 
205  for (unsigned int sid=0; sid<n_sides_; sid++)
206  delete side_ref_data[sid];
207 }
208 
209 
210 template<unsigned int spacedim>
212 {
213  ASSERT_EQ( dim_, cell.dim() );
214  data.cell = cell;
215 
216  // calculate Jacobian of mapping, JxW, inverse Jacobian
217  switch (dim_)
218  {
219  case 1:
220  fill_data<1>();
221  break;
222  case 2:
223  fill_data<2>();
224  break;
225  case 3:
226  fill_data<3>();
227  break;
228  default:
229  ASSERT_PERMANENT(false)(dim_).error("Unsupported dimension.\n");
230  break;
231  }
232 }
233 
234 
235 template<unsigned int spacedim>
236 void ElementValues<spacedim>::reinit(const Side & cell_side)
237 {
238  ASSERT_EQ( dim_, cell_side.dim()+1 );
239  data.side = cell_side;
240 
241  // calculate Jacobian of mapping, JxW, inverse Jacobian, normal vector(s)
242  switch (dim_)
243  {
244  case 1:
245  fill_data<1>();
246  fill_side_data<1>();
247  break;
248  case 2:
249  fill_data<2>();
250  fill_side_data<2>();
251  break;
252  case 3:
253  fill_data<3>();
254  fill_side_data<3>();
255  break;
256  default:
257  ASSERT_PERMANENT(false)(dim_).error("Unsupported dimension.\n");
258  break;
259  }
260 }
261 
262 
263 template<unsigned int spacedim>
264 template<unsigned int dim>
266 {
267  typename MappingP1<dim,spacedim>::ElementMap coords;
268 
269  if ((data.update_flags & update_jacobians) |
270  (data.update_flags & update_volume_elements) |
271  (data.update_flags & update_JxW_values) |
272  (data.update_flags & update_inverse_jacobians) |
273  (data.update_flags & update_normal_vectors) |
274  (data.update_flags & update_quadrature_points))
275  {
276  if (cell().is_valid())
277  coords = MappingP1<dim,spacedim>::element_map(cell());
278  else
279  coords = MappingP1<dim,spacedim>::element_map(side().element());
280  }
281 
282  // calculation of Jacobian dependent data
283  if ((data.update_flags & update_jacobians) |
284  (data.update_flags & update_volume_elements) |
285  (data.update_flags & update_JxW_values) |
286  (data.update_flags & update_inverse_jacobians) |
287  (data.update_flags & update_normal_vectors))
288  {
289  arma::mat::fixed<spacedim,dim> jac = MappingP1<dim,spacedim>::jacobian(coords);
290 
291  // update Jacobians
292  if (data.update_flags & update_jacobians)
293  for (unsigned int i=0; i<n_points_; i++)
294  data.jacobians.set(i) = Armor::mat<spacedim,dim>( jac );
295 
296  // calculation of determinant dependent data
297  if ((data.update_flags & update_volume_elements) |
298  (data.update_flags & update_JxW_values))
299  {
300  double det = fabs(::determinant(jac));
301 
302  // update determinants
303  if (data.update_flags & update_volume_elements)
304  for (unsigned int i=0; i<n_points_; i++)
305  data.determinants[i] = det;
306 
307  // update JxW values
308  if (data.update_flags & update_JxW_values)
309  for (unsigned int i=0; i<n_points_; i++)
310  data.JxW_values[i] = det*ref_data->weights[i];
311  }
312 
313  // update inverse Jacobians
314  if (data.update_flags & update_inverse_jacobians)
315  {
316  arma::mat::fixed<dim,spacedim> ijac;
317  if (dim==spacedim)
318  {
319  ijac = inv(jac);
320  }
321  else
322  {
323  ijac = pinv(jac);
324  }
325  for (unsigned int i=0; i<n_points_; i++)
326  data.inverse_jacobians.set(i) = Armor::mat<dim,spacedim>( ijac );
327  }
328  }
329 
330  // quadrature points in the actual cell coordinate system
331  if (cell().is_valid() && data.update_flags & update_quadrature_points)
332  {
333  for (unsigned int i=0; i<n_points_; i++)
334  data.points.set(i) = Armor::vec<spacedim>( coords*ref_data->bar_coords[i] );
335  }
336 }
337 
338 
339 template<unsigned int spacedim>
340 template<unsigned int dim>
342 {
343  const unsigned int side_idx = side().side_idx();
344 
345  // calculation of normal vectors to the side
346  if (data.update_flags & update_normal_vectors)
347  {
348  arma::vec::fixed<spacedim> n_cell;
349  n_cell = trans(data.inverse_jacobians.template mat<dim,spacedim>(0))*RefElement<dim>::normal_vector(side_idx);
350  n_cell = n_cell/norm(n_cell,2);
351  for (unsigned int i=0; i<n_points_; i++)
352  data.normal_vectors.set(i) = Armor::vec<spacedim>( n_cell );
353  }
354 
355  // Quadrature points in the actual cell coordinate system.
356  // The points location can vary from side to side.
357  if (data.update_flags & update_quadrature_points)
358  {
360  for (unsigned int i=0; i<n_points_; i++)
361  data.points.set(i) = Armor::vec<spacedim>( coords*side_ref_data[side_idx]->bar_coords[i] );
362  }
363 
364  if (data.update_flags & update_side_JxW_values)
365  {
366  double side_det;
367  if (dim <= 1)
368  {
369  side_det = 1;
370  }
371  else
372  {
373  arma::mat::fixed<spacedim,dim> side_coords;
374  arma::mat::fixed<spacedim, MatrixSizes<dim>::dim_minus_one > side_jac; // some compilers complain for case dim==0
375 
376  // calculation of side Jacobian
377  for (unsigned int n=0; n<dim; n++)
378  for (unsigned int c=0; c<spacedim; c++)
379  side_coords(c,n) = (*data.side.node(n))[c];
380  side_jac = MappingP1<MatrixSizes<dim>::dim_minus_one,spacedim>::jacobian(side_coords);
381 
382  // calculation of JxW
383  side_det = fabs(::determinant(side_jac));
384  }
385  for (unsigned int i=0; i<n_points_; i++)
386  data.side_JxW_values[i] = side_det*side_ref_data[side_idx]->weights[i];
387  }
388 }
389 
390 
391 
392 
393 
394 template class ElementValues<3>;
395 
ElementAccessor::dim
unsigned int dim() const
Definition: accessors.hh:188
RefElement::normal_vector
static LocalPoint normal_vector(unsigned int sid)
ElementValues::update_each
UpdateFlags update_each(UpdateFlags flags)
Determine quantities to be recomputed on each cell.
Definition: element_values.cc:129
MappingP1::ElementMap
arma::mat::fixed< spacedim, dim+1 > ElementMap
Definition: mapping_p1.hh:64
update_volume_elements
@ update_volume_elements
Determinant of the Jacobian.
Definition: update_flags.hh:147
RefElementData::weights
std::vector< double > weights
Quadrature weights.
Definition: element_values.hh:53
Quadrature::point
Armor::ArmaVec< double, point_dim > point(unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:91
Armor::vec
ArmaVec< double, N > vec
Definition: armor.hh:885
RefElementData::RefElementData
RefElementData(unsigned int np)
Resize vectors to size np.
Definition: element_values.cc:37
ElementValues::init_ref_data
RefElementData * init_ref_data(const Quadrature &q)
Precompute data on reference element.
Definition: element_values.cc:97
RefElementData
Structure for storing the precomputed element data.
Definition: element_values.hh:42
ASSERT
#define ASSERT(expr)
Definition: asserts.hh:351
ElementValues::side_ref_data
std::vector< RefElementData * > side_ref_data
Data on reference element (for each side ).
Definition: element_values.hh:297
MappingP1::update_each
static UpdateFlags update_each(UpdateFlags flags)
Determines which additional quantities have to be computed.
Definition: mapping_p1.cc:30
Quadrature::size
unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:86
Side::dim
unsigned int dim() const
Returns dimension of the side, that is dimension of the element minus one.
Definition: accessors_impl.hh:198
ElementValues::ref_data
RefElementData * ref_data
Data on reference element.
Definition: element_values.hh:294
ElementAccessor
Definition: dh_cell_accessor.hh:32
update_quadrature_points
@ update_quadrature_points
Transformed quadrature points.
Definition: update_flags.hh:102
arma
Definition: doxy_dummy_defs.hh:15
update_normal_vectors
@ update_normal_vectors
Normal vectors.
Definition: update_flags.hh:124
ElementValues::fill_data
void fill_data()
Compute data from reference cell and using MappingP1.
Definition: element_values.cc:265
ElementValues::~ElementValues
~ElementValues()
Correct deallocation of objects created by 'initialize' methods.
Definition: element_values.cc:201
update_inverse_jacobians
@ update_inverse_jacobians
Volume element.
Definition: update_flags.hh:141
ElementValues::ElementValues
ElementValues(Quadrature &_quadrature, UpdateFlags _flags, unsigned int dim)
Constructor.
Definition: element_values.cc:155
ElementValues::reinit
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: element_values.cc:211
Armor::mat
ArmaMat< double, N, M > mat
Definition: armor.hh:888
determinant
double determinant(const T &M)
Calculates determinant of a rectangular matrix.
ElementValues::fill_side_data
void fill_side_data()
Calculates the mapping data on a side of a cell.
Definition: element_values.cc:341
ASSERT_PERMANENT
#define ASSERT_PERMANENT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
Definition: asserts.hh:348
Quadrature::weight
double weight(unsigned int i) const
Returns the ith weight.
Definition: quadrature.hh:107
Side
Definition: accessors.hh:390
ASSERT_EQ
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:333
Quadrature::dim
unsigned int dim() const
Definition: quadrature.hh:72
ElementValues
Class for computation of data on cell and side.
Definition: element_values.hh:128
ElementValues::n_sides_
const unsigned int n_sides_
Number of sides in reference cell.
Definition: element_values.hh:291
MappingP1::jacobian
static arma::mat::fixed< spacedim, dim > jacobian(const ElementMap &coords)
Definition: mapping_p1.cc:58
MappingP1::element_map
static ElementMap element_map(ElementAccessor< 3 > elm)
Definition: mapping_p1.cc:48
ElementData::print
void print()
Print calculated data.
Definition: element_values.cc:62
ElementData::ElementData
ElementData(unsigned int size, UpdateFlags flags, unsigned int dim)
Resize the data arrays.
Definition: element_values.cc:46
fmt::printf
void printf(BasicWriter< Char > &w, BasicCStringRef< Char > format, ArgList args)
Definition: printf.h:444
std
Definition: doxy_dummy_defs.hh:5
MappingP1
Affine mapping between reference and actual cell.
Definition: mapping_p1.hh:58
quadrature.hh
Basic definitions of numerical quadrature rules.
update_JxW_values
@ update_JxW_values
Transformed quadrature weights.
Definition: update_flags.hh:114
element_values.hh
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
UpdateFlags
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
Definition: update_flags.hh:67
mapping_p1.hh
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell.
update_jacobians
@ update_jacobians
Volume element.
Definition: update_flags.hh:132
RefElement::local_to_bary
static BaryPoint local_to_bary(const LocalPoint &lp)
Converts from local to barycentric coordinates.
Definition: ref_element.cc:187
update_side_JxW_values
@ update_side_JxW_values
Transformed quadrature weight for cell sides.
Definition: update_flags.hh:153
Quadrature
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
RefElementData::bar_coords
std::vector< arma::vec > bar_coords
Barycentric coordinates of quadrature points.
Definition: element_values.hh:50
Quadrature::make_from_side
Quadrature make_from_side(unsigned int sid) const
Definition: quadrature.cc:47