Flow123d  DF_patch_fe_data_tables-1a16803
fe_values.hh
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.hh
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 #ifndef FE_VALUES_HH_
22 #define FE_VALUES_HH_
23 
24 #include <string.h> // for memcpy
25 #include <algorithm> // for swap
26 #include <new> // for operator new[]
27 #include <string> // for operator<<
28 #include <vector> // for vector
29 #include "fem/element_values.hh" // for ElementValues
30 #include "fem/fe_values_views.hh" // for FEValuesViews
31 #include "mesh/ref_element.hh" // for RefElement
32 #include "mesh/accessors.hh"
33 #include "fem/update_flags.hh" // for UpdateFlags
34 #include "tools/mixed.hh"
36 #include "fields/eval_subset.hh"
37 
38 class Quadrature;
39 template<unsigned int dim> class FiniteElement;
40 
41 template<class FV, unsigned int dim> class MapScalar;
42 template<class FV, unsigned int dim> class MapPiola;
43 template<class FV, unsigned int dim> class MapContravariant;
44 template<class FV, unsigned int dim> class MapVector;
45 template<class FV, unsigned int dim> class MapTensor;
46 template<class FV, unsigned int dim> class MapSystem;
47 
48 template<unsigned int spcedim> class FEValues;
49 template<unsigned int spcedim> class PatchFEValues_TEMP;
50 
51 
52 
53 
54 
55 
56 
58 
59 
60 /// Structure for storing the precomputed finite element data.
62 {
63 public:
64 
65  FEInternalData(unsigned int np, unsigned int nd);
66 
67  /// Create a new instance of FEInternalData for a FESystem component or subvector.
68  FEInternalData(const FEInternalData &fe_system_data,
69  const std::vector<unsigned int> &dof_indices,
70  unsigned int first_component_idx,
71  unsigned int ncomps = 1);
72 
73  /**
74  * @brief Precomputed values of basis functions at the quadrature points.
75  *
76  * Dimensions: (no. of quadrature points)
77  * x (no. of dofs)
78  * x (no. of components in ref. cell)
79  */
81 
82  /**
83  * @brief Precomputed gradients of basis functions at the quadrature points.
84  *
85  * Dimensions: (no. of quadrature points)
86  * x (no. of dofs)
87  * x ((dim_ of. ref. cell)x(no. of components in ref. cell))
88  */
90 
91  /// Number of quadrature points.
92  unsigned int n_points;
93 
94  /// Number of dofs (shape functions).
95  unsigned int n_dofs;
96 };
97 
98 
99 template<class FV, unsigned int spacedim = 3>
101 {
102 protected:
103  // internal structure that stores all possible views
104  // for scalar and vector-valued components of the FE
105  struct ViewsCache {
109 
110  template<unsigned int DIM>
111  void initialize(const FV &fv, const FiniteElement<DIM> &fe);
112  };
113 
114 
115 public:
116  /// Default constructor with postponed initialization.
117  FEValuesBase();
118 
119  /**
120  * @brief Initialize structures and calculates cell-independent data.
121  *
122  * @param _quadrature The quadrature rule for the cell associated
123  * to given finite element or for the cell side.
124  * @param _fe The finite element.
125  * @param _flags The update flags.
126  */
127  template<unsigned int DIM>
128  void initialize(Quadrature &_quadrature,
129  FiniteElement<DIM> &_fe,
130  UpdateFlags _flags);
131 
132  /**
133  * @brief Allocates space for computed data.
134  *
135  * @param n_points Number of quadrature points.
136  * @param _fe The finite element.
137  * @param flags The update flags.
138  */
139  template<unsigned int DIM>
140  void allocate(Quadrature &_quadrature,
141  FiniteElement<DIM> &_fe,
142  UpdateFlags flags);
143 
144  /**
145  * @brief Returns the number of quadrature points.
146  */
147  inline unsigned int n_points() const
148  { return n_points_; }
149 
150  /**
151  * @brief Returns the number of shape functions.
152  */
153  inline unsigned int n_dofs() const
154  {
155  return n_dofs_;
156  }
157 
158  /// Return dimension of reference space.
159  inline unsigned int dim() const
160  { return dim_; }
161 
162  /**
163  * @brief Accessor to scalar values of multicomponent FE.
164  * @param i Index of scalar component.
165  */
167  {
168  ASSERT_LT(i, views_cache_.scalars.size());
169  return views_cache_.scalars[i];
170  }
171 
172  /**
173  * @brief Accessor to vector values of multicomponent FE.
174  * @param i Index of first vector component.
175  */
177  {
178  ASSERT_LT(i, views_cache_.vectors.size());
179  return views_cache_.vectors[i];
180  }
181 
182  /**
183  * @brief Accessor to tensor values of multicomponent FE.
184  * @param i Index of first tensor component.
185  */
187  {
188  ASSERT_LT(i, views_cache_.tensors.size());
189  return views_cache_.tensors[i];
190  }
191 
192 protected:
193  /**
194  * @brief Computes the shape function values and gradients on the actual cell
195  * and fills the FEValues structure.
196  *
197  * @param fe_data Precomputed finite element data.
198  */
199  void fill_data(const ElementValues<spacedim> &elm_values, const FEInternalData &fe_data);
200 
201  /**
202  * @brief Computes the shape function values and gradients on the actual cell
203  * and fills the FEValues structure. Specialized variant of previous method for
204  * different FETypes given by template parameter.
205  */
206  template<class MapType>
207  void fill_data_specialized(const ElementValues<spacedim> &elm_values, const FEInternalData &fe_data);
208 
209 
210  /// Initialize vectors declared separately in descendants.
211  virtual void allocate_in(unsigned int) =0;
212 
213  /// Initialize ElementValues separately in descendants.
214  virtual void initialize_in(Quadrature &, unsigned int) =0;
215 
216  /// Initialize @p fe_values_vec only in PatchFEValues_TEMP.
217  virtual void init_fe_val_vec() =0;
218 
219  /// Precompute finite element data on reference element.
220  template<unsigned int DIM>
221  std::shared_ptr<FEInternalData> init_fe_data(const FiniteElement<DIM> &fe, const Quadrature &q);
222 
223  /// Dimension of reference space.
224  unsigned int dim_;
225 
226  /// Number of integration points.
227  unsigned int n_points_;
228 
229  /// Number of finite element dofs.
230  unsigned int n_dofs_;
231 
232  /// Type of finite element (scalar, vector, tensor).
234 
235  /// Dof indices of FESystem sub-elements.
237 
238  /// Numbers of components of FESystem sub-elements in reference space.
240 
241  /// Numbers of components of FESystem sub-elements in real space.
243 
244  /// Flags that indicate which finite element quantities are to be computed.
246 
247  /// Vector of FEValues for sub-elements of FESystem.
249 
250  /// Number of components of the FE.
251  unsigned int n_components_;
252 
253  /// Auxiliary storage of FEValuesViews accessors.
255 
256  /// Precomputed finite element data.
257  std::shared_ptr<FEInternalData> fe_data_;
258 
259  /// Precomputed FE data (shape functions on reference element) for all side quadrature points.
261 
262  /// Helper object, we need its for ViewsCache initialization.
263  FV *fv_;
264 
265  friend class MapScalar<FV, spacedim>;
266  friend class MapPiola<FV, spacedim>;
267  friend class MapContravariant<FV, spacedim>;
268  friend class MapVector<FV, spacedim>;
269  friend class MapTensor<FV, spacedim>;
270  friend class MapSystem<FV, spacedim>;
271 };
272 
273 
274 
275 
276 
277 /**
278  * @brief Calculates finite element data on the actual cell.
279  *
280  * FEValues takes care of the calculation of finite element data on
281  * the actual cell or on its side, (values of shape functions at quadrature
282  * points, gradients of shape functions, Jacobians of the mapping
283  * from the reference cell etc.).
284  * @param spacedim Dimension of the Euclidean space where the actual
285  * cell lives.
286  */
287 template<unsigned int spacedim = 3>
288 class FEValues : public FEValuesBase<FEValues<spacedim>, spacedim>
289 {
290 public:
291 
292  /// Default constructor with postponed initialization.
293  FEValues();
294 
295  /// Constructor with initialization of data structures
296  /// (see initialize() for description of parameters).
297  template<unsigned int DIM>
298  FEValues(Quadrature &_quadrature,
299  FiniteElement<DIM> &_fe,
300  UpdateFlags _flags)
301  : FEValues()
302  {
303  this->initialize(_quadrature, _fe, _flags);
304  }
305 
306 
307  /// Correct deallocation of objects created by 'initialize' methods.
308  ~FEValues();
309 
310 
311  /**
312  * @brief Update cell-dependent data (gradients, Jacobians etc.)
313  *
314  * @param cell The actual cell.
315  */
316  void reinit(const ElementAccessor<spacedim> &cell);
317 
318  /**
319  * @brief Update cell side-dependent FE data (values, gradients).
320  *
321  * @param cell_side Accessor to cell side.
322  */
323  void reinit(const Side &cell_side);
324 
325  /**
326  * @brief Return the value of the @p function_no-th shape function at
327  * the @p point_no-th quadrature point.
328  *
329  * @param function_no Number of the shape function.
330  * @param point_no Number of the quadrature point.
331  */
332  inline double shape_value(const unsigned int function_no, const unsigned int point_no) const
333  {
334  ASSERT_LT(function_no, this->n_dofs_);
335  ASSERT_LT(point_no, this->n_points_);
336  return shape_values_[point_no][function_no];
337  }
338 
339 
340  /**
341  * @brief Return the gradient of the @p function_no-th shape function at
342  * the @p point_no-th quadrature point.
343  *
344  * @param function_no Number of the shape function.
345  * @param point_no Number of the quadrature point.
346  */
347  inline arma::vec::fixed<spacedim> shape_grad(const unsigned int function_no, const unsigned int point_no) const
348  {
349  ASSERT_LT(function_no, this->n_dofs_);
350  ASSERT_LT(point_no, this->n_points_);
351  return shape_gradients_[point_no][function_no];
352  }
353 
354  /**
355  * @brief Return the value of the @p function_no-th shape function at
356  * the @p point_no-th quadrature point.
357  *
358  * For vectorial finite elements.
359  *
360  * @param function_no Number of the shape function.
361  * @param point_no Number of the quadrature point.
362  */
363  inline double shape_value_component(const unsigned int function_no,
364  const unsigned int point_no,
365  const unsigned int comp) const
366  {
367  ASSERT_LT(function_no, this->n_dofs_);
368  ASSERT_LT(point_no, this->n_points_);
369  ASSERT_LT(comp, this->n_components_);
370  return shape_values_[point_no][function_no*this->n_components_+comp];
371  }
372 
373  /**
374  * @brief Return the gradient of the @p function_no-th shape function at
375  * the @p point_no-th quadrature point.
376  *
377  * For vectorial finite elements.
378  *
379  * @param function_no Number of the shape function.
380  * @param point_no Number of the quadrature point.
381  */
382  arma::vec::fixed<spacedim> shape_grad_component(const unsigned int function_no,
383  const unsigned int point_no,
384  const unsigned int comp) const;
385 
386  /**
387  * Set shape value @p val of the @p i_point and @p i_func_comp.
388  */
389  inline void set_shape_value(unsigned int i_point, unsigned int i_func_comp, double val)
390  {
391  shape_values_[i_point][i_func_comp] = val;
392  }
393 
394  /**
395  * Set shape gradient @p val of the @p i_point and @p i_func_comp.
396  */
397  inline void set_shape_gradient(unsigned int i_point, unsigned int i_func_comp, arma::vec::fixed<spacedim> val)
398  {
399  shape_gradients_[i_point][i_func_comp] = val;
400  }
401 
402  /**
403  * @brief Return the relative volume change of the cell (Jacobian determinant).
404  *
405  * If dim_==spacedim then the sign may be negative, otherwise the
406  * result is a positive number.
407  *
408  * @param point_no Number of the quadrature point.
409  */
410  inline double determinant(const unsigned int point_no)
411  {
412  ASSERT_LT(point_no, this->n_points_);
413  return elm_values_->determinant(point_no);
414  }
415 
416  /**
417  * @brief Return the product of Jacobian determinant and the quadrature
418  * weight at given quadrature point.
419  *
420  * @param point_no Number of the quadrature point.
421  */
422  inline double JxW(const unsigned int point_no)
423  {
424  ASSERT_LT(point_no, this->n_points_);
425  // TODO: This is temporary solution to distinguish JxW on element and side_JxW on side.
426  // In future we should call the appropriate method in elm_values_.
427  return (elm_values_->cell().is_valid()) ? elm_values_->JxW(point_no) : elm_values_->side_JxW(point_no);
428  }
429 
430  /**
431  * @brief Return coordinates of the quadrature point in the actual cell system.
432  *
433  * @param point_no Number of the quadrature point.
434  */
435  inline arma::vec::fixed<spacedim> point(const unsigned int point_no)
436  {
437  ASSERT_LT(point_no, this->n_points_);
438  return elm_values_->point(point_no);
439  }
440 
441  /**
442  * @brief Return coordinates of all quadrature points in the actual cell system.
443  *
444  */
445  inline const Armor::array &point_list() const
446  {
447  return elm_values_->point_list();
448  }
449 
450 
451  /**
452  * @brief Returns the normal vector to a side at given quadrature point.
453  *
454  * @param point_no Number of the quadrature point.
455  */
456  inline arma::vec::fixed<spacedim> normal_vector(unsigned int point_no)
457  {
458  ASSERT_LT(point_no, this->n_points_);
459  return elm_values_->normal_vector(point_no);
460  }
461 
462 
463 protected:
464 
465  /// Implement @p FEValuesBase::allocate_in
466  void allocate_in(unsigned int q_dim) override;
467 
468  /// Implement @p FEValuesBase::initialize_in
469  void initialize_in(Quadrature &q, unsigned int dim) override;
470 
471  /// Implement @p FEValuesBase::initialize_in
472  void init_fe_val_vec() override
473  {}
474 
475  /// Shape functions evaluated at the quadrature points.
477 
478  /// Gradients of shape functions evaluated at the quadrature points.
479  /// Each row of the matrix contains the gradient of one shape function.
481 
482  /// Auxiliary object for calculation of element-dependent data.
483  std::shared_ptr<ElementValues<spacedim> > elm_values_;
484 
485  friend class MapScalar<FEValues<spacedim>, spacedim>;
486  friend class MapPiola<FEValues<spacedim>, spacedim>;
487  friend class MapContravariant<FEValues<spacedim>, spacedim>;
488  friend class MapVector<FEValues<spacedim>, spacedim>;
489  friend class MapTensor<FEValues<spacedim>, spacedim>;
490  friend class MapSystem<FEValues<spacedim>, spacedim>;
491 };
492 
493 
494 
495 
496 
497 
498 template<unsigned int spacedim = 3>
499 class PatchFEValues_TEMP : public FEValuesBase<PatchFEValues_TEMP<spacedim>, spacedim> {
500 public:
501  /// Constructor, set maximal number of elements on patch
502  PatchFEValues_TEMP(unsigned int max_size=0);
503 
504  /// Reinit data.
505  void reinit(PatchElementsList patch_elements);
506 
507  inline unsigned int used_size() const {
508  return used_size_;
509  }
510 
511  inline unsigned int max_size() const {
512  return element_data_.size();
513  }
514 
515  /// Set element that is selected for processing. Element is given by index on patch.
516  inline void get_cell(const unsigned int patch_cell_idx) {
517  patch_data_idx_ = element_patch_map_.find(patch_cell_idx)->second;
518  }
519 
520  /// Set element and side that are selected for processing. Element is given by index on patch.
521  inline void get_side(unsigned int patch_cell_idx, unsigned int side_idx) {
522  patch_data_idx_ = element_patch_map_.find(patch_cell_idx)->second + side_idx;
523  }
524 
525  /**
526  * @brief Return the value of the @p function_no-th shape function at
527  * the @p point_no-th quadrature point.
528  *
529  * @param function_no Number of the shape function.
530  * @param point_no Number of the quadrature point.
531  */
532  inline double shape_value(const unsigned int function_no, const unsigned int point_no) const
533  {
534  // TODO: obsolete method, will be replaced with following two methods.
535  ASSERT_LT(function_no, this->n_dofs_);
536  ASSERT_LT(point_no, this->n_points_);
537  return element_data_[patch_data_idx_].shape_values_[point_no][function_no];
538  }
539 
540 
541  /**
542  * @brief Return the value of the @p function_no-th shape function at
543  * the @p p quadrature point.
544  *
545  * @param function_no Number of the shape function.
546  * @param p BulkPoint corresponds to the quadrature point.
547  */
548  inline double shape_value(const unsigned int function_no, const BulkPoint &p) const
549  {
550  ASSERT_LT(function_no, this->n_dofs_);
551  unsigned int patch_data_idx = element_patch_map_.find(p.elem_patch_idx())->second;
552  return element_data_[patch_data_idx].shape_values_[p.eval_point_idx()][function_no];
553  }
554 
555 
556  /**
557  * @brief Return the value of the @p function_no-th shape function at
558  * the @p p quadrature point.
559  *
560  * @param function_no Number of the shape function.
561  * @param p SidePoint corresponds to the quadrature point.
562  */
563  inline double shape_value(const unsigned int function_no, const SidePoint &p) const
564  {
565  ASSERT_LT(function_no, this->n_dofs_);
566  unsigned int patch_data_idx = element_patch_map_.find(p.elem_patch_idx())->second + p.side_idx();
567  return element_data_[patch_data_idx].shape_values_[p.local_point_idx()][function_no];
568  }
569 
570 
571  /**
572  * @brief Return the gradient of the @p function_no-th shape function at
573  * the @p point_no-th quadrature point.
574  *
575  * @param function_no Number of the shape function.
576  * @param point_no Number of the quadrature point.
577  */
578  inline arma::vec::fixed<spacedim> shape_grad(const unsigned int function_no, const unsigned int point_no) const
579  {
580  // TODO: obsolete method, will be replaced with following two methods.
581  ASSERT_LT(function_no, this->n_dofs_);
582  ASSERT_LT(point_no, this->n_points_);
583  return element_data_[patch_data_idx_].shape_gradients_[point_no][function_no];
584  }
585 
586  /**
587  * @brief Return the gradient of the @p function_no-th shape function at
588  * the @p p quadrature point.
589  *
590  * @param function_no Number of the shape function.
591  * @param p BulkPoint corresponds to the quadrature point.
592  */
593  inline arma::vec::fixed<spacedim> shape_grad(const unsigned int function_no, const BulkPoint &p) const
594  {
595  ASSERT_LT(function_no, this->n_dofs_);
596  unsigned int patch_data_idx = element_patch_map_.find(p.elem_patch_idx())->second;
597  return element_data_[patch_data_idx].shape_gradients_[p.eval_point_idx()][function_no];;
598  }
599 
600  /**
601  * @brief Return the gradient of the @p function_no-th shape function at
602  * the @p p quadrature point.
603  *
604  * @param function_no Number of the shape function.
605  * @param p SidePoint corresponds to the quadrature point.
606  */
607  inline arma::vec::fixed<spacedim> shape_grad(const unsigned int function_no, const SidePoint &p) const
608  {
609  ASSERT_LT(function_no, this->n_dofs_);
610  unsigned int patch_data_idx = element_patch_map_.find(p.elem_patch_idx())->second + p.side_idx();
611  return element_data_[patch_data_idx].shape_gradients_[p.local_point_idx()][function_no];;
612  }
613 
614  /**
615  * @brief Return the value of the @p function_no-th shape function at
616  * the @p point_no-th quadrature point.
617  *
618  * For vectorial finite elements.
619  *
620  * @param function_no Number of the shape function.
621  * @param point_no Number of the quadrature point.
622  */
623  inline double shape_value_component(const unsigned int function_no,
624  const unsigned int point_no,
625  const unsigned int comp) const
626  {
627  ASSERT_LT(function_no, this->n_dofs_);
628  ASSERT_LT(point_no, this->n_points_);
629  ASSERT_LT(comp, this->n_components_);
630  return element_data_[patch_data_idx_].shape_values_[point_no][function_no*this->n_components_+comp];
631  }
632 
633  /**
634  * @brief Return the gradient of the @p function_no-th shape function at
635  * the @p point_no-th quadrature point.
636  *
637  * For vectorial finite elements.
638  *
639  * @param function_no Number of the shape function.
640  * @param point_no Number of the quadrature point.
641  */
642  arma::vec::fixed<spacedim> shape_grad_component(const unsigned int function_no,
643  const unsigned int point_no,
644  const unsigned int comp) const;
645 
646 
647  /**
648  * Set shape value @p val of the @p i_point and @p i_func_comp.
649  */
650  inline void set_shape_value(unsigned int i_point, unsigned int i_func_comp, double val)
651  {
652  element_data_[patch_data_idx_].shape_values_[i_point][i_func_comp] = val;
653  }
654 
655  /**
656  * Set shape gradient @p val of the @p i_point and @p i_func_comp.
657  */
658  inline void set_shape_gradient(unsigned int i_point, unsigned int i_func_comp, arma::vec::fixed<spacedim> val)
659  {
660  element_data_[patch_data_idx_].shape_gradients_[i_point][i_func_comp] = val;
661  }
662 
663  /**
664  * @brief Return the relative volume change of the cell (Jacobian determinant).
665  *
666  * If dim_==spacedim then the sign may be negative, otherwise the
667  * result is a positive number.
668  *
669  * @param point_no Number of the quadrature point.
670  */
671  inline double determinant(const unsigned int point_no)
672  {
673  ASSERT_LT(point_no, this->n_points_);
674  return element_data_[patch_data_idx_].elm_values_->determinant(point_no);
675  }
676 
677  /**
678  * @brief Return the product of Jacobian determinant and the quadrature
679  * weight at given quadrature point.
680  *
681  * @param point_no Number of the quadrature point.
682  */
683  inline double JxW(const unsigned int point_no)
684  {
685  // TODO: obsolete method, will be replaced with following two methods.
686  ASSERT_LT(point_no, this->n_points_);
687  return (object_type_==ElementFE) ? element_data_[patch_data_idx_].elm_values_->JxW(point_no)
688  : element_data_[patch_data_idx_].elm_values_->side_JxW(point_no);
689  }
690 
691  /**
692  * @brief Return the product of Jacobian determinant and the quadrature
693  * weight at given quadrature point.
694  *
695  * @param p BulkPoint corresponds to the quadrature point.
696  */
697  inline double JxW(const BulkPoint &p)
698  {
699  unsigned int patch_data_idx = element_patch_map_.find(p.elem_patch_idx())->second;
700  return element_data_[patch_data_idx].elm_values_->JxW(p.eval_point_idx());
701  }
702 
703  /**
704  * @brief Return the product of Jacobian determinant and the quadrature
705  * weight at given quadrature point.
706  *
707  * @param p SidePoint corresponds to the quadrature point.
708  */
709  inline double JxW(const SidePoint &p)
710  {
711  unsigned int patch_data_idx = element_patch_map_.find(p.elem_patch_idx())->second + p.side_idx();
712  return element_data_[patch_data_idx].elm_values_->side_JxW(p.local_point_idx());
713  }
714 
715  /**
716  * @brief Returns the normal vector to a side at given quadrature point.
717  *
718  * @param point_no Number of the quadrature point.
719  */
720  inline arma::vec::fixed<spacedim> normal_vector(unsigned int point_no)
721  {
722  // TODO: obsolete method, will be replaced with following method.
723  ASSERT_LT(point_no, this->n_points_);
724  return element_data_[patch_data_idx_].elm_values_->normal_vector(point_no);
725  }
726 
727  /**
728  * @brief Returns the normal vector to a side at given quadrature point.
729  *
730  * @param p SidePoint corresponds to the quadrature point.
731  */
732  inline arma::vec::fixed<spacedim> normal_vector(const SidePoint &p)
733  {
734  unsigned int patch_data_idx = element_patch_map_.find(p.elem_patch_idx())->second + p.side_idx();
735  return element_data_[patch_data_idx].elm_values_->normal_vector(p.local_point_idx());
736  }
737 
738 
739 protected:
742  SideFE = 1
743  };
744 
746  {
747  public:
749 
750  /// Shape functions evaluated at the quadrature points.
752 
753  /// Gradients of shape functions evaluated at the quadrature points.
754  /// Each row of the matrix contains the gradient of one shape function.
756 
757  /// Auxiliary object for calculation of element-dependent data.
758  std::shared_ptr<ElementValues<spacedim> > elm_values_;
759 
760  };
761 
762  /// Set size of ElementFEData. Important: Use only during the initialization of FESystem !
763  void resize(unsigned int max_size) {
764  element_data_.resize(max_size);
765  }
766 
767  /// Implement @p FEValuesBase::allocate_in
768  void allocate_in(unsigned int q_dim) override;
769 
770  /// Implement @p FEValuesBase::initialize_in
771  void initialize_in(Quadrature &q, unsigned int dim) override;
772 
773  /// Implement @p FEValuesBase::initialize_in
774  void init_fe_val_vec() override;
775 
776  /// Patch index of processed element / side.
777  unsigned int patch_data_idx_;
778 
779  /// Map of element patch indexes to element_data_.
781 
782  /// Data of elements / sides on patch
784 
785  /// Number of elements / sides on patch. Must be less or equal to size of element_data vector
786  unsigned int used_size_;
787 
788  /// Maximal number of elements on patch.
789  unsigned int max_n_elem_;
790 
791  /// Distinguishes using of PatchFEValues_TEMP for storing data of elements or sides.
793 
794 
795  friend class MapScalar<PatchFEValues_TEMP<spacedim>, spacedim>;
796  friend class MapPiola<PatchFEValues_TEMP<spacedim>, spacedim>;
797  friend class MapContravariant<PatchFEValues_TEMP<spacedim>, spacedim>;
798  friend class MapVector<PatchFEValues_TEMP<spacedim>, spacedim>;
799  friend class MapTensor<PatchFEValues_TEMP<spacedim>, spacedim>;
800  friend class MapSystem<PatchFEValues_TEMP<spacedim>, spacedim>;
801 };
802 
803 
804 
806  QGauss::array &quadrature,
808  UpdateFlags flags);
809 
810 
811 
812 
813 
814 
815 
816 #endif /* FE_VALUES_HH_ */
#define ASSERT_LT(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:301
Base point accessor class.
Definition: eval_subset.hh:55
unsigned int elem_patch_idx() const
Definition: eval_subset.hh:79
unsigned int eval_point_idx() const
Return index in EvalPoints object.
Definition: eval_subset.hh:84
Class for computation of data on cell and side.
Structure for storing the precomputed finite element data.
Definition: fe_values.hh:62
std::vector< std::vector< arma::mat > > ref_shape_grads
Precomputed gradients of basis functions at the quadrature points.
Definition: fe_values.hh:89
unsigned int n_points
Number of quadrature points.
Definition: fe_values.hh:92
FEInternalData(unsigned int np, unsigned int nd)
Definition: fe_values.cc:47
unsigned int n_dofs
Number of dofs (shape functions).
Definition: fe_values.hh:95
std::vector< std::vector< arma::vec > > ref_shape_values
Precomputed values of basis functions at the quadrature points.
Definition: fe_values.hh:80
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
Definition: fe_values.hh:245
void allocate(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags flags)
Allocates space for computed data.
Definition: fe_values.cc:173
unsigned int n_points_
Number of integration points.
Definition: fe_values.hh:227
FV * fv_
Helper object, we need its for ViewsCache initialization.
Definition: fe_values.hh:263
std::shared_ptr< FEInternalData > init_fe_data(const FiniteElement< DIM > &fe, const Quadrature &q)
Precompute finite element data on reference element.
Definition: fe_values.cc:218
std::vector< shared_ptr< FEInternalData > > side_fe_data_
Precomputed FE data (shape functions on reference element) for all side quadrature points.
Definition: fe_values.hh:260
std::vector< unsigned int > fe_sys_n_space_components_
Numbers of components of FESystem sub-elements in real space.
Definition: fe_values.hh:242
FEType fe_type_
Type of finite element (scalar, vector, tensor).
Definition: fe_values.hh:233
unsigned int n_dofs_
Number of finite element dofs.
Definition: fe_values.hh:230
const FEValuesViews::Scalar< FV, spacedim > & scalar_view(unsigned int i) const
Accessor to scalar values of multicomponent FE.
Definition: fe_values.hh:166
unsigned int n_points() const
Returns the number of quadrature points.
Definition: fe_values.hh:147
virtual void init_fe_val_vec()=0
Initialize fe_values_vec only in PatchFEValues_TEMP.
ViewsCache views_cache_
Auxiliary storage of FEValuesViews accessors.
Definition: fe_values.hh:254
virtual void allocate_in(unsigned int)=0
Initialize vectors declared separately in descendants.
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:254
std::vector< FV > fe_values_vec
Vector of FEValues for sub-elements of FESystem.
Definition: fe_values.hh:248
std::vector< unsigned int > fe_sys_n_components_
Numbers of components of FESystem sub-elements in reference space.
Definition: fe_values.hh:239
const FEValuesViews::Tensor< FV, spacedim > & tensor_view(unsigned int i) const
Accessor to tensor values of multicomponent FE.
Definition: fe_values.hh:186
std::vector< std::vector< unsigned int > > fe_sys_dofs_
Dof indices of FESystem sub-elements.
Definition: fe_values.hh:236
unsigned int n_components_
Number of components of the FE.
Definition: fe_values.hh:251
unsigned int dim() const
Return dimension of reference space.
Definition: fe_values.hh:159
FEValuesBase()
Default constructor with postponed initialization.
Definition: fe_values.cc:41
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....
Definition: fe_values.cc:284
unsigned int n_dofs() const
Returns the number of shape functions.
Definition: fe_values.hh:153
void initialize(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
Definition: fe_values.cc:129
virtual void initialize_in(Quadrature &, unsigned int)=0
Initialize ElementValues separately in descendants.
std::shared_ptr< FEInternalData > fe_data_
Precomputed finite element data.
Definition: fe_values.hh:257
const FEValuesViews::Vector< FV, spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
Definition: fe_values.hh:176
unsigned int dim_
Dimension of reference space.
Definition: fe_values.hh:224
Calculates finite element data on the actual cell.
Definition: fe_values.hh:289
double determinant(const unsigned int point_no)
Return the relative volume change of the cell (Jacobian determinant).
Definition: fe_values.hh:410
~FEValues()
Correct deallocation of objects created by 'initialize' methods.
Definition: fe_values.cc:314
void set_shape_value(unsigned int i_point, unsigned int i_func_comp, double val)
Definition: fe_values.hh:389
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
Definition: fe_values.hh:422
FEValues()
Default constructor with postponed initialization.
Definition: fe_values.cc:309
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:343
double shape_value(const unsigned int function_no, const unsigned int point_no) const
Return the value of the function_no-th shape function at the point_no-th quadrature point.
Definition: fe_values.hh:332
const Armor::array & point_list() const
Return coordinates of all quadrature points in the actual cell system.
Definition: fe_values.hh:445
FEValues(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Definition: fe_values.hh:298
std::vector< std::vector< double > > shape_values_
Shape functions evaluated at the quadrature points.
Definition: fe_values.hh:476
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:560
void initialize_in(Quadrature &q, unsigned int dim) override
Implement FEValuesBase::initialize_in.
Definition: fe_values.cc:319
void init_fe_val_vec() override
Implement FEValuesBase::initialize_in.
Definition: fe_values.hh:472
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
Definition: fe_values.hh:456
void set_shape_gradient(unsigned int i_point, unsigned int i_func_comp, arma::vec::fixed< spacedim > val)
Definition: fe_values.hh:397
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.hh:363
std::shared_ptr< ElementValues< spacedim > > elm_values_
Auxiliary object for calculation of element-dependent data.
Definition: fe_values.hh:483
arma::vec::fixed< spacedim > point(const unsigned int point_no)
Return coordinates of the quadrature point in the actual cell system.
Definition: fe_values.hh:435
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const unsigned int point_no) const
Return the gradient of the function_no-th shape function at the point_no-th quadrature point.
Definition: fe_values.hh:347
std::vector< std::vector< arma::vec::fixed< spacedim > > > shape_gradients_
Definition: fe_values.hh:480
void allocate_in(unsigned int q_dim) override
Implement FEValuesBase::allocate_in.
Definition: fe_values.cc:329
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Helper class allows update values and gradients of FEValues of FEVectorContravariant type.
Helper class allows update values and gradients of FEValues of FEVectorPiola type.
Helper class allows update values and gradients of FEValues of FEScalar type.
Helper class allows update values and gradients of FEValues of FEMixedSystem type.
Helper class allows update values and gradients of FEValues of FETensor type.
Helper class allows update values and gradients of FEValues of FEVector type.
std::shared_ptr< ElementValues< spacedim > > elm_values_
Auxiliary object for calculation of element-dependent data.
Definition: fe_values.hh:758
std::vector< std::vector< arma::vec::fixed< spacedim > > > shape_gradients_
Definition: fe_values.hh:755
std::vector< std::vector< double > > shape_values_
Shape functions evaluated at the quadrature points.
Definition: fe_values.hh:751
void get_cell(const unsigned int patch_cell_idx)
Set element that is selected for processing. Element is given by index on patch.
Definition: fe_values.hh:516
PatchFEValues_TEMP(unsigned int max_size=0)
Constructor, set maximal number of elements on patch.
double shape_value(const unsigned int function_no, const BulkPoint &p) const
Return the value of the function_no-th shape function at the p quadrature point.
Definition: fe_values.hh:548
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
Definition: fe_values.hh:720
unsigned int used_size() const
Definition: fe_values.hh:507
arma::vec::fixed< spacedim > normal_vector(const SidePoint &p)
Returns the normal vector to a side at given quadrature point.
Definition: fe_values.hh:732
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const SidePoint &p) const
Return the gradient of the function_no-th shape function at the p quadrature point.
Definition: fe_values.hh:607
void allocate_in(unsigned int q_dim) override
Implement FEValuesBase::allocate_in.
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
Definition: fe_values.hh:683
unsigned int max_n_elem_
Maximal number of elements on patch.
Definition: fe_values.hh:789
void set_shape_value(unsigned int i_point, unsigned int i_func_comp, double val)
Definition: fe_values.hh:650
MeshObjectType object_type_
Distinguishes using of PatchFEValues_TEMP for storing data of elements or sides.
Definition: fe_values.hh:792
std::map< unsigned int, unsigned int > element_patch_map_
Map of element patch indexes to element_data_.
Definition: fe_values.hh:780
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.
void get_side(unsigned int patch_cell_idx, unsigned int side_idx)
Set element and side that are selected for processing. Element is given by index on patch.
Definition: fe_values.hh:521
double shape_value(const unsigned int function_no, const SidePoint &p) const
Return the value of the function_no-th shape function at the p quadrature point.
Definition: fe_values.hh:563
double determinant(const unsigned int point_no)
Return the relative volume change of the cell (Jacobian determinant).
Definition: fe_values.hh:671
unsigned int max_size() const
Definition: fe_values.hh:511
std::vector< ElementFEData > element_data_
Data of elements / sides on patch.
Definition: fe_values.hh:783
unsigned int patch_data_idx_
Patch index of processed element / side.
Definition: fe_values.hh:777
void initialize_in(Quadrature &q, unsigned int dim) override
Implement FEValuesBase::initialize_in.
double JxW(const SidePoint &p)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
Definition: fe_values.hh:709
void init_fe_val_vec() override
Implement FEValuesBase::initialize_in.
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const BulkPoint &p) const
Return the gradient of the function_no-th shape function at the p quadrature point.
Definition: fe_values.hh:593
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.hh:623
double JxW(const BulkPoint &p)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
Definition: fe_values.hh:697
double shape_value(const unsigned int function_no, const unsigned int point_no) const
Return the value of the function_no-th shape function at the point_no-th quadrature point.
Definition: fe_values.hh:532
unsigned int used_size_
Number of elements / sides on patch. Must be less or equal to size of element_data vector.
Definition: fe_values.hh:786
void set_shape_gradient(unsigned int i_point, unsigned int i_func_comp, arma::vec::fixed< spacedim > val)
Definition: fe_values.hh:658
void resize(unsigned int max_size)
Set size of ElementFEData. Important: Use only during the initialization of FESystem !
Definition: fe_values.hh:763
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const unsigned int point_no) const
Return the gradient of the function_no-th shape function at the point_no-th quadrature point.
Definition: fe_values.hh:578
void reinit(PatchElementsList patch_elements)
Reinit data.
std::array< QGauss, 4 > array
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
General point a+ side_begin_ + ccessor allow iterate over quadrature points of given side defined in ...
Definition: eval_subset.hh:116
virtual unsigned int side_idx() const =0
Intermediate step in implementation of PatcFEValues.
unsigned int local_point_idx() const
Return local index in quadrature. Temporary method - intermediate step in implementation of PatcFEVal...
Definition: eval_subset.hh:135
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
std::vector< FEValues< 3 > > mixed_fe_values(QGauss::array &quadrature, MixedPtr< FiniteElement > fe, UpdateFlags flags)
Definition: fe_values.cc:685
std::list< std::pair< ElementAccessor< 3 >, unsigned int > > PatchElementsList
Definition: fe_values.hh:49
FEType
Definitions of particular quadrature rules on simplices.
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc.
vector< FEValuesViews::Scalar< FV, spacedim > > scalars
Definition: fe_values.hh:106
vector< FEValuesViews::Vector< FV, spacedim > > vectors
Definition: fe_values.hh:107
void initialize(const FV &fv, const FiniteElement< DIM > &fe)
Definition: fe_values.cc:74
vector< FEValuesViews::Tensor< FV, spacedim > > tensors
Definition: fe_values.hh:108
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
Definition: update_flags.hh:68