Flow123d  release_3.0.0-1165-gffc6a5f
fe_values.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 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 FEValues calculates finite element data on the actual
16  * cells such as shape function values, gradients, Jacobian of
17  * the mapping from the reference cell etc.
18  * @author Jan Stebel
19  */
20 
21 #include "fem/mapping.hh"
22 #include "quadrature/quadrature.hh"
23 #include "fem/finite_element.hh"
24 #include "fem/fe_values.hh"
25 #include "fem/fe_system.hh"
26 
27 
28 
29 using namespace arma;
30 using namespace std;
31 
32 
33 
34 
35 
36 
37 
38 
39 FEInternalData::FEInternalData(unsigned int np, unsigned int nd)
40  : n_points(np),
41  n_dofs(nd)
42 {
43  ref_shape_values.resize(np, vector<arma::vec>(nd));
44  ref_shape_grads.resize(np, vector<arma::mat>(nd));
45 }
46 
47 
48 
49 template<unsigned int dim, unsigned int spacedim>
50 void FEValuesData<dim,spacedim>::allocate(unsigned int size, UpdateFlags flags, unsigned int n_comp)
51 {
52  update_flags = flags;
53 
54  // resize the arrays of computed quantities
55  if (update_flags & update_jacobians)
56  jacobians.resize(size);
57 
58  if (update_flags & update_volume_elements)
59  determinants.resize(size);
60 
61  if ((update_flags & update_JxW_values) |
62  (update_flags & update_side_JxW_values))
63  JxW_values.resize(size);
64 
65  if (update_flags & update_inverse_jacobians)
66  inverse_jacobians.resize(size);
67 
68  if (update_flags & update_values)
69  {
70  shape_values.resize(size, vector<double>(n_comp));
71  }
72 
73  if (update_flags & update_gradients)
74  {
75  shape_gradients.resize(size, vector<arma::vec::fixed<spacedim> >(n_comp));
76  }
77 
78  if (update_flags & update_quadrature_points)
79  points.resize(size);
80 
81  if (update_flags & update_normal_vectors)
82  normal_vectors.resize(size);
83 }
84 
85 
86 template<unsigned int dim, unsigned int spacedim>
88 {
89  scalars.clear();
90  vectors.clear();
91  tensors.clear();
92  switch (fv.get_fe()->type_) {
93  case FEType::FEScalar:
94  scalars.push_back(FEValuesViews::Scalar<dim,spacedim>(fv, 0));
95  break;
96  case FEType::FEVector:
99  vectors.push_back(FEValuesViews::Vector<dim,spacedim>(fv, 0));
100  break;
101  case FEType::FETensor:
102  tensors.push_back(FEValuesViews::Tensor<dim,spacedim>(fv, 0));
103  break;
105  FESystem<dim> *fe_sys = dynamic_cast<FESystem<dim>*>(fv.get_fe());
106  ASSERT_DBG(fe_sys != nullptr).error("Mixed system must be represented by FESystem.");
107 
108  // Loop through sub-elements and add views according to their types.
109  // Note that the component index is calculated using fe->n_space_components(),
110  // not fe->n_components()!
111  unsigned int comp_offset = 0;
112  for (auto fe : fe_sys->fe())
113  {
114  switch (fe->type_)
115  {
116  case FEType::FEScalar:
117  scalars.push_back(FEValuesViews::Scalar<dim,spacedim>(fv,comp_offset));
118  break;
119  case FEType::FEVector:
122  vectors.push_back(FEValuesViews::Vector<dim,spacedim>(fv,comp_offset));
123  break;
124  case FEType::FETensor:
125  tensors.push_back(FEValuesViews::Tensor<dim,spacedim>(fv,comp_offset));
126  break;
127  default:
128  ASSERT_DBG(false).error("Not implemented.");
129  break;
130  }
131 
132  comp_offset += fe->n_space_components(spacedim);
133  }
134  break;
135  }
136 }
137 
138 
139 
140 template<unsigned int dim,unsigned int spacedim>
142 : mapping(NULL), n_points_(0), fe(NULL), mapping_data(NULL), fe_data(NULL)
143 {
144 }
145 
146 
147 
148 template<unsigned int dim,unsigned int spacedim>
150  if (mapping_data) delete mapping_data;
151  if (fe_data) delete fe_data;
152 }
153 
154 
155 
156 template<unsigned int dim, unsigned int spacedim>
158  unsigned int n_points,
159  FiniteElement<dim> & _fe,
160  UpdateFlags _flags)
161 {
162  // For FEVector and FETensor check number of components.
163  // This cannot be done in FiniteElement since it does not know spacedim.
164  if (_fe.type_ == FEVector) {
165  ASSERT_DBG(_fe.n_components() == spacedim).error("FEVector must have spacedim components.");
166  } else if (_fe.type_ == FETensor) {
167  ASSERT_DBG(_fe.n_components() == spacedim*spacedim).error("FETensor must have spacedim*spacedim components.");
168  }
169 
170  mapping = &_mapping;
172  fe = &_fe;
173  n_components_ = fe->n_space_components(spacedim);
174 
175  // add flags required by the finite element or mapping
176  data.allocate(n_points_, update_each(_flags), fe->n_dofs()*n_components_);
177 
178  views_cache_.initialize(*this);
179 }
180 
181 
182 
183 template<unsigned int dim, unsigned int spacedim>
185 {
186  ASSERT_DBG( q->dim() == dim );
187  FEInternalData *data = new FEInternalData(q->size(), fe->n_dofs());
188 
189  arma::mat shape_values(fe->n_dofs(), fe->n_components());
190  for (unsigned int i=0; i<q->size(); i++)
191  {
192  for (unsigned int j=0; j<fe->n_dofs(); j++)
193  {
194  for (unsigned int c=0; c<fe->n_components(); c++)
195  shape_values(j,c) = fe->shape_value(j, q->point<dim>(i).arma(), c);
196 
197  data->ref_shape_values[i][j] = trans(shape_values.row(j));
198  }
199  }
200 
201  arma::mat grad(dim, fe->n_components());
202  for (unsigned int i=0; i<q->size(); i++)
203  {
204  for (unsigned int j=0; j<fe->n_dofs(); j++)
205  {
206  grad.zeros();
207  for (unsigned int c=0; c<fe->n_components(); c++)
208  grad.col(c) += fe->shape_grad(j, q->point<dim>(i).arma(), c);
209 
210  data->ref_shape_grads[i][j] = grad;
211  }
212  }
213 
214  return data;
215 }
216 
217 
218 
219 template<unsigned int dim, unsigned int spacedim>
221 {
222  UpdateFlags f = flags | fe->update_each(flags);
223  f |= mapping->update_each(f);
224  return f;
225 }
226 
227 
228 template<unsigned int dim, unsigned int spacedim>
229 double FEValuesBase<dim,spacedim>::shape_value(const unsigned int function_no, const unsigned int point_no)
230 {
231  ASSERT_LT_DBG(function_no, fe->n_dofs());
232  ASSERT_LT_DBG(point_no, n_points_);
233  return data.shape_values[point_no][function_no];
234 }
235 
236 
237 template<unsigned int dim, unsigned int spacedim>
238 arma::vec::fixed<spacedim> FEValuesBase<dim,spacedim>::shape_grad(const unsigned int function_no, const unsigned int point_no)
239 {
240  ASSERT_LT_DBG(function_no, fe->n_dofs());
241  ASSERT_LT_DBG(point_no, n_points_);
242  return data.shape_gradients[point_no][function_no];
243 }
244 
245 
246 template<unsigned int dim, unsigned int spacedim>
247 double FEValuesBase<dim,spacedim>::shape_value_component(const unsigned int function_no,
248  const unsigned int point_no,
249  const unsigned int comp) const
250 {
251  ASSERT_LT_DBG(function_no, fe->n_dofs());
252  ASSERT_LT_DBG(point_no, n_points_);
254  return data.shape_values[point_no][function_no*n_components_+comp];
255 }
256 
257 
258 template<unsigned int dim, unsigned int spacedim>
259 arma::vec::fixed<spacedim> FEValuesBase<dim,spacedim>::shape_grad_component(const unsigned int function_no,
260  const unsigned int point_no,
261  const unsigned int comp) const
262 {
263  ASSERT_LT_DBG(function_no, fe->n_dofs());
264  ASSERT_LT_DBG(point_no, n_points_);
266  return data.shape_gradients[point_no][function_no*n_components_+comp];
267 }
268 
269 
270 template<unsigned int dim, unsigned int spacedim>
272 {
273  ASSERT_DBG(fe->type_ == FEScalar);
274 
275  // shape values
276  if (data.update_flags & update_values)
277  for (unsigned int i = 0; i < fe_data.n_points; i++)
278  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
279  data.shape_values[i][j] = fe_data.ref_shape_values[i][j][0];
280 
281  // shape gradients
282  if (data.update_flags & update_gradients)
283  for (unsigned int i = 0; i < fe_data.n_points; i++)
284  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
285  data.shape_gradients[i][j] = trans(data.inverse_jacobians[i]) * fe_data.ref_shape_grads[i][j];
286 }
287 
288 
289 template<unsigned int dim, unsigned int spacedim>
291 {
292  ASSERT_DBG(fe->type_ == FEVector);
293 
294  // shape values
295  if (data.update_flags & update_values)
296  {
297  for (unsigned int i = 0; i < fe_data.n_points; i++)
298  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
299  {
300  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
301  for (unsigned int c=0; c<spacedim; c++)
302  data.shape_values[i][j*spacedim+c] = fv_vec[c];
303  }
304  }
305 
306  // shape gradients
307  if (data.update_flags & update_gradients)
308  {
309  for (unsigned int i = 0; i < fe_data.n_points; i++)
310  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
311  {
312  arma::mat grads = trans(data.inverse_jacobians[i]) * fe_data.ref_shape_grads[i][j];
313  for (unsigned int c=0; c<spacedim; c++)
314  data.shape_gradients[i][j*spacedim+c] = grads.col(c);
315  }
316  }
317 }
318 
319 
320 template<unsigned int dim, unsigned int spacedim>
322 {
324 
325  // shape values
326  if (data.update_flags & update_values)
327  {
328  for (unsigned int i = 0; i < fe_data.n_points; i++)
329  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
330  {
331  arma::vec fv_vec = data.jacobians[i] * fe_data.ref_shape_values[i][j];
332  for (unsigned int c=0; c<spacedim; c++)
333  data.shape_values[i][j*spacedim+c] = fv_vec[c];
334  }
335  }
336 
337  // shape gradients
338  if (data.update_flags & update_gradients)
339  {
340  for (unsigned int i = 0; i < fe_data.n_points; i++)
341  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
342  {
343  arma::mat grads = trans(data.inverse_jacobians[i]) * fe_data.ref_shape_grads[i][j] * trans(data.jacobians[i]);
344  for (unsigned int c=0; c<spacedim; c++)
345  data.shape_gradients[i][j*spacedim+c] = grads.col(c);
346  }
347  }
348 }
349 
350 
351 template<unsigned int dim, unsigned int spacedim>
353 {
354  ASSERT_DBG(fe->type_ == FEVectorPiola);
355 
356  // shape values
357  if (data.update_flags & update_values)
358  {
359  for (unsigned int i = 0; i < fe_data.n_points; i++)
360  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
361  {
362  arma::vec fv_vec = data.jacobians[i]*fe_data.ref_shape_values[i][j]/data.determinants[i];
363  for (unsigned int c=0; c<spacedim; c++)
364  data.shape_values[i][j*spacedim+c] = fv_vec(c);
365  }
366  }
367 
368  // shape gradients
369  if (data.update_flags & update_gradients)
370  {
371  for (unsigned int i = 0; i < fe_data.n_points; i++)
372  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
373  {
374  arma::mat grads = trans(data.inverse_jacobians[i]) * fe_data.ref_shape_grads[i][j] * trans(data.jacobians[i])
375  / data.determinants[i];
376  for (unsigned int c=0; c<spacedim; c++)
377  data.shape_gradients[i][j*spacedim+c] = grads.col(c);
378  }
379  }
380 }
381 
382 
383 template<unsigned int dim, unsigned int spacedim>
385 {
386  ASSERT_DBG(fe->type_ == FETensor);
387 
388  // shape values
389  if (data.update_flags & update_values)
390  {
391  for (unsigned int i = 0; i < fe_data.n_points; i++)
392  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
393  {
394  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
395  for (unsigned int c=0; c<spacedim*spacedim; c++)
396  data.shape_values[i][j*spacedim*spacedim+c] = fv_vec[c];
397  }
398  }
399 
400  // shape gradients
401  if (data.update_flags & update_gradients)
402  {
403  for (unsigned int i = 0; i < fe_data.n_points; i++)
404  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
405  {
406  arma::mat grads = trans(data.inverse_jacobians[i]) * fe_data.ref_shape_grads[i][j];
407  for (unsigned int c=0; c<spacedim*spacedim; c++)
408  data.shape_gradients[i][j*spacedim*spacedim+c] = grads.col(c);
409  }
410  }
411 }
412 
413 
414 template<unsigned int dim, unsigned int spacedim>
416 {
417  ASSERT_DBG(fe->type_ == FEMixedSystem);
418 
419  // for mixed system we first fill data in sub-elements
420  FESystem<dim> *fe_sys = dynamic_cast<FESystem<dim>*>(fe);
421  ASSERT_DBG(fe_sys != nullptr).error("Mixed system must be represented by FESystem.");
422  for (unsigned int f=0; f<fe_sys->fe().size(); f++)
423  {
424  // fill fe_values for base FE
425  fe_values_vec[f]->data.jacobians = data.jacobians;
426  fe_values_vec[f]->data.inverse_jacobians = data.inverse_jacobians;
427  fe_values_vec[f]->data.determinants = data.determinants;
428  fe_values_vec[f]->fill_data(*fe_values_vec[f]->fe_data);
429  }
430 
431  unsigned int n_space_components = fe->n_space_components(spacedim);
432 
433  // shape values
434  if (data.update_flags & update_values)
435  {
436  arma::vec fv_vec;
437  unsigned int comp_offset = 0;
438  unsigned int shape_offset = 0;
439  for (unsigned int f=0; f<fe_sys->fe().size(); f++)
440  {
441  unsigned int n_sub_space_components = fe_sys->fe()[f]->n_space_components(spacedim);
442 
443  // gather fe_values in vectors for FESystem
444  for (unsigned int i=0; i<fe_data.n_points; i++)
445  for (unsigned int n=0; n<fe_sys->fe()[f]->n_dofs(); n++)
446  for (unsigned int c=0; c<n_sub_space_components; c++)
447  data.shape_values[i][shape_offset+n_space_components*n+comp_offset+c] = fe_values_vec[f]->data.shape_values[i][n*n_sub_space_components+c];
448 
449  comp_offset += n_sub_space_components;
450  shape_offset += fe_sys->fe()[f]->n_dofs()*n_space_components;
451  }
452  }
453 
454  // shape gradients
455  if (data.update_flags & update_gradients)
456  {
457  arma::mat grads;
458  unsigned int comp_offset = 0;
459  unsigned int shape_offset = 0;
460  for (unsigned int f=0; f<fe_sys->fe().size(); f++)
461  {
462  unsigned int n_sub_space_components = fe_sys->fe()[f]->n_space_components(spacedim);
463 
464  // gather fe_values in vectors for FESystem
465  for (unsigned int i=0; i<fe_data.n_points; i++)
466  for (unsigned int n=0; n<fe_sys->fe()[f]->n_dofs(); n++)
467  for (unsigned int c=0; c<n_sub_space_components; c++)
468  data.shape_gradients[i][shape_offset+n_space_components*n+comp_offset+c] = fe_values_vec[f]->data.shape_gradients[i][n*n_sub_space_components+c];
469 
470  comp_offset += n_sub_space_components;
471  shape_offset += fe_sys->fe()[f]->n_dofs()*n_space_components;
472  }
473  }
474 
475 }
476 
477 
478 template<unsigned int dim, unsigned int spacedim>
480 {
481  switch (fe->type_) {
482  case FEScalar:
483  fill_scalar_data(fe_data);
484  break;
485  case FEVector:
486  fill_vec_data(fe_data);
487  break;
490  break;
491  case FEVectorPiola:
492  fill_vec_piola_data(fe_data);
493  break;
494  case FETensor:
495  fill_tensor_data(fe_data);
496  break;
497  case FEMixedSystem:
498  fill_system_data(fe_data);
499  break;
500  default:
501  ASSERT(false).error("Not implemented.");
502  }
503 }
504 
505 
506 
507 
508 
509 
510 
511 
512 template<unsigned int dim, unsigned int spacedim>
514  Quadrature &_quadrature,
515  FiniteElement<dim> &_fe,
516  UpdateFlags _flags)
517 : FEValuesBase<dim, spacedim>(),
518  quadrature(&_quadrature)
519 {
520  ASSERT_DBG( _quadrature.dim() == dim );
521  this->allocate(_mapping, _quadrature.size(), _fe, _flags);
522 
523  // precompute the maping data and finite element data
524  this->mapping_data = this->mapping->initialize(*quadrature, this->data.update_flags);
525  this->fe_data = this->init_fe_data(quadrature);
526 
527  // In case of mixed system allocate data for sub-elements.
528  if (this->fe->type_ == FEMixedSystem)
529  {
530  FESystem<dim> *fe = dynamic_cast<FESystem<dim>*>(this->fe);
531  ASSERT_DBG(fe != nullptr).error("Mixed system must be represented by FESystem.");
532 
533  for (auto fe_sub : fe->fe())
534  this->fe_values_vec.push_back(make_shared<FEValues<dim,spacedim> >(_mapping, _quadrature, *fe_sub, _flags));
535  }
536 }
537 
538 
539 
540 template<unsigned int dim,unsigned int spacedim>
542 {
543  OLD_ASSERT_EQUAL( dim, cell->dim() );
544  this->data.present_cell = &cell;
545 
546  // calculate Jacobian of mapping, JxW, inverse Jacobian, normal vector(s)
547  this->mapping->fill_fe_values(cell,
548  *this->quadrature,
549  *this->mapping_data,
550  this->data);
551 
552  this->fill_data(*this->fe_data);
553 }
554 
555 
556 
557 
558 
559 
560 
561 
562 
563 template<unsigned int dim,unsigned int spacedim>
565  Quadrature & _sub_quadrature,
566  FiniteElement<dim> & _fe,
567  const UpdateFlags _flags)
568 : FEValuesBase<dim,spacedim>(),
569  side_quadrature(RefElement<dim>::n_sides, std::vector<Quadrature>(RefElement<dim>::n_side_permutations, Quadrature(dim)))
570 {
571  ASSERT_DBG( _sub_quadrature.dim() + 1 == dim );
572  sub_quadrature = &_sub_quadrature;
573 
574  this->allocate(_mapping, _sub_quadrature.size(), _fe, _flags);
575 
576  for (unsigned int sid = 0; sid < RefElement<dim>::n_sides; sid++)
577  {
578  for (unsigned int pid = 0; pid < RefElement<dim>::n_side_permutations; pid++)
579  {
580  // transform the side quadrature points to the cell quadrature points
581  side_quadrature[sid][pid] = _sub_quadrature.make_from_side<dim>(sid, pid);
582  side_mapping_data[sid][pid] = this->mapping->initialize(side_quadrature[sid][pid], this->data.update_flags);
583  side_fe_data[sid][pid] = this->init_fe_data(&side_quadrature[sid][pid]);
584  }
585  }
586 
587  // In case of mixed system allocate data for sub-elements.
588  if (this->fe->type_ == FEMixedSystem)
589  {
590  FESystem<dim> *fe = dynamic_cast<FESystem<dim>*>(this->fe);
591  ASSERT_DBG(fe != nullptr).error("Mixed system must be represented by FESystem.");
592 
593  for (auto fe_sub : fe->fe())
594  this->fe_values_vec.push_back(make_shared<FESideValues<dim,spacedim> >(_mapping, _sub_quadrature, *fe_sub, _flags));
595  }
596 }
597 
598 
599 
600 template<unsigned int dim,unsigned int spacedim>
602 {
603  for (unsigned int sid=0; sid<RefElement<dim>::n_sides; sid++)
604  {
605  for (unsigned int pid=0; pid<RefElement<dim>::n_side_permutations; pid++)
606  {
607  delete side_mapping_data[sid][pid];
608  delete side_fe_data[sid][pid];
609  }
610  }
611 }
612 
613 
614 template<unsigned int dim,unsigned int spacedim>
616  unsigned int sid)
617 {
618  ASSERT_LT_DBG( sid, cell->n_sides());
619  ASSERT_EQ_DBG(dim, cell->dim());
620  this->data.present_cell = &cell;
621 
622  side_idx_ = sid;
623  side_perm_ = cell->permutation_idx(sid);
625  // calculate Jacobian of mapping, JxW, inverse Jacobian, normal vector(s)
626  this->mapping->fill_fe_side_values(cell,
627  sid,
629  *side_mapping_data[sid][side_perm_],
630  this->data);
631 
632  // calculation of finite element data
633  this->fill_data(*side_fe_data[sid][side_perm_]);
634 }
635 
636 
637 
638 
639 
640 
641 
642 template class FEValuesBase<0,3>;
643 template class FEValuesBase<1,3>;
644 template class FEValuesBase<2,3>;
645 template class FEValuesBase<3,3>;
646 
647 template class FEValues<0,3>;
648 template class FEValues<1,3>;
649 template class FEValues<2,3>;
650 template class FEValues<3,3>;
651 
652 template class FESideValues<1,3>;
653 template class FESideValues<2,3>;
654 template class FESideValues<3,3>;
655 
Shape function values.
Definition: update_flags.hh:87
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Definition: update_flags.hh:67
void fill_vec_contravariant_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for vectorial FE.
Definition: fe_values.cc:321
const std::vector< std::shared_ptr< FiniteElement< dim > > > & fe()
Definition: fe_system.hh:142
Transformed quadrature weight for cell sides.
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:331
Mat< double, N, M > mat
Definition: armor.hh:214
unsigned int dim() const
Definition: quadrature.hh:71
void reinit(ElementAccessor< 3 > &cell, unsigned int sid)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:615
void fill_scalar_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for scalar FE.
Definition: fe_values.cc:271
arma::vec::fixed< spacedim > shape_grad_component(const unsigned int function_no, const unsigned int point_no, const unsigned int comp) const
Return the gradient of the function_no-th shape function at the point_no-th quadrature point...
Definition: fe_values.cc:259
Determinant of the Jacobian.
void fill_data(const FEInternalData &fe_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure...
Definition: fe_values.cc:479
std::vector< std::vector< Quadrature > > side_quadrature
Definition: fe_values.hh:626
MappingInternalData * mapping_data
Precomputed mapping data.
Definition: fe_values.hh:488
Quadrature make_from_side(unsigned int sid, unsigned int pid)
Definition: quadrature.cc:52
void initialize(FEValuesBase &fv)
Definition: fe_values.cc:87
UpdateFlags update_each(UpdateFlags flags)
Determine quantities to be recomputed on each cell.
Definition: fe_values.cc:220
unsigned int side_perm_
Permutation index of current side.
Definition: fe_values.hh:636
void allocate(Mapping< dim, spacedim > &_mapping, unsigned int n_points, FiniteElement< dim > &_fe, UpdateFlags flags)
Allocates space for computed data.
Definition: fe_values.cc:157
virtual ~FEValuesBase()
Definition: fe_values.cc:149
FEValuesBase()
Default constructor.
Definition: fe_values.cc:141
Class FEValues calculates finite element data on the actual cells such as shape function values...
Mat< double, N, 1 > vec
Definition: armor.hh:211
FEValuesData< dim, spacedim > data
Data computed by the mapping and finite element.
Definition: fe_values.hh:498
Class FESystem for compound finite elements.
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:346
FEInternalData * fe_data
Precomputed finite element data.
Definition: fe_values.hh:493
FiniteElement< dim > * get_fe() const
Returns the finite element in use.
Definition: fe_values.hh:427
unsigned int n_points()
Returns the number of quadrature points.
Definition: fe_values.hh:407
void allocate(unsigned int size, UpdateFlags flags, unsigned n_comp)
Resize the data arrays.
Definition: fe_values.cc:50
unsigned int side_idx_
Current side on which FESideValues was recomputed.
Definition: fe_values.hh:633
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
Volume element.
unsigned int n_components_
Number of components of the FE.
Definition: fe_values.hh:504
virtual ~FESideValues()
Destructor.
Definition: fe_values.cc:601
Quadrature * quadrature
The quadrature rule used to calculate integrals.
Definition: fe_values.hh:562
Armor::vec< point_dim > point(const unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:90
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const unsigned int point_no) override
Return the gradient of the function_no-th shape function at the point_no-th quadrature point...
Definition: fe_values.cc:238
unsigned int dim() const
Definition: elements.h:124
Transformed quadrature points.
Abstract class for the mapping between reference and actual cell.
Definition: fe_values.hh:38
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:101
const unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:85
void fill_tensor_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for tensorial FE.
Definition: fe_values.cc:384
unsigned int n_sides() const
Definition: elements.h:135
void fill_system_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for mixed system of FE.
Definition: fe_values.cc:415
void fill_vec_piola_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for Raviart-Thomas FE.
Definition: fe_values.cc:352
unsigned int n_components() const
Returns numer of components of the basis function.
FEValues(Mapping< dim, spacedim > &_mapping, Quadrature &_quadrature, FiniteElement< dim > &_fe, UpdateFlags _flags)
Constructor.
Definition: fe_values.cc:513
Basic definitions of numerical quadrature rules.
Shape function gradients.
Definition: update_flags.hh:95
MappingInternalData * side_mapping_data[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
Definition: fe_values.hh:628
FEInternalData(unsigned int np, unsigned int nd)
Definition: fe_values.cc:39
Class Mapping calculates data related to the mapping of the reference cell to the actual cell...
Normal vectors.
FESideValues(Mapping< dim, spacedim > &_mapping, Quadrature &_sub_quadrature, FiniteElement< dim > &_fe, UpdateFlags flags)
Constructor.
Definition: fe_values.cc:564
void fill_vec_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for vectorial FE.
Definition: fe_values.cc:290
std::vector< std::vector< arma::mat > > ref_shape_grads
Precomputed gradients of basis functions at the quadrature points.
Definition: fe_values.hh:69
unsigned int n_points_
Number of integration points.
Definition: fe_values.hh:478
FEInternalData * side_fe_data[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
Definition: fe_values.hh:630
ViewsCache views_cache_
Auxiliary storage of FEValuesViews accessors.
Definition: fe_values.hh:507
const Quadrature * sub_quadrature
Quadrature for the integration on the element sides.
Definition: fe_values.hh:624
double shape_value(const unsigned int function_no, const unsigned int point_no) override
Return the value of the function_no-th shape function at the point_no-th quadrature point...
Definition: fe_values.cc:229
Mapping< dim, spacedim > * mapping
The mapping from the reference cell to the actual cell.
Definition: fe_values.hh:475
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int comp) const
Return the value of the function_no-th shape function at the point_no-th quadrature point...
Definition: fe_values.cc:247
FiniteElement< dim > * fe
The used finite element.
Definition: fe_values.hh:483
#define ASSERT_DBG(expr)
Definition: asserts.hh:349
Calculates finite element data on the actual cell.
Definition: fe_values.hh:525
unsigned int n_points
Number of quadrature points.
Definition: fe_values.hh:72
FEInternalData * init_fe_data(const Quadrature *q)
Precompute finite element data on reference element.
Definition: fe_values.cc:184
#define OLD_ASSERT_EQUAL(a, b)
Definition: global_defs.h:133
std::vector< std::shared_ptr< FEValuesBase< dim, spacedim > > > fe_values_vec
Vector of FEValues for sub-elements of FESystem.
Definition: fe_values.hh:501
Structure for storing the precomputed finite element data.
Definition: fe_values.hh:47
std::vector< std::vector< arma::vec > > ref_shape_values
Precomputed values of basis functions at the quadrature points.
Definition: fe_values.hh:60
Abstract class for description of finite elements.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Base class for FEValues and FESideValues.
Definition: fe_values.hh:37
unsigned int n_dofs
Number of dofs (shape functions).
Definition: fe_values.hh:75
unsigned int permutation_idx(unsigned int prm_idx) const
Return permutation_idx of given index.
Definition: elements.h:144
void reinit(ElementAccessor< 3 > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:541
Transformed quadrature weights.
FEType type_
Type of FiniteElement.
Calculates finite element data on a side.
Definition: fe_values.hh:582
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:299