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