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