Flow123d  JS_before_hm-1711-g43b724b37
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_DBG(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_DBG(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_DBG(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].resize(RefElement<DIM>::n_side_permutations);
169 
170  // For each side transform the side quadrature points to the cell quadrature points
171  // and then precompute side_fe_data.
172  for (unsigned int pid = 0; pid < RefElement<DIM>::n_side_permutations; pid++)
173  side_fe_data[sid][pid] = init_fe_data(_fe, q.make_from_side<DIM>(sid,pid));
174  }
175  }
176  else
177  ASSERT_DBG(false)(q.dim())(DIM).error("Dimension mismatch in FEValues::initialize().");
178 }
179 
180 
181 
182 template<unsigned int spacedim>
183 template<unsigned int DIM>
185  unsigned int n_points,
186  FiniteElement<DIM> & _fe,
187  UpdateFlags _flags)
188 {
189  // For FEVector and FETensor check number of components.
190  // This cannot be done in FiniteElement since it does not know spacedim.
191  if (_fe.type_ == FEVector) {
192  ASSERT_DBG(_fe.n_components() == spacedim).error("FEVector must have spacedim components.");
193  } else if (_fe.type_ == FETensor) {
194  ASSERT_DBG(_fe.n_components() == spacedim*spacedim).error("FETensor must have spacedim*spacedim components.");
195  }
196 
197  fe_sys_dofs_.clear();
198  fe_sys_n_components_.clear();
199  fe_sys_n_space_components_.clear();
200 
201  dim_ = DIM;
202  n_points_ = n_points;
203  n_dofs_ = _fe.n_dofs();
204  n_components_ = _fe.n_space_components(spacedim);
205  fe_type_ = _fe.type_;
206  FESystem<DIM> *fe_sys = dynamic_cast<FESystem<DIM>*>(&_fe);
207  if (fe_sys != nullptr)
208  {
209  for (unsigned int f=0; f<fe_sys->fe().size(); f++)
210  {
211  fe_sys_dofs_.push_back(fe_sys->fe_dofs(f));
212  fe_sys_n_components_.push_back(fe_sys->fe()[f]->n_components());
213  fe_sys_n_space_components_.push_back(fe_sys->fe()[f]->n_space_components(spacedim));
214  }
215  }
216 
217  // add flags required by the finite element or mapping
218  update_flags = _flags | _fe.update_each(_flags);
219  update_flags |= MappingP1<DIM,spacedim>::update_each(update_flags);
220  if (update_flags & update_values)
221  shape_values.resize(n_points_, vector<double>(n_dofs_*n_components_));
222 
223  if (update_flags & update_gradients)
224  shape_gradients.resize(n_points_, vector<arma::vec::fixed<spacedim> >(n_dofs_*n_components_));
225 
226  views_cache_.initialize(*this, _fe);
227 }
228 
229 
230 
231 template<unsigned int spacedim>
232 template<unsigned int DIM>
233 std::shared_ptr<typename FEValues<spacedim>::FEInternalData> FEValues<spacedim>::init_fe_data(const FiniteElement<DIM> &fe, const Quadrature &q)
234 {
235  ASSERT_DBG( DIM == dim_ );
236  ASSERT_DBG( q.dim() == DIM );
237  std::shared_ptr<FEInternalData> data = std::make_shared<FEInternalData>(q.size(), n_dofs_);
238 
239  arma::mat shape_values(n_dofs_, fe.n_components());
240  for (unsigned int i=0; i<q.size(); i++)
241  {
242  for (unsigned int j=0; j<n_dofs_; j++)
243  {
244  for (unsigned int c=0; c<fe.n_components(); c++)
245  shape_values(j,c) = fe.shape_value(j, q.point<DIM>(i), c);
246 
247  data->ref_shape_values[i][j] = trans(shape_values.row(j));
248  }
249  }
250 
251  arma::mat grad(DIM, fe.n_components());
252  for (unsigned int i=0; i<q.size(); i++)
253  {
254  for (unsigned int j=0; j<n_dofs_; j++)
255  {
256  grad.zeros();
257  for (unsigned int c=0; c<fe.n_components(); c++)
258  grad.col(c) += fe.shape_grad(j, q.point<DIM>(i), c);
259 
260  data->ref_shape_grads[i][j] = grad;
261  }
262  }
263 
264  return data;
265 }
266 
267 
268 /*template<unsigned int spacedim>
269 double FEValues<spacedim>::shape_value_component(const unsigned int function_no,
270  const unsigned int point_no,
271  const unsigned int comp) const
272 {
273  ASSERT_LT_DBG(function_no, n_dofs_);
274  ASSERT_LT_DBG(point_no, n_points_);
275  ASSERT_LT_DBG(comp, n_components_);
276  return shape_values[point_no][function_no*n_components_+comp];
277 }*/
278 
279 
280 template<unsigned int spacedim>
281 arma::vec::fixed<spacedim> FEValues<spacedim>::shape_grad_component(const unsigned int function_no,
282  const unsigned int point_no,
283  const unsigned int comp) const
284 {
285  ASSERT_LT_DBG(function_no, n_dofs_);
286  ASSERT_LT_DBG(point_no, n_points_);
287  ASSERT_LT_DBG(comp, n_components_);
288  return shape_gradients[point_no][function_no*n_components_+comp];
289 }
290 
291 
292 /*template<unsigned int spacedim>
293 void FEValues<spacedim>::fill_scalar_data(const ElementValues<spacedim> &elm_values, const FEInternalData &fe_data)
294 {
295  ASSERT_DBG(fe_type_ == FEScalar);
296 
297  // shape values
298  if (update_flags & update_values)
299  for (unsigned int i = 0; i < fe_data.n_points; i++)
300  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
301  shape_values[i][j] = fe_data.ref_shape_values[i][j][0];
302 
303  // shape gradients
304  if (update_flags & update_gradients)
305  for (unsigned int i = 0; i < fe_data.n_points; i++)
306  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
307  shape_gradients[i][j] = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
308 }
309 
310 
311 template<unsigned int spacedim>
312 void FEValues<spacedim>::fill_vec_data(const ElementValues<spacedim> &elm_values,
313  const FEInternalData &fe_data)
314 {
315  ASSERT_DBG(fe_type_ == FEVector);
316 
317  // shape values
318  if (update_flags & update_values)
319  {
320  for (unsigned int i = 0; i < fe_data.n_points; i++)
321  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
322  {
323  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
324  for (unsigned int c=0; c<spacedim; c++)
325  shape_values[i][j*spacedim+c] = fv_vec[c];
326  }
327  }
328 
329  // shape gradients
330  if (update_flags & update_gradients)
331  {
332  for (unsigned int i = 0; i < fe_data.n_points; i++)
333  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
334  {
335  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
336  for (unsigned int c=0; c<spacedim; c++)
337  shape_gradients[i][j*spacedim+c] = grads.col(c);
338  }
339  }
340 }
341 
342 
343 template<unsigned int spacedim>
344 void FEValues<spacedim>::fill_vec_contravariant_data(const ElementValues<spacedim> &elm_values,
345  const FEInternalData &fe_data)
346 {
347  ASSERT_DBG(fe_type_ == FEVectorContravariant);
348 
349  // shape values
350  if (update_flags & update_values)
351  {
352  for (unsigned int i = 0; i < fe_data.n_points; i++)
353  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
354  {
355  arma::vec fv_vec = elm_values.jacobian(i) * fe_data.ref_shape_values[i][j];
356  for (unsigned int c=0; c<spacedim; c++)
357  shape_values[i][j*spacedim+c] = fv_vec[c];
358  }
359  }
360 
361  // shape gradients
362  if (update_flags & update_gradients)
363  {
364  for (unsigned int i = 0; i < fe_data.n_points; i++)
365  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
366  {
367  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j] * trans(elm_values.jacobian(i));
368  for (unsigned int c=0; c<spacedim; c++)
369  shape_gradients[i][j*spacedim+c] = grads.col(c);
370  }
371  }
372 }
373 
374 
375 template<unsigned int spacedim>
376 void FEValues<spacedim>::fill_vec_piola_data(const ElementValues<spacedim> &elm_values,
377  const FEInternalData &fe_data)
378 {
379  ASSERT_DBG(fe_type_ == FEVectorPiola);
380 
381  // shape values
382  if (update_flags & update_values)
383  {
384  for (unsigned int i = 0; i < fe_data.n_points; i++)
385  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
386  {
387  arma::vec fv_vec = elm_values.jacobian(i)*fe_data.ref_shape_values[i][j]/elm_values.determinant(i);
388  for (unsigned int c=0; c<spacedim; c++)
389  shape_values[i][j*spacedim+c] = fv_vec(c);
390  }
391  }
392 
393  // shape gradients
394  if (update_flags & update_gradients)
395  {
396  for (unsigned int i = 0; i < fe_data.n_points; i++)
397  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
398  {
399  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j] * trans(elm_values.jacobian(i))
400  / elm_values.determinant(i);
401  for (unsigned int c=0; c<spacedim; c++)
402  shape_gradients[i][j*spacedim+c] = grads.col(c);
403  }
404  }
405 }
406 
407 
408 template<unsigned int spacedim>
409 void FEValues<spacedim>::fill_tensor_data(const ElementValues<spacedim> &elm_values,
410  const FEInternalData &fe_data)
411 {
412  ASSERT_DBG(fe_type_ == FETensor);
413 
414  // shape values
415  if (update_flags & update_values)
416  {
417  for (unsigned int i = 0; i < fe_data.n_points; i++)
418  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
419  {
420  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
421  for (unsigned int c=0; c<spacedim*spacedim; c++)
422  shape_values[i][j*spacedim*spacedim+c] = fv_vec[c];
423  }
424  }
425 
426  // shape gradients
427  if (update_flags & update_gradients)
428  {
429  for (unsigned int i = 0; i < fe_data.n_points; i++)
430  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
431  {
432  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
433  for (unsigned int c=0; c<spacedim*spacedim; c++)
434  shape_gradients[i][j*spacedim*spacedim+c] = grads.col(c);
435  }
436  }
437 }
438 
439 
440 template<unsigned int spacedim>
441 void FEValues<spacedim>::fill_system_data(const ElementValues<spacedim> &elm_values, const FEInternalData &fe_data)
442 {
443  ASSERT_DBG(fe_type_ == FEMixedSystem);
444 
445  // for mixed system we first fill data in sub-elements
446  unsigned int comp_offset = 0;
447  for (unsigned int f=0; f<fe_sys_dofs_.size(); f++)
448  {
449  // fill fe_values for base FE
450  FEInternalData vec_fe_data(fe_data, fe_sys_dofs_[f], comp_offset, fe_sys_n_components_[f]);
451  fe_values_vec[f].fill_data(elm_values, vec_fe_data); // fe_values.fe_values_vec
452 
453  comp_offset += fe_sys_n_components_[f];
454  }
455 
456  // shape values
457  if (update_flags & update_values)
458  {
459  arma::vec fv_vec;
460  unsigned int comp_offset = 0;
461  unsigned int shape_offset = 0;
462  for (unsigned int f=0; f<fe_sys_dofs_.size(); f++)
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_dofs_[f].size(); n++)
467  for (unsigned int c=0; c<fe_sys_n_space_components_[f]; c++)
468  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];
469 
470  comp_offset += fe_sys_n_space_components_[f];
471  shape_offset += fe_sys_dofs_[f].size()*n_components_;
472  }
473  }
474 
475  // shape gradients
476  if (update_flags & update_gradients)
477  {
478  arma::mat grads;
479  unsigned int comp_offset = 0;
480  unsigned int shape_offset = 0;
481  for (unsigned int f=0; f<fe_sys_dofs_.size(); f++)
482  {
483  // gather fe_values in vectors for FESystem
484  for (unsigned int i=0; i<fe_data.n_points; i++)
485  for (unsigned int n=0; n<fe_sys_dofs_[f].size(); n++)
486  for (unsigned int c=0; c<fe_sys_n_space_components_[f]; c++)
487  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];
488 
489  comp_offset += fe_sys_n_space_components_[f];
490  shape_offset += fe_sys_dofs_[f].size()*n_components_;
491  }
492  }
493 
494 }*/
495 
496 
497 template<unsigned int spacedim>
498 void FEValues<spacedim>::fill_data(const ElementValues<spacedim> &elm_values, const FEInternalData &fe_data)
499 {
500  switch (fe_type_) {
501  case FEScalar:
502  this->fill_data_specialized<MapScalar<spacedim>>(elm_values, fe_data);
503  break;
504  case FEVector:
505  this->fill_data_specialized<MapVector<spacedim>>(elm_values, fe_data);
506  break;
508  this->fill_data_specialized<MapContravariant<spacedim>>(elm_values, fe_data);
509  break;
510  case FEVectorPiola:
511  this->fill_data_specialized<MapPiola<spacedim>>(elm_values, fe_data);
512  break;
513  case FETensor:
514  this->fill_data_specialized<MapTensor<spacedim>>(elm_values, fe_data);
515  break;
516  case FEMixedSystem:
517  this->fill_data_specialized<MapSystem<spacedim>>(elm_values, fe_data);
518  break;
519  default:
520  ASSERT(false).error("Not implemented.");
521  }
522 }
523 
524 
525 
526 template<unsigned int spacedim>
527 template<class MapType>
529  MapType map_type;
530  map_type.fill_values_vec(*this, elm_values, fe_data);
531  if (update_flags & update_values)
532  map_type.update_values(*this, elm_values, fe_data);
533  if (update_flags & update_gradients)
534  map_type.update_gradients(*this, elm_values, fe_data);
535 }
536 
537 
538 
539 
540 template<unsigned int spacedim>
542 {
543  ASSERT_EQ_DBG( dim_, cell.dim() );
544 
545  if (!elm_values->cell().is_valid() ||
546  elm_values->cell() != cell)
547  {
548  elm_values->reinit(cell);
549  }
550 
551  fill_data(*elm_values, *fe_data);
552 }
553 
554 
555 template<unsigned int spacedim>
556 void FEValues<spacedim>::reinit(const Side &cell_side)
557 {
558  ASSERT_EQ_DBG( dim_, cell_side.dim() );
559 
560  if (!elm_values->side().is_valid() ||
561  elm_values->side() != cell_side)
562  {
563  elm_values->reinit(cell_side);
564  }
565 
566  const LongIdx sid = cell_side.side_idx();
567  const unsigned int pid = elm_values->side().element()->permutation_idx(sid);
568 
569  // calculation of finite element data
570  fill_data(*elm_values, *side_fe_data[sid][pid]);
571 }
572 
573 
574 
576  QGauss::array &quadrature,
578  UpdateFlags flags)
579 {
581  fv[0].initialize(quadrature[0], *fe[0_d], flags);
582  fv[1].initialize(quadrature[1], *fe[1_d], flags);
583  fv[2].initialize(quadrature[2], *fe[2_d], flags);
584  fv[3].initialize(quadrature[3], *fe[3_d], flags);
585  return fv;
586 }
587 
588 
589 
590 // explicit instantiation
601 
602 
603 
604 
605 
606 
607 
608 
609 
610 
611 
612 
613 
614 
615 
616 template class FEValues<3>;
ElementAccessor::dim
unsigned int dim() const
Definition: accessors.hh:157
Quadrature::point
Armor::ArmaVec< double, point_dim > point(unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:92
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:163
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)
Allow use shorter versions of macro names if these names is not used with external library.
Definition: asserts.hh:347
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,...
ASSERT_DBG
#define ASSERT_DBG(expr)
Definition: include_fadbad.hh:28
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:87
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:213
std::vector< arma::vec >
ElementAccessor
Definition: dh_cell_accessor.hh:32
ASSERT_EQ_DBG
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:332
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.
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:575
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:498
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
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:361
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
Quadrature::dim
unsigned int dim() const
Definition: quadrature.hh:73
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:415
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:184
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....
ASSERT_LT_DBG
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:300
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:281
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:541
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