Flow123d  JS_before_hm-913-g695b665
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), elm_values(nullptr), fe_data(nullptr)
124 {
125 }
126 
127 
128 
129 template<unsigned int spacedim>
131  if (elm_values != nullptr) delete elm_values;
132  if (fe_data) delete fe_data;
133  for (unsigned int sid=0; sid<side_fe_data.size(); sid++)
134  for (unsigned int pid=0; pid<side_fe_data[sid].size(); pid++)
135  delete side_fe_data[sid][pid];
136 }
137 
138 
139 
140 template<unsigned int spacedim>
141 template<unsigned int DIM>
143  unsigned int n_points,
144  FiniteElement<DIM> & _fe,
145  UpdateFlags _flags)
146 {
147  // For FEVector and FETensor check number of components.
148  // This cannot be done in FiniteElement since it does not know spacedim.
149  if (_fe.type_ == FEVector) {
150  ASSERT_DBG(_fe.n_components() == spacedim).error("FEVector must have spacedim components.");
151  } else if (_fe.type_ == FETensor) {
152  ASSERT_DBG(_fe.n_components() == spacedim*spacedim).error("FETensor must have spacedim*spacedim components.");
153  }
154 
155  dim_ = DIM;
157  n_dofs_ = _fe.n_dofs();
158  n_components_ = _fe.n_space_components(spacedim);
159  fe_type_ = _fe.type_;
160  FESystem<DIM> *fe_sys = dynamic_cast<FESystem<DIM>*>(&_fe);
161  if (fe_sys != nullptr)
162  {
163  for (unsigned int f=0; f<fe_sys->fe().size(); f++)
164  {
165  fe_sys_dofs_.push_back(fe_sys->fe_dofs(f));
166  fe_sys_n_components_.push_back(fe_sys->fe()[f]->n_components());
167  fe_sys_n_space_components_.push_back(fe_sys->fe()[f]->n_space_components(spacedim));
168  }
169  }
170 
171  // add flags required by the finite element or mapping
172  update_flags = _flags | _fe.update_each(_flags);
176 
178  shape_gradients.resize(n_points_, vector<arma::vec::fixed<spacedim> >(n_dofs_*n_components_));
179 
180  views_cache_.initialize(*this, _fe);
181 }
182 
183 
184 
185 template<unsigned int spacedim>
186 template<unsigned int DIM>
188 {
189  ASSERT_DBG( DIM == dim_ );
190  ASSERT_DBG( q.dim() == DIM );
191  FEInternalData *data = new FEInternalData(q.size(), n_dofs_);
192 
194  for (unsigned int i=0; i<q.size(); i++)
195  {
196  for (unsigned int j=0; j<n_dofs_; j++)
197  {
198  for (unsigned int c=0; c<fe.n_components(); c++)
199  shape_values(j,c) = fe.shape_value(j, q.point<DIM>(i), c);
200 
201  data->ref_shape_values[i][j] = trans(shape_values.row(j));
202  }
203  }
204 
205  arma::mat grad(DIM, fe.n_components());
206  for (unsigned int i=0; i<q.size(); i++)
207  {
208  for (unsigned int j=0; j<n_dofs_; j++)
209  {
210  grad.zeros();
211  for (unsigned int c=0; c<fe.n_components(); c++)
212  grad.col(c) += fe.shape_grad(j, q.point<DIM>(i), c);
213 
214  data->ref_shape_grads[i][j] = grad;
215  }
216  }
217 
218  return data;
219 }
220 
221 
222 template<unsigned int spacedim>
223 double FEValues<spacedim>::shape_value(const unsigned int function_no, const unsigned int point_no)
224 {
225  ASSERT_LT_DBG(function_no, n_dofs_);
226  ASSERT_LT_DBG(point_no, n_points_);
227  return shape_values[point_no][function_no];
228 }
229 
230 
231 template<unsigned int spacedim>
232 arma::vec::fixed<spacedim> FEValues<spacedim>::shape_grad(const unsigned int function_no, const unsigned int point_no)
233 {
234  ASSERT_LT_DBG(function_no, n_dofs_);
235  ASSERT_LT_DBG(point_no, n_points_);
236  return shape_gradients[point_no][function_no];
237 }
238 
239 
240 template<unsigned int spacedim>
241 double FEValues<spacedim>::shape_value_component(const unsigned int function_no,
242  const unsigned int point_no,
243  const unsigned int comp) const
244 {
245  ASSERT_LT_DBG(function_no, n_dofs_);
246  ASSERT_LT_DBG(point_no, n_points_);
248  return shape_values[point_no][function_no*n_components_+comp];
249 }
250 
251 
252 template<unsigned int spacedim>
253 arma::vec::fixed<spacedim> FEValues<spacedim>::shape_grad_component(const unsigned int function_no,
254  const unsigned int point_no,
255  const unsigned int comp) const
256 {
257  ASSERT_LT_DBG(function_no, n_dofs_);
258  ASSERT_LT_DBG(point_no, n_points_);
260  return shape_gradients[point_no][function_no*n_components_+comp];
261 }
262 
263 
264 template<unsigned int spacedim>
266 {
268 
269  // shape values
271  for (unsigned int i = 0; i < fe_data.n_points; i++)
272  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
273  shape_values[i][j] = fe_data.ref_shape_values[i][j][0];
274 
275  // shape gradients
277  for (unsigned int i = 0; i < fe_data.n_points; i++)
278  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
279  shape_gradients[i][j] = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
280 }
281 
282 
283 template<unsigned int spacedim>
285  const FEInternalData &fe_data)
286 {
288 
289  // shape values
291  {
292  for (unsigned int i = 0; i < fe_data.n_points; i++)
293  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
294  {
295  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
296  for (unsigned int c=0; c<spacedim; c++)
297  shape_values[i][j*spacedim+c] = fv_vec[c];
298  }
299  }
300 
301  // shape gradients
303  {
304  for (unsigned int i = 0; i < fe_data.n_points; i++)
305  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
306  {
307  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
308  for (unsigned int c=0; c<spacedim; c++)
309  shape_gradients[i][j*spacedim+c] = grads.col(c);
310  }
311  }
312 }
313 
314 
315 template<unsigned int spacedim>
317  const FEInternalData &fe_data)
318 {
320 
321  // shape values
323  {
324  for (unsigned int i = 0; i < fe_data.n_points; i++)
325  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
326  {
327  arma::vec fv_vec = elm_values.jacobian(i) * fe_data.ref_shape_values[i][j];
328  for (unsigned int c=0; c<spacedim; c++)
329  shape_values[i][j*spacedim+c] = fv_vec[c];
330  }
331  }
332 
333  // shape gradients
335  {
336  for (unsigned int i = 0; i < fe_data.n_points; i++)
337  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
338  {
339  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j] * trans(elm_values.jacobian(i));
340  for (unsigned int c=0; c<spacedim; c++)
341  shape_gradients[i][j*spacedim+c] = grads.col(c);
342  }
343  }
344 }
345 
346 
347 template<unsigned int spacedim>
349  const FEInternalData &fe_data)
350 {
352 
353  // shape values
355  {
356  for (unsigned int i = 0; i < fe_data.n_points; i++)
357  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
358  {
359  arma::vec fv_vec = elm_values.jacobian(i)*fe_data.ref_shape_values[i][j]/elm_values.determinant(i);
360  for (unsigned int c=0; c<spacedim; c++)
361  shape_values[i][j*spacedim+c] = fv_vec(c);
362  }
363  }
364 
365  // shape gradients
367  {
368  for (unsigned int i = 0; i < fe_data.n_points; i++)
369  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
370  {
371  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j] * trans(elm_values.jacobian(i))
372  / elm_values.determinant(i);
373  for (unsigned int c=0; c<spacedim; c++)
374  shape_gradients[i][j*spacedim+c] = grads.col(c);
375  }
376  }
377 }
378 
379 
380 template<unsigned int spacedim>
382  const FEInternalData &fe_data)
383 {
385 
386  // shape values
388  {
389  for (unsigned int i = 0; i < fe_data.n_points; i++)
390  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
391  {
392  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
393  for (unsigned int c=0; c<spacedim*spacedim; c++)
394  shape_values[i][j*spacedim*spacedim+c] = fv_vec[c];
395  }
396  }
397 
398  // shape gradients
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::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
405  for (unsigned int c=0; c<spacedim*spacedim; c++)
406  shape_gradients[i][j*spacedim*spacedim+c] = grads.col(c);
407  }
408  }
409 }
410 
411 
412 template<unsigned int spacedim>
414 {
416 
417  // for mixed system we first fill data in sub-elements
418  unsigned int comp_offset = 0;
419  for (unsigned int f=0; f<fe_sys_dofs_.size(); f++)
420  {
421  // fill fe_values for base FE
422  FEInternalData vec_fe_data(fe_data, fe_sys_dofs_[f], comp_offset, fe_sys_n_components_[f]);
423  fe_values_vec[f].fill_data(elm_values, vec_fe_data);
424 
425  comp_offset += fe_sys_n_components_[f];
426  }
427 
428  // shape values
430  {
431  arma::vec fv_vec;
432  unsigned int comp_offset = 0;
433  unsigned int shape_offset = 0;
434  for (unsigned int f=0; f<fe_sys_dofs_.size(); f++)
435  {
436  // gather fe_values in vectors for FESystem
437  for (unsigned int i=0; i<fe_data.n_points; i++)
438  for (unsigned int n=0; n<fe_sys_dofs_[f].size(); n++)
439  for (unsigned int c=0; c<fe_sys_n_space_components_[f]; c++)
440  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];
441 
442  comp_offset += fe_sys_n_space_components_[f];
443  shape_offset += fe_sys_dofs_[f].size()*n_components_;
444  }
445  }
446 
447  // shape gradients
449  {
450  arma::mat grads;
451  unsigned int comp_offset = 0;
452  unsigned int shape_offset = 0;
453  for (unsigned int f=0; f<fe_sys_dofs_.size(); f++)
454  {
455  // gather fe_values in vectors for FESystem
456  for (unsigned int i=0; i<fe_data.n_points; i++)
457  for (unsigned int n=0; n<fe_sys_dofs_[f].size(); n++)
458  for (unsigned int c=0; c<fe_sys_n_space_components_[f]; c++)
459  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];
460 
461  comp_offset += fe_sys_n_space_components_[f];
462  shape_offset += fe_sys_dofs_[f].size()*n_components_;
463  }
464  }
465 
466 }
467 
468 
469 template<unsigned int spacedim>
471 {
472  switch (fe_type_) {
473  case FEScalar:
474  fill_scalar_data(elm_values, fe_data);
475  break;
476  case FEVector:
477  fill_vec_data(elm_values, fe_data);
478  break;
480  fill_vec_contravariant_data(elm_values, fe_data);
481  break;
482  case FEVectorPiola:
483  fill_vec_piola_data(elm_values, fe_data);
484  break;
485  case FETensor:
486  fill_tensor_data(elm_values, fe_data);
487  break;
488  case FEMixedSystem:
489  fill_system_data(elm_values, fe_data);
490  break;
491  default:
492  ASSERT(false).error("Not implemented.");
493  }
494 }
495 
496 
497 
498 
499 
500 
501 
502 
503 template<unsigned int spacedim>
504 template<unsigned int DIM>
506  Quadrature &q,
507  FiniteElement<DIM> &_fe,
508  UpdateFlags _flags)
509 {
510  if (DIM == 0) return; // avoid unnecessary allocation of dummy 0 dimensional objects
511 
512  allocate( q.size(), _fe, _flags);
514 
515  // In case of mixed system allocate data for sub-elements.
516  if (fe_type_ == FEMixedSystem)
517  {
518  FESystem<DIM> *fe = dynamic_cast<FESystem<DIM>*>(&_fe);
519  ASSERT_DBG(fe != nullptr).error("Mixed system must be represented by FESystem.");
520 
521  fe_values_vec.resize(fe->fe().size());
522  for (unsigned int f=0; f<fe->fe().size(); f++)
523  fe_values_vec[f].initialize(q, *fe->fe()[f], update_flags);
524  }
525 
526  // precompute finite element data
527  if ( q.dim() == DIM )
528  {
529  fe_data = init_fe_data(_fe, q);
530  }
531  else if ( q.dim() + 1 == DIM )
532  {
534  for (unsigned int sid = 0; sid < RefElement<DIM>::n_sides; sid++)
535  {
537 
538  // For each side transform the side quadrature points to the cell quadrature points
539  // and then precompute side_fe_data.
540  for (unsigned int pid = 0; pid < RefElement<DIM>::n_side_permutations; pid++)
541  side_fe_data[sid][pid] = init_fe_data(_fe, q.make_from_side<DIM>(sid,pid));
542  }
543  }
544  else
545  ASSERT_DBG(false)(q.dim())(DIM).error("Dimension mismatch in FEValues::initialize().");
546 }
547 
548 
549 
550 
551 
552 template<unsigned int spacedim>
554 {
555  ASSERT_EQ_DBG( dim_, cell.dim() );
556 
557  if (!elm_values->cell().is_valid() ||
558  elm_values->cell() != cell)
559  {
560  elm_values->reinit(cell);
561  }
562 
564 }
565 
566 
567 template<unsigned int spacedim>
568 void FEValues<spacedim>::reinit(const Side &cell_side)
569 {
570  ASSERT_EQ_DBG( dim_, cell_side.dim() );
571 
572  if (!elm_values->side().is_valid() ||
573  elm_values->side() != cell_side)
574  {
575  elm_values->reinit(cell_side);
576  }
577 
578  const LongIdx sid = cell_side.side_idx();
579  const unsigned int pid = elm_values->side().element()->permutation_idx(sid);
580 
581  // calculation of finite element data
582  fill_data(*elm_values, *side_fe_data[sid][pid]);
583 }
584 
585 
586 
588  QGauss::array &quadrature,
590  UpdateFlags flags)
591 {
593  fv[0].initialize(quadrature[0], *fe[0_d], flags);
594  fv[1].initialize(quadrature[1], *fe[1_d], flags);
595  fv[2].initialize(quadrature[2], *fe[2_d], flags);
596  fv[3].initialize(quadrature[3], *fe[3_d], flags);
597  return fv;
598 }
599 
600 
601 
602 // explicit instantiation
607 
608 
609 
610 
611 
612 
613 
614 
615 
616 
617 
618 
619 
620 
621 
622 template class FEValues<3>;
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
std::vector< std::vector< typename FEValues< spacedim >::FEInternalData * > > side_fe_data
Precomputed FE data (shape functions on reference element) for all sides and permuted quadrature poin...
Definition: fe_values.hh:411
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:232
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::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
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:470
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:241
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:413
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:284
#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:553
void allocate(unsigned int n_points, FiniteElement< DIM > &_fe, UpdateFlags flags)
Allocates space for computed data.
Definition: fe_values.cc:142
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:223
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.
FEValues< spacedim >::FEInternalData * fe_data
Precomputed finite element data.
Definition: fe_values.hh:408
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:253
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:505
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
FEInternalData * init_fe_data(const FiniteElement< DIM > &fe, const Quadrature &q)
Precompute finite element data on reference element.
Definition: fe_values.cc:187
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:348
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:316
#define ASSERT_DBG(expr)
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
ElementValues< spacedim > * elm_values
Auxiliary object for calculation of element-dependent data.
Definition: fe_values.hh:396
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:381
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:587
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:265
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