Flow123d  3.9.0-68a662b83
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_p1.hh"
22 #include "quadrature/quadrature.hh"
23 #include "fem/element_values.hh"
24 #include "fem/finite_element.hh"
25 #include "fem/fe_values.hh"
26 #include "fem/fe_system.hh"
27 #include "fem/fe_values_map.hh"
28 
29 
30 
31 using namespace arma;
32 using namespace std;
33 
34 
35 
36 
37 
38 
39 
40 template<unsigned int spacedim>
41 FEValues<spacedim>::FEInternalData::FEInternalData(unsigned int np, unsigned int nd)
42  : n_points(np),
43  n_dofs(nd)
44 {
45  ref_shape_values.resize(np, vector<arma::vec>(nd));
46  ref_shape_grads.resize(np, vector<arma::mat>(nd));
47 }
48 
49 
50 template<unsigned int spacedim>
52  const std::vector<unsigned int> &dof_indices,
53  unsigned int first_component_idx,
54  unsigned int ncomps)
55  : FEInternalData(fe_system_data.n_points, dof_indices.size())
56 {
57  for (unsigned int ip=0; ip<n_points; ip++)
58  for (unsigned int id=0; id<dof_indices.size(); id++)
59  {
60  ref_shape_values[ip][id] = fe_system_data.ref_shape_values[ip][dof_indices[id]].subvec(first_component_idx, first_component_idx+ncomps-1);
61  ref_shape_grads[ip][id] = fe_system_data.ref_shape_grads[ip][dof_indices[id]].cols(first_component_idx, first_component_idx+ncomps-1);
62  }
63 }
64 
65 
66 
67 template<unsigned int spacedim>
68 template<unsigned int DIM>
70 {
71  scalars.clear();
72  vectors.clear();
73  tensors.clear();
74  switch (fe.type_) {
75  case FEType::FEScalar:
76  scalars.push_back(FEValuesViews::Scalar<spacedim>(fv, 0));
77  break;
78  case FEType::FEVector:
81  vectors.push_back(FEValuesViews::Vector<spacedim>(fv, 0));
82  break;
83  case FEType::FETensor:
84  tensors.push_back(FEValuesViews::Tensor<spacedim>(fv, 0));
85  break;
87  const FESystem<DIM> *fe_sys = dynamic_cast<const FESystem<DIM>*>(&fe);
88  ASSERT(fe_sys != nullptr).error("Mixed system must be represented by FESystem.");
89 
90  // Loop through sub-elements and add views according to their types.
91  // Note that the component index is calculated using fe->n_space_components(),
92  // not fe->n_components()!
93  unsigned int comp_offset = 0;
94  for (auto fe : fe_sys->fe())
95  {
96  switch (fe->type_)
97  {
98  case FEType::FEScalar:
99  scalars.push_back(FEValuesViews::Scalar<spacedim>(fv,comp_offset));
100  break;
101  case FEType::FEVector:
104  vectors.push_back(FEValuesViews::Vector<spacedim>(fv,comp_offset));
105  break;
106  case FEType::FETensor:
107  tensors.push_back(FEValuesViews::Tensor<spacedim>(fv,comp_offset));
108  break;
109  default:
110  ASSERT(false).error("Not implemented.");
111  break;
112  }
113 
114  comp_offset += fe->n_space_components(spacedim);
115  }
116  break;
117  }
118 }
119 
120 
121 
122 template<unsigned int spacedim>
124 : dim_(-1), n_points_(0), n_dofs_(0)
125 {
126 }
127 
128 
129 
130 template<unsigned int spacedim>
132 }
133 
134 
135 template<unsigned int spacedim>
136 template<unsigned int DIM>
138  Quadrature &q,
139  FiniteElement<DIM> &_fe,
140  UpdateFlags _flags)
141 {
142  if (DIM == 0) return; // avoid unnecessary allocation of dummy 0 dimensional objects
143 
144  allocate( q.size(), _fe, _flags);
145  elm_values = std::make_shared<ElementValues<spacedim> >(q, update_flags, DIM);
146 
147  // In case of mixed system allocate data for sub-elements.
148  if (fe_type_ == FEMixedSystem)
149  {
150  FESystem<DIM> *fe = dynamic_cast<FESystem<DIM>*>(&_fe);
151  ASSERT(fe != nullptr).error("Mixed system must be represented by FESystem.");
152 
153  fe_values_vec.resize(fe->fe().size());
154  for (unsigned int f=0; f<fe->fe().size(); f++)
155  fe_values_vec[f].initialize(q, *fe->fe()[f], update_flags);
156  }
157 
158  // precompute finite element data
159  if ( q.dim() == DIM )
160  {
161  fe_data = init_fe_data(_fe, q);
162  }
163  else if ( q.dim() + 1 == DIM )
164  {
165  side_fe_data.resize(RefElement<DIM>::n_sides);
166  for (unsigned int sid = 0; sid < RefElement<DIM>::n_sides; sid++)
167  {
168  side_fe_data[sid] = init_fe_data(_fe, q.make_from_side<DIM>(sid));
169  }
170  }
171  else
172  ASSERT(false)(q.dim())(DIM).error("Dimension mismatch in FEValues::initialize().");
173 }
174 
175 
176 
177 template<unsigned int spacedim>
178 template<unsigned int DIM>
180  unsigned int n_points,
181  FiniteElement<DIM> & _fe,
182  UpdateFlags _flags)
183 {
184  // For FEVector and FETensor check number of components.
185  // This cannot be done in FiniteElement since it does not know spacedim.
186  if (_fe.type_ == FEVector) {
187  ASSERT(_fe.n_components() == spacedim).error("FEVector must have spacedim components.");
188  } else if (_fe.type_ == FETensor) {
189  ASSERT(_fe.n_components() == spacedim*spacedim).error("FETensor must have spacedim*spacedim components.");
190  }
191 
192  fe_sys_dofs_.clear();
193  fe_sys_n_components_.clear();
194  fe_sys_n_space_components_.clear();
195 
196  dim_ = DIM;
197  n_points_ = n_points;
198  n_dofs_ = _fe.n_dofs();
199  n_components_ = _fe.n_space_components(spacedim);
200  fe_type_ = _fe.type_;
201  FESystem<DIM> *fe_sys = dynamic_cast<FESystem<DIM>*>(&_fe);
202  if (fe_sys != nullptr)
203  {
204  for (unsigned int f=0; f<fe_sys->fe().size(); f++)
205  {
206  fe_sys_dofs_.push_back(fe_sys->fe_dofs(f));
207  fe_sys_n_components_.push_back(fe_sys->fe()[f]->n_components());
208  fe_sys_n_space_components_.push_back(fe_sys->fe()[f]->n_space_components(spacedim));
209  }
210  }
211 
212  // add flags required by the finite element or mapping
213  update_flags = _flags | _fe.update_each(_flags);
214  update_flags |= MappingP1<DIM,spacedim>::update_each(update_flags);
215  if (update_flags & update_values)
216  shape_values.resize(n_points_, vector<double>(n_dofs_*n_components_));
217 
218  if (update_flags & update_gradients)
219  shape_gradients.resize(n_points_, vector<arma::vec::fixed<spacedim> >(n_dofs_*n_components_));
220 
221  views_cache_.initialize(*this, _fe);
222 }
223 
224 
225 
226 template<unsigned int spacedim>
227 template<unsigned int DIM>
228 std::shared_ptr<typename FEValues<spacedim>::FEInternalData> FEValues<spacedim>::init_fe_data(const FiniteElement<DIM> &fe, const Quadrature &q)
229 {
230  ASSERT( DIM == dim_ );
231  ASSERT( q.dim() == DIM );
232  std::shared_ptr<FEInternalData> data = std::make_shared<FEInternalData>(q.size(), n_dofs_);
233 
234  arma::mat shape_values(n_dofs_, fe.n_components());
235  for (unsigned int i=0; i<q.size(); i++)
236  {
237  for (unsigned int j=0; j<n_dofs_; j++)
238  {
239  for (unsigned int c=0; c<fe.n_components(); c++)
240  shape_values(j,c) = fe.shape_value(j, q.point<DIM>(i), c);
241 
242  data->ref_shape_values[i][j] = trans(shape_values.row(j));
243  }
244  }
245 
246  arma::mat grad(DIM, fe.n_components());
247  for (unsigned int i=0; i<q.size(); i++)
248  {
249  for (unsigned int j=0; j<n_dofs_; j++)
250  {
251  grad.zeros();
252  for (unsigned int c=0; c<fe.n_components(); c++)
253  grad.col(c) += fe.shape_grad(j, q.point<DIM>(i), c);
254 
255  data->ref_shape_grads[i][j] = grad;
256  }
257  }
258 
259  return data;
260 }
261 
262 
263 /*template<unsigned int spacedim>
264 double FEValues<spacedim>::shape_value_component(const unsigned int function_no,
265  const unsigned int point_no,
266  const unsigned int comp) const
267 {
268  ASSERT_LT(function_no, n_dofs_);
269  ASSERT_LT(point_no, n_points_);
270  ASSERT_LT(comp, n_components_);
271  return shape_values[point_no][function_no*n_components_+comp];
272 }*/
273 
274 
275 template<unsigned int spacedim>
276 arma::vec::fixed<spacedim> FEValues<spacedim>::shape_grad_component(const unsigned int function_no,
277  const unsigned int point_no,
278  const unsigned int comp) const
279 {
280  ASSERT_LT(function_no, n_dofs_);
281  ASSERT_LT(point_no, n_points_);
282  ASSERT_LT(comp, n_components_);
283  return shape_gradients[point_no][function_no*n_components_+comp];
284 }
285 
286 
287 /*template<unsigned int spacedim>
288 void FEValues<spacedim>::fill_scalar_data(const ElementValues<spacedim> &elm_values, const FEInternalData &fe_data)
289 {
290  ASSERT(fe_type_ == FEScalar);
291 
292  // shape values
293  if (update_flags & update_values)
294  for (unsigned int i = 0; i < fe_data.n_points; i++)
295  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
296  shape_values[i][j] = fe_data.ref_shape_values[i][j][0];
297 
298  // shape gradients
299  if (update_flags & update_gradients)
300  for (unsigned int i = 0; i < fe_data.n_points; i++)
301  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
302  shape_gradients[i][j] = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
303 }
304 
305 
306 template<unsigned int spacedim>
307 void FEValues<spacedim>::fill_vec_data(const ElementValues<spacedim> &elm_values,
308  const FEInternalData &fe_data)
309 {
310  ASSERT(fe_type_ == FEVector);
311 
312  // shape values
313  if (update_flags & update_values)
314  {
315  for (unsigned int i = 0; i < fe_data.n_points; i++)
316  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
317  {
318  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
319  for (unsigned int c=0; c<spacedim; c++)
320  shape_values[i][j*spacedim+c] = fv_vec[c];
321  }
322  }
323 
324  // shape gradients
325  if (update_flags & update_gradients)
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::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
331  for (unsigned int c=0; c<spacedim; c++)
332  shape_gradients[i][j*spacedim+c] = grads.col(c);
333  }
334  }
335 }
336 
337 
338 template<unsigned int spacedim>
339 void FEValues<spacedim>::fill_vec_contravariant_data(const ElementValues<spacedim> &elm_values,
340  const FEInternalData &fe_data)
341 {
342  ASSERT(fe_type_ == FEVectorContravariant);
343 
344  // shape values
345  if (update_flags & update_values)
346  {
347  for (unsigned int i = 0; i < fe_data.n_points; i++)
348  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
349  {
350  arma::vec fv_vec = elm_values.jacobian(i) * fe_data.ref_shape_values[i][j];
351  for (unsigned int c=0; c<spacedim; c++)
352  shape_values[i][j*spacedim+c] = fv_vec[c];
353  }
354  }
355 
356  // shape gradients
357  if (update_flags & update_gradients)
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::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j] * trans(elm_values.jacobian(i));
363  for (unsigned int c=0; c<spacedim; c++)
364  shape_gradients[i][j*spacedim+c] = grads.col(c);
365  }
366  }
367 }
368 
369 
370 template<unsigned int spacedim>
371 void FEValues<spacedim>::fill_vec_piola_data(const ElementValues<spacedim> &elm_values,
372  const FEInternalData &fe_data)
373 {
374  ASSERT(fe_type_ == FEVectorPiola);
375 
376  // shape values
377  if (update_flags & update_values)
378  {
379  for (unsigned int i = 0; i < fe_data.n_points; i++)
380  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
381  {
382  arma::vec fv_vec = elm_values.jacobian(i)*fe_data.ref_shape_values[i][j]/elm_values.determinant(i);
383  for (unsigned int c=0; c<spacedim; c++)
384  shape_values[i][j*spacedim+c] = fv_vec(c);
385  }
386  }
387 
388  // shape gradients
389  if (update_flags & update_gradients)
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::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j] * trans(elm_values.jacobian(i))
395  / elm_values.determinant(i);
396  for (unsigned int c=0; c<spacedim; c++)
397  shape_gradients[i][j*spacedim+c] = grads.col(c);
398  }
399  }
400 }
401 
402 
403 template<unsigned int spacedim>
404 void FEValues<spacedim>::fill_tensor_data(const ElementValues<spacedim> &elm_values,
405  const FEInternalData &fe_data)
406 {
407  ASSERT(fe_type_ == FETensor);
408 
409  // shape values
410  if (update_flags & update_values)
411  {
412  for (unsigned int i = 0; i < fe_data.n_points; i++)
413  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
414  {
415  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
416  for (unsigned int c=0; c<spacedim*spacedim; c++)
417  shape_values[i][j*spacedim*spacedim+c] = fv_vec[c];
418  }
419  }
420 
421  // shape gradients
422  if (update_flags & update_gradients)
423  {
424  for (unsigned int i = 0; i < fe_data.n_points; i++)
425  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
426  {
427  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
428  for (unsigned int c=0; c<spacedim*spacedim; c++)
429  shape_gradients[i][j*spacedim*spacedim+c] = grads.col(c);
430  }
431  }
432 }
433 
434 
435 template<unsigned int spacedim>
436 void FEValues<spacedim>::fill_system_data(const ElementValues<spacedim> &elm_values, const FEInternalData &fe_data)
437 {
438  ASSERT(fe_type_ == FEMixedSystem);
439 
440  // for mixed system we first fill data in sub-elements
441  unsigned int comp_offset = 0;
442  for (unsigned int f=0; f<fe_sys_dofs_.size(); f++)
443  {
444  // fill fe_values for base FE
445  FEInternalData vec_fe_data(fe_data, fe_sys_dofs_[f], comp_offset, fe_sys_n_components_[f]);
446  fe_values_vec[f].fill_data(elm_values, vec_fe_data); // fe_values.fe_values_vec
447 
448  comp_offset += fe_sys_n_components_[f];
449  }
450 
451  // shape values
452  if (update_flags & update_values)
453  {
454  arma::vec fv_vec;
455  unsigned int comp_offset = 0;
456  unsigned int shape_offset = 0;
457  for (unsigned int f=0; f<fe_sys_dofs_.size(); f++)
458  {
459  // gather fe_values in vectors for FESystem
460  for (unsigned int i=0; i<fe_data.n_points; i++)
461  for (unsigned int n=0; n<fe_sys_dofs_[f].size(); n++)
462  for (unsigned int c=0; c<fe_sys_n_space_components_[f]; c++)
463  shape_values[i][shape_offset+n_components_*n+comp_offset+c] = fe_values_vec[f].shape_values[i][n*fe_sys_n_space_components_[f]+c];
464 
465  comp_offset += fe_sys_n_space_components_[f];
466  shape_offset += fe_sys_dofs_[f].size()*n_components_;
467  }
468  }
469 
470  // shape gradients
471  if (update_flags & update_gradients)
472  {
473  arma::mat grads;
474  unsigned int comp_offset = 0;
475  unsigned int shape_offset = 0;
476  for (unsigned int f=0; f<fe_sys_dofs_.size(); f++)
477  {
478  // gather fe_values in vectors for FESystem
479  for (unsigned int i=0; i<fe_data.n_points; i++)
480  for (unsigned int n=0; n<fe_sys_dofs_[f].size(); n++)
481  for (unsigned int c=0; c<fe_sys_n_space_components_[f]; c++)
482  shape_gradients[i][shape_offset+n_components_*n+comp_offset+c] = fe_values_vec[f].shape_gradients[i][n*fe_sys_n_space_components_[f]+c];
483 
484  comp_offset += fe_sys_n_space_components_[f];
485  shape_offset += fe_sys_dofs_[f].size()*n_components_;
486  }
487  }
488 
489 }*/
490 
491 
492 template<unsigned int spacedim>
493 void FEValues<spacedim>::fill_data(const ElementValues<spacedim> &elm_values, const FEInternalData &fe_data)
494 {
495  switch (fe_type_) {
496  case FEScalar:
497  this->fill_data_specialized<MapScalar<spacedim>>(elm_values, fe_data);
498  break;
499  case FEVector:
500  this->fill_data_specialized<MapVector<spacedim>>(elm_values, fe_data);
501  break;
503  this->fill_data_specialized<MapContravariant<spacedim>>(elm_values, fe_data);
504  break;
505  case FEVectorPiola:
506  this->fill_data_specialized<MapPiola<spacedim>>(elm_values, fe_data);
507  break;
508  case FETensor:
509  this->fill_data_specialized<MapTensor<spacedim>>(elm_values, fe_data);
510  break;
511  case FEMixedSystem:
512  this->fill_data_specialized<MapSystem<spacedim>>(elm_values, fe_data);
513  break;
514  default:
515  ASSERT_PERMANENT(false).error("Not implemented.");
516  }
517 }
518 
519 
520 
521 template<unsigned int spacedim>
522 template<class MapType>
524  MapType map_type;
525  map_type.fill_values_vec(*this, elm_values, fe_data);
526  if (update_flags & update_values)
527  map_type.update_values(*this, elm_values, fe_data);
528  if (update_flags & update_gradients)
529  map_type.update_gradients(*this, elm_values, fe_data);
530 }
531 
532 
533 
534 
535 template<unsigned int spacedim>
537 {
538  ASSERT_EQ( dim_, cell.dim() );
539 
540  if (!elm_values->cell().is_valid() ||
541  elm_values->cell() != cell)
542  {
543  elm_values->reinit(cell);
544  }
545 
546  fill_data(*elm_values, *fe_data);
547 }
548 
549 
550 template<unsigned int spacedim>
551 void FEValues<spacedim>::reinit(const Side &cell_side)
552 {
553  ASSERT_EQ( dim_, cell_side.dim()+1 );
554 
555  if (!elm_values->side().is_valid() ||
556  elm_values->side() != cell_side)
557  {
558  elm_values->reinit(cell_side);
559  }
560 
561  const LongIdx sid = cell_side.side_idx();
562 
563  // calculation of finite element data
564  fill_data(*elm_values, *(side_fe_data[sid]) );
565 }
566 
567 
568 
570  QGauss::array &quadrature,
572  UpdateFlags flags)
573 {
575  fv[0].initialize(quadrature[0], *fe[0_d], flags);
576  fv[1].initialize(quadrature[1], *fe[1_d], flags);
577  fv[2].initialize(quadrature[2], *fe[2_d], flags);
578  fv[3].initialize(quadrature[3], *fe[3_d], flags);
579  return fv;
580 }
581 
582 
583 
584 // explicit instantiation
595 
596 
597 
598 
599 
600 
601 
602 
603 
604 
605 
606 
607 
608 
609 
610 template class FEValues<3>;
ElementAccessor::dim
unsigned int dim() const
Definition: accessors.hh:188
Quadrature::point
Armor::ArmaVec< double, point_dim > point(unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:91
FEValues::~FEValues
~FEValues()
Correct deallocation of objects created by 'initialize' methods.
Definition: fe_values.cc:131
update_gradients
@ update_gradients
Shape function gradients.
Definition: update_flags.hh:95
MapPiola
Helper class allows update values and gradients of FEValues of FEVectorPiola type.
Definition: fe_values.hh:41
FEValues::n_points
unsigned int n_points() const
Returns the number of quadrature points.
Definition: fe_values.hh:296
FEScalar
@ FEScalar
Definition: finite_element.hh:200
FEValuesViews::Scalar
Definition: fe_values_views.hh:38
RefElement
Definition: ref_element.hh:339
FEValuesViews::Vector
Definition: fe_values_views.hh:75
FiniteElement::shape_grad
arma::vec::fixed< dim > shape_grad(const unsigned int i, const arma::vec::fixed< dim > &p, const unsigned int comp=0) const
Calculates the comp-th component of the gradient of the i-th shape function at the point p on the ref...
Definition: finite_element.cc:127
ASSERT
#define ASSERT(expr)
Definition: asserts.hh:351
FEValues::FEValues
FEValues()
Default constructor with postponed initialization.
Definition: fe_values.cc:123
FEValuesViews::Tensor
Definition: fe_values_views.hh:126
MappingP1::update_each
static UpdateFlags update_each(UpdateFlags flags)
Determines which additional quantities have to be computed.
Definition: mapping_p1.cc:30
fe_values.hh
Class FEValues calculates finite element data on the actual cells such as shape function values,...
update_values
@ update_values
Shape function values.
Definition: update_flags.hh:87
FEValues::n_dofs_
unsigned int n_dofs_
Number of finite element dofs.
Definition: fe_values.hh:384
Quadrature::size
unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:86
FEValues::initialize
void initialize(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
Definition: fe_values.cc:137
Side::dim
unsigned int dim() const
Returns dimension of the side, that is dimension of the element minus one.
Definition: accessors_impl.hh:198
std::vector< arma::vec >
ElementAccessor
Definition: dh_cell_accessor.hh:32
FEMixedSystem
@ FEMixedSystem
Definition: finite_element.hh:205
arma
Definition: doxy_dummy_defs.hh:15
FEValues::ViewsCache::initialize
void initialize(const FEValues &fv, const FiniteElement< DIM > &fe)
Definition: fe_values.cc:69
fe_system.hh
Class FESystem for compound finite elements.
ASSERT_LT
#define ASSERT_LT(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:301
MapVector
Helper class allows update values and gradients of FEValues of FEVector type.
Definition: fe_values.hh:43
FEValues
Calculates finite element data on the actual cell.
Definition: fe_values.hh:66
FEValues::dim_
unsigned int dim_
Dimension of reference space.
Definition: fe_values.hh:378
FiniteElement::n_components
unsigned int n_components() const
Returns numer of components of the basis function.
Definition: finite_element.hh:292
mixed_fe_values
std::vector< FEValues< 3 > > mixed_fe_values(QGauss::array &quadrature, MixedPtr< FiniteElement > fe, UpdateFlags flags)
Definition: fe_values.cc:569
FEValues::fill_data
void fill_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure.
Definition: fe_values.cc:493
FESystem::fe
const std::vector< std::shared_ptr< FiniteElement< dim > > > & fe() const
Definition: fe_system.hh:147
FiniteElement::update_each
virtual UpdateFlags update_each(UpdateFlags flags)
Decides which additional quantities have to be computed for each cell.
Definition: finite_element.cc:146
FEValues::FEInternalData
Structure for storing the precomputed finite element data.
Definition: fe_values.hh:315
Armor::mat
ArmaMat< double, N, M > mat
Definition: armor.hh:888
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
FiniteElement::type_
FEType type_
Type of FiniteElement.
Definition: finite_element.hh:364
FEValues::FEInternalData::ref_shape_grads
std::vector< std::vector< arma::mat > > ref_shape_grads
Precomputed gradients of basis functions at the quadrature points.
Definition: fe_values.hh:343
finite_element.hh
Abstract class for description of finite elements.
FEVectorContravariant
@ FEVectorContravariant
Definition: finite_element.hh:202
Side
Definition: accessors.hh:390
MapTensor
Helper class allows update values and gradients of FEValues of FETensor type.
Definition: fe_values.hh:44
FETensor
@ FETensor
Definition: finite_element.hh:204
QGauss::array
std::array< QGauss, 4 > array
Definition: quadrature_lib.hh:36
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
FEVector
@ FEVector
Definition: finite_element.hh:201
Side::side_idx
unsigned int side_idx() const
Returns local index of the side on the element.
Definition: accessors.hh:444
MapSystem
Helper class allows update values and gradients of FEValues of FEMixedSystem type.
Definition: fe_values.hh:45
FiniteElement
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Definition: discrete_space.hh:26
FiniteElement::n_dofs
unsigned int n_dofs() const
Returns the number of degrees of freedom needed by the finite element.
Definition: finite_element.hh:262
ElementValues
Class for computation of data on cell and side.
Definition: element_values.hh:128
FEValues::FEInternalData::ref_shape_values
std::vector< std::vector< arma::vec > > ref_shape_values
Precomputed values of basis functions at the quadrature points.
Definition: fe_values.hh:334
FEValues::n_points_
unsigned int n_points_
Number of integration points.
Definition: fe_values.hh:381
LongIdx
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: index_types.hh:24
MapScalar
Helper class allows update values and gradients of FEValues of FEScalar type.
Definition: fe_values.hh:40
FESystem::fe_dofs
std::vector< unsigned int > fe_dofs(unsigned int fe_index)
Return dof indices belonging to given sub-FE.
Definition: fe_system.cc:267
std
Definition: doxy_dummy_defs.hh:5
FEValues::allocate
void allocate(unsigned int n_points, FiniteElement< DIM > &_fe, UpdateFlags flags)
Allocates space for computed data.
Definition: fe_values.cc:179
fe_values_map.hh
quadrature.hh
Basic definitions of numerical quadrature rules.
MixedPtr< FiniteElement >
FEValues::fill_data_specialized
void fill_data_specialized(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure....
FEValues::shape_grad_component
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:276
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.
FEValues::reinit
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:536
Quadrature
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
FEVectorPiola
@ FEVectorPiola
Definition: finite_element.hh:203
FESystem
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:101
FEValues::init_fe_data
std::shared_ptr< FEInternalData > init_fe_data(const FiniteElement< DIM > &fe, const Quadrature &q)
Precompute finite element data on reference element.
MapContravariant
Helper class allows update values and gradients of FEValues of FEVectorContravariant type.
Definition: fe_values.hh:42
FiniteElement::n_space_components
unsigned int n_space_components(unsigned int spacedim)
Number of components of FE in a mapped space with dimension spacedim.
Definition: finite_element.cc:178
FEValues::FEInternalData::FEInternalData
FEInternalData(unsigned int np, unsigned int nd)
Definition: fe_values.cc:41
FEValues::FEInternalData::n_points
unsigned int n_points
Number of quadrature points.
Definition: fe_values.hh:346
FiniteElement::shape_value
double shape_value(const unsigned int i, const arma::vec::fixed< dim > &p, const unsigned int comp=0) const
Calculates the value of the comp-th component of the i-th shape function at the point p on the refere...
Definition: finite_element.cc:112