Flow123d  release_3.0.0-875-gdc24e59
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.hh"
22 #include "quadrature/quadrature.hh"
23 #include "fem/finite_element.hh"
24 #include "fem/fe_values.hh"
25 #include "fem/fe_system.hh"
26 
27 
28 
29 using namespace arma;
30 using namespace std;
31 
32 
33 
34 
35 
36 
37 
38 
39 FEInternalData::FEInternalData(unsigned int np, unsigned int nd)
40  : n_points(np),
41  n_dofs(nd)
42 {
43  ref_shape_values.resize(np, vector<arma::vec>(nd));
44  ref_shape_grads.resize(np, vector<arma::mat>(nd));
45 }
46 
47 
48 
49 template<unsigned int dim, unsigned int spacedim>
50 void FEValuesData<dim,spacedim>::allocate(unsigned int size, UpdateFlags flags, unsigned int n_comp)
51 {
52  update_flags = flags;
53 
54  // resize the arrays of computed quantities
55  if (update_flags & update_jacobians)
56  jacobians.resize(size);
57 
58  if (update_flags & update_volume_elements)
59  determinants.resize(size);
60 
61  if ((update_flags & update_JxW_values) |
62  (update_flags & update_side_JxW_values))
63  JxW_values.resize(size);
64 
65  if (update_flags & update_inverse_jacobians)
66  inverse_jacobians.resize(size);
67 
68  if (update_flags & update_values)
69  {
70  shape_values.resize(size, vector<double>(n_comp));
71  }
72 
73  if (update_flags & update_gradients)
74  {
75  shape_gradients.resize(size, vector<arma::vec::fixed<spacedim> >(n_comp));
76  }
77 
78  if (update_flags & update_quadrature_points)
79  points.resize(size);
80 
81  if (update_flags & update_normal_vectors)
82  normal_vectors.resize(size);
83 }
84 
85 
86 template<unsigned int dim, unsigned int spacedim>
88 {
89  scalars.clear();
90  vectors.clear();
91  tensors.clear();
92  switch (fv.get_fe()->type_) {
93  case FEType::FEScalar:
94  scalars.push_back(FEValuesViews::Scalar<dim,spacedim>(fv, 0));
95  break;
96  case FEType::FEVector:
99  vectors.push_back(FEValuesViews::Vector<dim,spacedim>(fv, 0));
100  break;
101  case FEType::FETensor:
102  tensors.push_back(FEValuesViews::Tensor<dim,spacedim>(fv, 0));
103  break;
105  FESystem<dim> *fe = dynamic_cast<FESystem<dim>*>(fv.get_fe());
106  OLD_ASSERT(fe != nullptr, "Mixed system must be represented by FESystem.");
110  for (auto si : sc)
111  scalars.push_back(FEValuesViews::Scalar<dim,spacedim>(fv,si));
112  for (auto vi : vc)
113  vectors.push_back(FEValuesViews::Vector<dim,spacedim>(fv,vi));
114  for (auto ti : tc)
115  tensors.push_back(FEValuesViews::Tensor<dim,spacedim>(fv,ti));
116  break;
117  }
118 }
119 
120 
121 
122 template<unsigned int dim,unsigned int spacedim>
124 : mapping(NULL), quadrature(NULL), fe(NULL), mapping_data(NULL), fe_data(NULL)
125 {
126 }
127 
128 
129 
130 template<unsigned int dim,unsigned int spacedim>
132  if (mapping_data) delete mapping_data;
133  if (fe_data) delete fe_data;
134 }
135 
136 
137 
138 template<unsigned int dim, unsigned int spacedim>
140  Quadrature<dim> & _quadrature,
141  FiniteElement<dim> & _fe,
142  UpdateFlags _flags)
143 {
144  // For FEVector and FETensor check number of components.
145  // This cannot be done in FiniteElement since it does not know spacedim.
146  if (_fe.type_ == FEVector)
147  ASSERT_DBG(_fe.n_components() == spacedim).error("FEVector must have spacedim components.");
148  else if (_fe.type_ == FETensor)
149  ASSERT_DBG(_fe.n_components() == spacedim*spacedim).error("FETensor must have spacedim*spacedim components.");
150 
151  mapping = &_mapping;
152  quadrature = &_quadrature;
153  fe = &_fe;
154 
155  switch (fe->type_) {
156  case FEScalar:
157  n_components_ = 1;
158  break;
159  case FEVector:
161  case FEVectorPiola:
162  n_components_ = spacedim;
163  break;
164  case FETensor:
165  n_components_ = spacedim*spacedim;
166  break;
167  case FEMixedSystem:
168  n_components_ = fe->n_components();
169  break;
170  }
171  // add flags required by the finite element or mapping
172  data.allocate(quadrature->size(), update_each(_flags), fe->n_dofs()*n_components_);
173 
174  views_cache_.initialize(*this);
175 }
176 
177 
178 
179 template<unsigned int dim, unsigned int spacedim>
181 {
182  FEInternalData *data = new FEInternalData(q->size(), fe->n_dofs());
183 
184  arma::mat shape_values(fe->n_dofs(), fe->n_components());
185  for (unsigned int i=0; i<q->size(); i++)
186  {
187  for (unsigned int j=0; j<fe->n_dofs(); j++)
188  {
189  for (unsigned int c=0; c<fe->n_components(); c++)
190  shape_values(j,c) = fe->shape_value(j, q->point(i), c);
191 
192  data->ref_shape_values[i][j] = trans(shape_values.row(j));
193  }
194  }
195 
196  arma::mat grad(dim, fe->n_components());
197  for (unsigned int i=0; i<q->size(); i++)
198  {
199  for (unsigned int j=0; j<fe->n_dofs(); j++)
200  {
201  grad.zeros();
202  for (unsigned int c=0; c<fe->n_components(); c++)
203  grad.col(c) += fe->shape_grad(j, q->point(i), c);
204 
205  data->ref_shape_grads[i][j] = grad;
206  }
207  }
208 
209  return data;
210 }
211 
212 
213 
214 template<unsigned int dim, unsigned int spacedim>
216 {
217  UpdateFlags f = flags | fe->update_each(flags);
218  f |= mapping->update_each(f);
219  return f;
220 }
221 
222 
223 template<unsigned int dim, unsigned int spacedim>
224 double FEValuesBase<dim,spacedim>::shape_value(const unsigned int function_no, const unsigned int point_no)
225 {
226  ASSERT_LT_DBG(function_no, fe->n_dofs());
227  ASSERT_LT_DBG(point_no, quadrature->size());
228  return data.shape_values[point_no][function_no];
229 }
230 
231 
232 template<unsigned int dim, unsigned int spacedim>
233 arma::vec::fixed<spacedim> FEValuesBase<dim,spacedim>::shape_grad(const unsigned int function_no, const unsigned int point_no)
234 {
235  ASSERT_LT_DBG(function_no, fe->n_dofs());
236  ASSERT_LT_DBG(point_no, quadrature->size());
237  return data.shape_gradients[point_no][function_no];
238 }
239 
240 
241 template<unsigned int dim, unsigned int spacedim>
242 double FEValuesBase<dim,spacedim>::shape_value_component(const unsigned int function_no,
243  const unsigned int point_no,
244  const unsigned int comp) const
245 {
246  ASSERT_LT_DBG(function_no, fe->n_dofs());
247  ASSERT_LT_DBG(point_no, quadrature->size());
249  return data.shape_values[point_no][function_no*n_components_+comp];
250 }
251 
252 
253 template<unsigned int dim, unsigned int spacedim>
254 arma::vec::fixed<spacedim> FEValuesBase<dim,spacedim>::shape_grad_component(const unsigned int function_no,
255  const unsigned int point_no,
256  const unsigned int comp) const
257 {
258  ASSERT_LT_DBG(function_no, fe->n_dofs());
259  ASSERT_LT_DBG(point_no, quadrature->size());
261  return data.shape_gradients[point_no][function_no*n_components_+comp];
262 }
263 
264 
265 template<unsigned int dim, unsigned int spacedim>
267 {
268  ASSERT_DBG(fe->type_ == FEScalar);
269 
270  // shape values
271  if (data.update_flags & update_values)
272  for (unsigned int i = 0; i < fe_data.n_points; i++)
273  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
274  data.shape_values[i][j] = fe_data.ref_shape_values[i][j][0];
275 
276  // shape gradients
277  if (data.update_flags & update_gradients)
278  for (unsigned int i = 0; i < fe_data.n_points; i++)
279  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
280  data.shape_gradients[i][j] = trans(data.inverse_jacobians[i]) * fe_data.ref_shape_grads[i][j];
281 }
282 
283 
284 template<unsigned int dim, unsigned int spacedim>
286 {
287  ASSERT_DBG(fe->type_ == FEVector);
288 
289  // shape values
290  if (data.update_flags & update_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  data.shape_values[i][j*spacedim+c] = fv_vec[c];
298  }
299  }
300 
301  // shape gradients
302  if (data.update_flags & update_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(data.inverse_jacobians[i]) * fe_data.ref_shape_grads[i][j];
308  for (unsigned int c=0; c<spacedim; c++)
309  data.shape_gradients[i][j*spacedim+c] = grads.col(c);
310  }
311  }
312 }
313 
314 
315 template<unsigned int dim, unsigned int spacedim>
317 {
319 
320  // shape values
321  if (data.update_flags & update_values)
322  {
323  for (unsigned int i = 0; i < fe_data.n_points; i++)
324  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
325  {
326  arma::vec fv_vec = data.jacobians[i] * fe_data.ref_shape_values[i][j];
327  for (unsigned int c=0; c<spacedim; c++)
328  data.shape_values[i][j*spacedim+c] = fv_vec[c];
329  }
330  }
331 
332  // shape gradients
333  if (data.update_flags & update_gradients)
334  {
335  for (unsigned int i = 0; i < fe_data.n_points; i++)
336  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
337  {
338  arma::mat grads = trans(data.inverse_jacobians[i]) * fe_data.ref_shape_grads[i][j] * trans(data.jacobians[i]);
339  for (unsigned int c=0; c<spacedim; c++)
340  data.shape_gradients[i][j*spacedim+c] = grads.col(c);
341  }
342  }
343 }
344 
345 
346 template<unsigned int dim, unsigned int spacedim>
348 {
349  ASSERT_DBG(fe->type_ == FEVectorPiola);
350 
351  // shape values
352  if (data.update_flags & update_values)
353  {
354  for (unsigned int i = 0; i < fe_data.n_points; i++)
355  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
356  {
357  arma::vec fv_vec = data.jacobians[i]*fe_data.ref_shape_values[i][j]/data.determinants[i];
358  for (unsigned int c=0; c<spacedim; c++)
359  data.shape_values[i][j*spacedim+c] = fv_vec(c);
360  }
361  }
362 
363  // shape gradients
364  if (data.update_flags & update_gradients)
365  {
366  for (unsigned int i = 0; i < fe_data.n_points; i++)
367  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
368  {
369  arma::mat grads = trans(data.inverse_jacobians[i]) * fe_data.ref_shape_grads[i][j] * trans(data.jacobians[i])
370  / data.determinants[i];
371  for (unsigned int c=0; c<spacedim; c++)
372  data.shape_gradients[i][j*spacedim+c] = grads.col(c);
373  }
374  }
375 }
376 
377 
378 template<unsigned int dim, unsigned int spacedim>
380 {
381  ASSERT_DBG(fe->type_ == FETensor);
382 
383  // shape values
384  if (data.update_flags & update_values)
385  {
386  for (unsigned int i = 0; i < fe_data.n_points; i++)
387  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
388  {
389  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
390  for (unsigned int c=0; c<spacedim*spacedim; c++)
391  data.shape_values[i][j*spacedim*spacedim+c] = fv_vec[c];
392  }
393  }
394 
395  // shape gradients
396  if (data.update_flags & update_gradients)
397  {
398  for (unsigned int i = 0; i < fe_data.n_points; i++)
399  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
400  {
401  arma::mat grads = trans(data.inverse_jacobians[i]) * fe_data.ref_shape_grads[i][j];
402  for (unsigned int c=0; c<spacedim*spacedim; c++)
403  data.shape_gradients[i][j*spacedim*spacedim+c] = grads.col(c);
404  }
405  }
406 }
407 
408 
409 template<unsigned int dim, unsigned int spacedim>
411 {
412  ASSERT_DBG(fe->type_ == FEMixedSystem);
413 
414  // for mixed system we first fill data in sub-elements
415  FESystem<dim> *fe_sys = dynamic_cast<FESystem<dim>*>(fe);
416  ASSERT_DBG(fe_sys != nullptr).error("Mixed system must be represented by FESystem.");
417  for (unsigned int f=0; f<fe_sys->fe().size(); f++)
418  {
419  // fill fe_values for base FE
420  fe_values_vec[f]->data.jacobians = data.jacobians;
421  fe_values_vec[f]->data.inverse_jacobians = data.inverse_jacobians;
422  fe_values_vec[f]->data.determinants = data.determinants;
423  fe_values_vec[f]->fill_data(*fe_values_vec[f]->fe_data);
424  }
425 
426  // shape values
427  if (data.update_flags & update_values)
428  {
429  arma::vec fv_vec;
430  unsigned int comp_offset = 0;
431  unsigned int shape_offset = 0;
432  for (unsigned int f=0; f<fe_sys->fe().size(); f++)
433  {
434  // gather fe_values in vectors for FESystem
435  for (unsigned int i=0; i<fe_data.n_points; i++)
436  for (unsigned int n=0; n<fe_sys->fe()[f]->n_dofs(); n++)
437  for (unsigned int c=0; c<fe_sys->fe()[f]->n_components(); c++)
438  data.shape_values[i][shape_offset+fe_sys->function_space_->n_components()*n+comp_offset+c] = fe_values_vec[f]->data.shape_values[i][n*fe_sys->fe()[f]->n_components()+c];
439 
440  comp_offset += fe_sys->fe()[f]->n_components();
441  shape_offset += fe_sys->fe()[f]->n_dofs()*fe_sys->function_space_->n_components();
442  }
443  }
444 
445  // shape gradients
446  if (data.update_flags & update_gradients)
447  {
448  arma::mat grads;
449  unsigned int comp_offset = 0;
450  unsigned int shape_offset = 0;
451  for (unsigned int f=0; f<fe_sys->fe().size(); f++)
452  {
453  // gather fe_values in vectors for FESystem
454  for (unsigned int i=0; i<fe_data.n_points; i++)
455  for (unsigned int n=0; n<fe_sys->fe()[f]->n_dofs(); n++)
456  for (unsigned int c=0; c<fe_sys->fe()[f]->n_components(); c++)
457  data.shape_gradients[i][shape_offset+fe_sys->function_space_->n_components()*n+comp_offset+c] = fe_values_vec[f]->data.shape_gradients[i][n*fe_sys->fe()[f]->n_components()+c];
458 
459  comp_offset += fe_sys->fe()[f]->n_components();
460  shape_offset += fe_sys->fe()[f]->n_dofs()*fe_sys->function_space_->n_components();
461  }
462  }
463 
464 }
465 
466 
467 template<unsigned int dim, unsigned int spacedim>
469 {
470  switch (fe->type_) {
471  case FEScalar:
472  fill_scalar_data(fe_data);
473  break;
474  case FEVector:
475  fill_vec_data(fe_data);
476  break;
479  break;
480  case FEVectorPiola:
481  fill_vec_piola_data(fe_data);
482  break;
483  case FETensor:
484  fill_tensor_data(fe_data);
485  break;
486  case FEMixedSystem:
487  fill_system_data(fe_data);
488  break;
489  default:
490  ASSERT(false).error("Not implemented.");
491  }
492 }
493 
494 
495 
496 
497 
498 
499 
500 
501 template<unsigned int dim, unsigned int spacedim>
503  Quadrature<dim> &_quadrature,
504  FiniteElement<dim> &_fe,
505  UpdateFlags _flags)
506 :FEValuesBase<dim, spacedim>()
507 {
508  this->allocate(_mapping, _quadrature, _fe, _flags);
509 
510  // precompute the maping data and finite element data
511  this->mapping_data = this->mapping->initialize(*this->quadrature, this->data.update_flags);
512  this->fe_data = this->init_fe_data(this->quadrature);
513 
514  // In case of mixed system allocate data for sub-elements.
515  if (this->fe->type_ == FEMixedSystem)
516  {
517  FESystem<dim> *fe = dynamic_cast<FESystem<dim>*>(this->fe);
518  ASSERT_DBG(fe != nullptr).error("Mixed system must be represented by FESystem.");
519 
520  for (auto fe_sub : fe->fe())
521  this->fe_values_vec.push_back(make_shared<FEValues<dim,spacedim> >(_mapping, _quadrature, *fe_sub, _flags));
522  }
523 }
524 
525 
526 
527 template<unsigned int dim,unsigned int spacedim>
529 {
530  OLD_ASSERT_EQUAL( dim, cell->dim() );
531  this->data.present_cell = &cell;
532 
533  // calculate Jacobian of mapping, JxW, inverse Jacobian, normal vector(s)
534  this->mapping->fill_fe_values(cell,
535  *this->quadrature,
536  *this->mapping_data,
537  this->data);
538 
539  this->fill_data(*this->fe_data);
540 }
541 
542 
543 
544 
545 
546 
547 
548 
549 
550 template<unsigned int dim,unsigned int spacedim>
552  Quadrature<dim-1> & _sub_quadrature,
553  FiniteElement<dim> & _fe,
554  const UpdateFlags _flags)
555 :FEValuesBase<dim,spacedim>()
556 {
557  sub_quadrature = &_sub_quadrature;
558  Quadrature<dim> *q = new Quadrature<dim>(_sub_quadrature.size());
559  this->allocate(_mapping, *q, _fe, _flags);
560 
561  for (unsigned int sid = 0; sid < RefElement<dim>::n_sides; sid++)
562  {
563  for (unsigned int pid = 0; pid < RefElement<dim>::n_side_permutations; pid++)
564  {
565  // transform the side quadrature points to the cell quadrature points
566  side_quadrature[sid][pid] = Quadrature<dim>(_sub_quadrature, sid, pid);
567  side_mapping_data[sid][pid] = this->mapping->initialize(side_quadrature[sid][pid], this->data.update_flags);
568  side_fe_data[sid][pid] = this->init_fe_data(&side_quadrature[sid][pid]);
569  }
570  }
571 
572  // In case of mixed system allocate data for sub-elements.
573  if (this->fe->type_ == FEMixedSystem)
574  {
575  FESystem<dim> *fe = dynamic_cast<FESystem<dim>*>(this->fe);
576  ASSERT_DBG(fe != nullptr).error("Mixed system must be represented by FESystem.");
577 
578  for (auto fe_sub : fe->fe())
579  this->fe_values_vec.push_back(make_shared<FESideValues<dim,spacedim> >(_mapping, _sub_quadrature, *fe_sub, _flags));
580  }
581 }
582 
583 
584 
585 template<unsigned int dim,unsigned int spacedim>
587 {
588  for (unsigned int sid=0; sid<RefElement<dim>::n_sides; sid++)
589  {
590  for (unsigned int pid=0; pid<RefElement<dim>::n_side_permutations; pid++)
591  {
592  delete side_mapping_data[sid][pid];
593  delete side_fe_data[sid][pid];
594  }
595  }
596 
597  // Since quadrature is an auxiliary internal variable allocated
598  // by the constructor, it must be destroyed here.
599  delete this->quadrature;
600 }
601 
602 
603 template<unsigned int dim,unsigned int spacedim>
605  unsigned int sid)
606 {
607  ASSERT_LT_DBG( sid, cell->n_sides());
608  ASSERT_EQ_DBG(dim, cell->dim());
609  this->data.present_cell = &cell;
610 
611  unsigned int pid = cell->permutation_idx(sid);
613  // calculate Jacobian of mapping, JxW, inverse Jacobian, normal vector(s)
614  this->mapping->fill_fe_side_values(cell,
615  sid,
616  side_quadrature[sid][pid],
617  *side_mapping_data[sid][pid],
618  this->data);
619 
620  // calculation of finite element data
621  this->fill_data(*side_fe_data[sid][pid]);
622 }
623 
624 
625 
626 
627 
628 
629 
630 template class FEValuesBase<0,3>;
631 template class FEValuesBase<1,3>;
632 template class FEValuesBase<2,3>;
633 template class FEValuesBase<3,3>;
634 
635 template class FEValues<0,3>;
636 template class FEValues<1,3>;
637 template class FEValues<2,3>;
638 template class FEValues<3,3>;
639 
640 template class FESideValues<1,3>;
641 template class FESideValues<2,3>;
642 template class FESideValues<3,3>;
643 
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
void fill_vec_contravariant_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for vectorial FE.
Definition: fe_values.cc:316
const std::vector< std::shared_ptr< FiniteElement< dim > > > & fe()
Definition: fe_system.hh:140
Transformed quadrature weight for cell sides.
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:331
std::vector< unsigned int > get_tensor_components() const
Definition: fe_system.hh:135
void reinit(ElementAccessor< 3 > &cell, unsigned int sid)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:604
void fill_scalar_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for scalar FE.
Definition: fe_values.cc:266
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:254
Determinant of the Jacobian.
void fill_data(const FEInternalData &fe_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure...
Definition: fe_values.cc:468
void allocate(Mapping< dim, spacedim > &_mapping, Quadrature< dim > &_quadrature, FiniteElement< dim > &_fe, UpdateFlags flags)
Allocates space for computed data.
Definition: fe_values.cc:139
MappingInternalData * mapping_data
Precomputed mapping data.
Definition: fe_values.hh:495
void initialize(FEValuesBase &fv)
Definition: fe_values.cc:87
UpdateFlags update_each(UpdateFlags flags)
Determine quantities to be recomputed on each cell.
Definition: fe_values.cc:215
virtual ~FEValuesBase()
Definition: fe_values.cc:131
FEValuesBase()
Default constructor.
Definition: fe_values.cc:123
Class FEValues calculates finite element data on the actual cells such as shape function values...
FEValuesData< dim, spacedim > data
Data computed by the mapping and finite element.
Definition: fe_values.hh:505
Class FESystem for compound finite elements.
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:346
std::vector< unsigned int > get_scalar_components() const
Definition: fe_system.hh:129
FEInternalData * fe_data
Precomputed finite element data.
Definition: fe_values.hh:500
FiniteElement< dim > * get_fe() const
Returns the finite element in use.
Definition: fe_values.hh:432
void allocate(unsigned int size, UpdateFlags flags, unsigned n_comp)
Resize the data arrays.
Definition: fe_values.cc:50
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: fe_values.hh:35
Volume element.
unsigned int n_components_
Number of components of the FE.
Definition: fe_values.hh:511
virtual ~FESideValues()
Destructor.
Definition: fe_values.cc:586
FEValues(Mapping< dim, spacedim > &_mapping, Quadrature< dim > &_quadrature, FiniteElement< dim > &_fe, UpdateFlags _flags)
Constructor.
Definition: fe_values.cc:502
unsigned int dim() const
Definition: elements.h:124
#define OLD_ASSERT(...)
Definition: global_defs.h:131
Transformed quadrature points.
std::vector< unsigned int > get_vector_components() const
Definition: fe_system.hh:132
Abstract class for the mapping between reference and actual cell.
Definition: fe_values.hh:38
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:99
void fill_tensor_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for tensorial FE.
Definition: fe_values.cc:379
unsigned int n_sides() const
Definition: elements.h:135
void fill_system_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for mixed system of FE.
Definition: fe_values.cc:410
void fill_vec_piola_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for Raviart-Thomas FE.
Definition: fe_values.cc:347
unsigned int n_components() const
Returns numer of components of the basis function.
Basic definitions of numerical quadrature rules.
const Quadrature< dim-1 > * sub_quadrature
Quadrature for the integration on the element sides.
Definition: fe_values.hh:617
Shape function gradients.
Definition: update_flags.hh:95
MappingInternalData * side_mapping_data[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
Definition: fe_values.hh:621
Quadrature< dim > side_quadrature[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
Definition: fe_values.hh:619
FEInternalData * init_fe_data(const Quadrature< dim > *q)
Precompute finite element data on reference element.
Definition: fe_values.cc:180
FEInternalData(unsigned int np, unsigned int nd)
Definition: fe_values.cc:39
Quadrature< dim > * quadrature
The quadrature rule used to calculate integrals.
Definition: fe_values.hh:485
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:224
Class Mapping calculates data related to the mapping of the reference cell to the actual cell...
Normal vectors.
void fill_vec_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for vectorial FE.
Definition: fe_values.cc:285
std::vector< std::vector< arma::mat > > ref_shape_grads
Precomputed gradients of basis functions at the quadrature points.
Definition: fe_values.hh:69
FEInternalData * side_fe_data[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
Definition: fe_values.hh:623
FESideValues(Mapping< dim, spacedim > &_mapping, Quadrature< dim-1 > &_sub_quadrature, FiniteElement< dim > &_fe, UpdateFlags flags)
Constructor.
Definition: fe_values.cc:551
ViewsCache views_cache_
Auxiliary storage of FEValuesViews accessors.
Definition: fe_values.hh:514
std::shared_ptr< FunctionSpace > function_space_
Function space defining the FE.
Mapping< dim, spacedim > * mapping
The mapping from the reference cell to the actual cell.
Definition: fe_values.hh:480
const unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:139
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:242
FiniteElement< dim > * fe
The used finite element.
Definition: fe_values.hh:490
#define ASSERT_DBG(expr)
Definition: asserts.hh:349
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:233
const arma::vec::fixed< dim > & point(const unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:144
Calculates finite element data on the actual cell.
Definition: fe_values.hh:532
unsigned int n_points
Number of quadrature points.
Definition: fe_values.hh:72
#define OLD_ASSERT_EQUAL(a, b)
Definition: global_defs.h:133
std::vector< std::shared_ptr< FEValuesBase< dim, spacedim > > > fe_values_vec
Vector of FEValues for sub-elements of FESystem.
Definition: fe_values.hh:508
Structure for storing the precomputed finite element data.
Definition: fe_values.hh:47
std::vector< std::vector< arma::vec > > ref_shape_values
Precomputed values of basis functions at the quadrature points.
Definition: fe_values.hh:60
Abstract class for description of finite elements.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Base class for FEValues and FESideValues.
Definition: fe_values.hh:37
unsigned int n_dofs
Number of dofs (shape functions).
Definition: fe_values.hh:75
unsigned int permutation_idx(unsigned int prm_idx) const
Return permutation_idx of given index.
Definition: elements.h:144
void reinit(ElementAccessor< 3 > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:528
Transformed quadrature weights.
FEType type_
Type of FiniteElement.
Calculates finite element data on a side.
Definition: fe_values.hh:577
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:299