Flow123d  release_3.0.0-1111-g6aae175
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  FEInternalData *data = new FEInternalData(q->size(), fe->n_dofs());
187 
188  arma::mat shape_values(fe->n_dofs(), fe->n_components());
189  for (unsigned int i=0; i<q->size(); i++)
190  {
191  for (unsigned int j=0; j<fe->n_dofs(); j++)
192  {
193  for (unsigned int c=0; c<fe->n_components(); c++)
194  shape_values(j,c) = fe->shape_value(j, q->point(i), c);
195 
196  data->ref_shape_values[i][j] = trans(shape_values.row(j));
197  }
198  }
199 
200  arma::mat grad(dim, fe->n_components());
201  for (unsigned int i=0; i<q->size(); i++)
202  {
203  for (unsigned int j=0; j<fe->n_dofs(); j++)
204  {
205  grad.zeros();
206  for (unsigned int c=0; c<fe->n_components(); c++)
207  grad.col(c) += fe->shape_grad(j, q->point(i), c);
208 
209  data->ref_shape_grads[i][j] = grad;
210  }
211  }
212 
213  return data;
214 }
215 
216 
217 
218 template<unsigned int dim, unsigned int spacedim>
220 {
221  UpdateFlags f = flags | fe->update_each(flags);
222  f |= mapping->update_each(f);
223  return f;
224 }
225 
226 
227 template<unsigned int dim, unsigned int spacedim>
228 double FEValuesBase<dim,spacedim>::shape_value(const unsigned int function_no, const unsigned int point_no)
229 {
230  ASSERT_LT_DBG(function_no, fe->n_dofs());
231  ASSERT_LT_DBG(point_no, n_points_);
232  return data.shape_values[point_no][function_no];
233 }
234 
235 
236 template<unsigned int dim, unsigned int spacedim>
237 arma::vec::fixed<spacedim> FEValuesBase<dim,spacedim>::shape_grad(const unsigned int function_no, const unsigned int point_no)
238 {
239  ASSERT_LT_DBG(function_no, fe->n_dofs());
240  ASSERT_LT_DBG(point_no, n_points_);
241  return data.shape_gradients[point_no][function_no];
242 }
243 
244 
245 template<unsigned int dim, unsigned int spacedim>
246 double FEValuesBase<dim,spacedim>::shape_value_component(const unsigned int function_no,
247  const unsigned int point_no,
248  const unsigned int comp) const
249 {
250  ASSERT_LT_DBG(function_no, fe->n_dofs());
251  ASSERT_LT_DBG(point_no, n_points_);
253  return data.shape_values[point_no][function_no*n_components_+comp];
254 }
255 
256 
257 template<unsigned int dim, unsigned int spacedim>
258 arma::vec::fixed<spacedim> FEValuesBase<dim,spacedim>::shape_grad_component(const unsigned int function_no,
259  const unsigned int point_no,
260  const unsigned int comp) const
261 {
262  ASSERT_LT_DBG(function_no, fe->n_dofs());
263  ASSERT_LT_DBG(point_no, n_points_);
265  return data.shape_gradients[point_no][function_no*n_components_+comp];
266 }
267 
268 
269 template<unsigned int dim, unsigned int spacedim>
271 {
272  ASSERT_DBG(fe->type_ == FEScalar);
273 
274  // shape values
275  if (data.update_flags & update_values)
276  for (unsigned int i = 0; i < fe_data.n_points; i++)
277  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
278  data.shape_values[i][j] = fe_data.ref_shape_values[i][j][0];
279 
280  // shape gradients
281  if (data.update_flags & update_gradients)
282  for (unsigned int i = 0; i < fe_data.n_points; i++)
283  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
284  data.shape_gradients[i][j] = trans(data.inverse_jacobians[i]) * fe_data.ref_shape_grads[i][j];
285 }
286 
287 
288 template<unsigned int dim, unsigned int spacedim>
290 {
291  ASSERT_DBG(fe->type_ == FEVector);
292 
293  // shape values
294  if (data.update_flags & update_values)
295  {
296  for (unsigned int i = 0; i < fe_data.n_points; i++)
297  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
298  {
299  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
300  for (unsigned int c=0; c<spacedim; c++)
301  data.shape_values[i][j*spacedim+c] = fv_vec[c];
302  }
303  }
304 
305  // shape gradients
306  if (data.update_flags & update_gradients)
307  {
308  for (unsigned int i = 0; i < fe_data.n_points; i++)
309  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
310  {
311  arma::mat grads = trans(data.inverse_jacobians[i]) * fe_data.ref_shape_grads[i][j];
312  for (unsigned int c=0; c<spacedim; c++)
313  data.shape_gradients[i][j*spacedim+c] = grads.col(c);
314  }
315  }
316 }
317 
318 
319 template<unsigned int dim, unsigned int spacedim>
321 {
323 
324  // shape values
325  if (data.update_flags & update_values)
326  {
327  for (unsigned int i = 0; i < fe_data.n_points; i++)
328  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
329  {
330  arma::vec fv_vec = data.jacobians[i] * fe_data.ref_shape_values[i][j];
331  for (unsigned int c=0; c<spacedim; c++)
332  data.shape_values[i][j*spacedim+c] = fv_vec[c];
333  }
334  }
335 
336  // shape gradients
337  if (data.update_flags & update_gradients)
338  {
339  for (unsigned int i = 0; i < fe_data.n_points; i++)
340  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
341  {
342  arma::mat grads = trans(data.inverse_jacobians[i]) * fe_data.ref_shape_grads[i][j] * trans(data.jacobians[i]);
343  for (unsigned int c=0; c<spacedim; c++)
344  data.shape_gradients[i][j*spacedim+c] = grads.col(c);
345  }
346  }
347 }
348 
349 
350 template<unsigned int dim, unsigned int spacedim>
352 {
353  ASSERT_DBG(fe->type_ == FEVectorPiola);
354 
355  // shape values
356  if (data.update_flags & update_values)
357  {
358  for (unsigned int i = 0; i < fe_data.n_points; i++)
359  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
360  {
361  arma::vec fv_vec = data.jacobians[i]*fe_data.ref_shape_values[i][j]/data.determinants[i];
362  for (unsigned int c=0; c<spacedim; c++)
363  data.shape_values[i][j*spacedim+c] = fv_vec(c);
364  }
365  }
366 
367  // shape gradients
368  if (data.update_flags & update_gradients)
369  {
370  for (unsigned int i = 0; i < fe_data.n_points; i++)
371  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
372  {
373  arma::mat grads = trans(data.inverse_jacobians[i]) * fe_data.ref_shape_grads[i][j] * trans(data.jacobians[i])
374  / data.determinants[i];
375  for (unsigned int c=0; c<spacedim; c++)
376  data.shape_gradients[i][j*spacedim+c] = grads.col(c);
377  }
378  }
379 }
380 
381 
382 template<unsigned int dim, unsigned int spacedim>
384 {
385  ASSERT_DBG(fe->type_ == FETensor);
386 
387  // shape values
388  if (data.update_flags & update_values)
389  {
390  for (unsigned int i = 0; i < fe_data.n_points; i++)
391  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
392  {
393  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
394  for (unsigned int c=0; c<spacedim*spacedim; c++)
395  data.shape_values[i][j*spacedim*spacedim+c] = fv_vec[c];
396  }
397  }
398 
399  // shape gradients
400  if (data.update_flags & update_gradients)
401  {
402  for (unsigned int i = 0; i < fe_data.n_points; i++)
403  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
404  {
405  arma::mat grads = trans(data.inverse_jacobians[i]) * fe_data.ref_shape_grads[i][j];
406  for (unsigned int c=0; c<spacedim*spacedim; c++)
407  data.shape_gradients[i][j*spacedim*spacedim+c] = grads.col(c);
408  }
409  }
410 }
411 
412 
413 template<unsigned int dim, unsigned int spacedim>
415 {
416  ASSERT_DBG(fe->type_ == FEMixedSystem);
417 
418  // for mixed system we first fill data in sub-elements
419  FESystem<dim> *fe_sys = dynamic_cast<FESystem<dim>*>(fe);
420  ASSERT_DBG(fe_sys != nullptr).error("Mixed system must be represented by FESystem.");
421  for (unsigned int f=0; f<fe_sys->fe().size(); f++)
422  {
423  // fill fe_values for base FE
424  fe_values_vec[f]->data.jacobians = data.jacobians;
425  fe_values_vec[f]->data.inverse_jacobians = data.inverse_jacobians;
426  fe_values_vec[f]->data.determinants = data.determinants;
427  fe_values_vec[f]->fill_data(*fe_values_vec[f]->fe_data);
428  }
429 
430  unsigned int n_space_components = fe->n_space_components(spacedim);
431 
432  // shape values
433  if (data.update_flags & update_values)
434  {
435  arma::vec fv_vec;
436  unsigned int comp_offset = 0;
437  unsigned int shape_offset = 0;
438  for (unsigned int f=0; f<fe_sys->fe().size(); f++)
439  {
440  unsigned int n_sub_space_components = fe_sys->fe()[f]->n_space_components(spacedim);
441 
442  // gather fe_values in vectors for FESystem
443  for (unsigned int i=0; i<fe_data.n_points; i++)
444  for (unsigned int n=0; n<fe_sys->fe()[f]->n_dofs(); n++)
445  for (unsigned int c=0; c<n_sub_space_components; c++)
446  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];
447 
448  comp_offset += n_sub_space_components;
449  shape_offset += fe_sys->fe()[f]->n_dofs()*n_space_components;
450  }
451  }
452 
453  // shape gradients
454  if (data.update_flags & update_gradients)
455  {
456  arma::mat grads;
457  unsigned int comp_offset = 0;
458  unsigned int shape_offset = 0;
459  for (unsigned int f=0; f<fe_sys->fe().size(); f++)
460  {
461  unsigned int n_sub_space_components = fe_sys->fe()[f]->n_space_components(spacedim);
462 
463  // gather fe_values in vectors for FESystem
464  for (unsigned int i=0; i<fe_data.n_points; i++)
465  for (unsigned int n=0; n<fe_sys->fe()[f]->n_dofs(); n++)
466  for (unsigned int c=0; c<n_sub_space_components; c++)
467  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];
468 
469  comp_offset += n_sub_space_components;
470  shape_offset += fe_sys->fe()[f]->n_dofs()*n_space_components;
471  }
472  }
473 
474 }
475 
476 
477 template<unsigned int dim, unsigned int spacedim>
479 {
480  switch (fe->type_) {
481  case FEScalar:
482  fill_scalar_data(fe_data);
483  break;
484  case FEVector:
485  fill_vec_data(fe_data);
486  break;
489  break;
490  case FEVectorPiola:
491  fill_vec_piola_data(fe_data);
492  break;
493  case FETensor:
494  fill_tensor_data(fe_data);
495  break;
496  case FEMixedSystem:
497  fill_system_data(fe_data);
498  break;
499  default:
500  ASSERT(false).error("Not implemented.");
501  }
502 }
503 
504 
505 
506 
507 
508 
509 
510 
511 template<unsigned int dim, unsigned int spacedim>
513  Quadrature<dim> &_quadrature,
514  FiniteElement<dim> &_fe,
515  UpdateFlags _flags)
516 : FEValuesBase<dim, spacedim>(),
517  quadrature(&_quadrature)
518 {
519  this->allocate(_mapping, _quadrature.size(), _fe, _flags);
520 
521  // precompute the maping data and finite element data
522  this->mapping_data = this->mapping->initialize(*quadrature, this->data.update_flags);
523  this->fe_data = this->init_fe_data(quadrature);
524 
525  // In case of mixed system allocate data for sub-elements.
526  if (this->fe->type_ == FEMixedSystem)
527  {
528  FESystem<dim> *fe = dynamic_cast<FESystem<dim>*>(this->fe);
529  ASSERT_DBG(fe != nullptr).error("Mixed system must be represented by FESystem.");
530 
531  for (auto fe_sub : fe->fe())
532  this->fe_values_vec.push_back(make_shared<FEValues<dim,spacedim> >(_mapping, _quadrature, *fe_sub, _flags));
533  }
534 }
535 
536 
537 
538 template<unsigned int dim,unsigned int spacedim>
540 {
541  OLD_ASSERT_EQUAL( dim, cell->dim() );
542  this->data.present_cell = &cell;
543 
544  // calculate Jacobian of mapping, JxW, inverse Jacobian, normal vector(s)
545  this->mapping->fill_fe_values(cell,
546  *this->quadrature,
547  *this->mapping_data,
548  this->data);
549 
550  this->fill_data(*this->fe_data);
551 }
552 
553 
554 
555 
556 
557 
558 
559 
560 
561 template<unsigned int dim,unsigned int spacedim>
563  Quadrature<dim-1> & _sub_quadrature,
564  FiniteElement<dim> & _fe,
565  const UpdateFlags _flags)
566 :FEValuesBase<dim,spacedim>()
567 {
568  sub_quadrature = &_sub_quadrature;
569  this->allocate(_mapping, _sub_quadrature.size(), _fe, _flags);
570 
571  for (unsigned int sid = 0; sid < RefElement<dim>::n_sides; sid++)
572  {
573  for (unsigned int pid = 0; pid < RefElement<dim>::n_side_permutations; pid++)
574  {
575  // transform the side quadrature points to the cell quadrature points
576  side_quadrature[sid][pid] = Quadrature<dim>(_sub_quadrature, sid, pid);
577  side_mapping_data[sid][pid] = this->mapping->initialize(side_quadrature[sid][pid], this->data.update_flags);
578  side_fe_data[sid][pid] = this->init_fe_data(&side_quadrature[sid][pid]);
579  }
580  }
581 
582  // In case of mixed system allocate data for sub-elements.
583  if (this->fe->type_ == FEMixedSystem)
584  {
585  FESystem<dim> *fe = dynamic_cast<FESystem<dim>*>(this->fe);
586  ASSERT_DBG(fe != nullptr).error("Mixed system must be represented by FESystem.");
587 
588  for (auto fe_sub : fe->fe())
589  this->fe_values_vec.push_back(make_shared<FESideValues<dim,spacedim> >(_mapping, _sub_quadrature, *fe_sub, _flags));
590  }
591 }
592 
593 
594 
595 template<unsigned int dim,unsigned int spacedim>
597 {
598  for (unsigned int sid=0; sid<RefElement<dim>::n_sides; sid++)
599  {
600  for (unsigned int pid=0; pid<RefElement<dim>::n_side_permutations; pid++)
601  {
602  delete side_mapping_data[sid][pid];
603  delete side_fe_data[sid][pid];
604  }
605  }
606 }
607 
608 
609 template<unsigned int dim,unsigned int spacedim>
611  unsigned int sid)
612 {
613  ASSERT_LT_DBG( sid, cell->n_sides());
614  ASSERT_EQ_DBG(dim, cell->dim());
615  this->data.present_cell = &cell;
616 
617  side_idx_ = sid;
618  side_perm_ = cell->permutation_idx(sid);
620  // calculate Jacobian of mapping, JxW, inverse Jacobian, normal vector(s)
621  this->mapping->fill_fe_side_values(cell,
622  sid,
624  *side_mapping_data[sid][side_perm_],
625  this->data);
626 
627  // calculation of finite element data
628  this->fill_data(*side_fe_data[sid][side_perm_]);
629 }
630 
631 
632 
633 
634 
635 
636 
637 template class FEValuesBase<0,3>;
638 template class FEValuesBase<1,3>;
639 template class FEValuesBase<2,3>;
640 template class FEValuesBase<3,3>;
641 
642 template class FEValues<0,3>;
643 template class FEValues<1,3>;
644 template class FEValues<2,3>;
645 template class FEValues<3,3>;
646 
647 template class FESideValues<1,3>;
648 template class FESideValues<2,3>;
649 template class FESideValues<3,3>;
650 
Shape function values.
Definition: update_flags.hh:87
Quadrature< dim > * quadrature
The quadrature rule used to calculate integrals.
Definition: fe_values.hh:562
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:320
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
void reinit(ElementAccessor< 3 > &cell, unsigned int sid)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:610
void fill_scalar_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for scalar FE.
Definition: fe_values.cc:270
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:258
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:478
MappingInternalData * mapping_data
Precomputed mapping data.
Definition: fe_values.hh:488
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:219
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...
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: fe_values.hh:35
Volume element.
unsigned int n_components_
Number of components of the FE.
Definition: fe_values.hh:504
virtual ~FESideValues()
Destructor.
Definition: fe_values.cc:596
FEValues(Mapping< dim, spacedim > &_mapping, Quadrature< dim > &_quadrature, FiniteElement< dim > &_fe, UpdateFlags _flags)
Constructor.
Definition: fe_values.cc:512
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:237
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
void fill_tensor_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for tensorial FE.
Definition: fe_values.cc:383
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:414
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:351
unsigned int n_components() const
Returns numer of components of the basis function.
Basic definitions of numerical quadrature rules.
const Quadrature< dim-1 > * sub_quadrature
Quadrature for the integration on the element sides.
Definition: fe_values.hh:624
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
Quadrature< dim > side_quadrature[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
Definition: fe_values.hh:626
FEInternalData * init_fe_data(const Quadrature< dim > *q)
Precompute finite element data on reference element.
Definition: fe_values.cc:184
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.
void fill_vec_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for vectorial FE.
Definition: fe_values.cc:289
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
FESideValues(Mapping< dim, spacedim > &_mapping, Quadrature< dim-1 > &_sub_quadrature, FiniteElement< dim > &_fe, UpdateFlags flags)
Constructor.
Definition: fe_values.cc:562
ViewsCache views_cache_
Auxiliary storage of FEValuesViews accessors.
Definition: fe_values.hh:507
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:228
Mapping< dim, spacedim > * mapping
The mapping from the reference cell to the actual cell.
Definition: fe_values.hh:475
const unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:139
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:246
FiniteElement< dim > * fe
The used finite element.
Definition: fe_values.hh:483
#define ASSERT_DBG(expr)
Definition: asserts.hh:349
const arma::vec::fixed< dim > & point(const unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:144
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
#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:539
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