Flow123d  DF_patch_fe_data_tables-956a7bc
patch_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 patch_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, David Flanderka
19  */
20 
21 #ifndef PATCH_FE_VALUES_HH_
22 #define PATCH_FE_VALUES_HH_
23 
24 
25 #include <string.h> // for memcpy
26 #include <algorithm> // for swap
27 #include <new> // for operator new[]
28 #include <string> // for operator<<
29 #include <vector> // for vector
30 #include "fem/element_values.hh" // for ElementValues
31 #include "fem/fe_values.hh" // for FEValuesBase
32 #include "fem/fe_values_views.hh" // for FEValuesViews
33 #include "fem/fe_system.hh" // for FESystem
34 #include "fem/eigen_tools.hh"
36 #include "fem/mapping_p1.hh"
37 #include "mesh/ref_element.hh" // for RefElement
38 #include "mesh/accessors.hh"
39 #include "fem/update_flags.hh" // for UpdateFlags
41 #include "fields/eval_subset.hh"
42 #include "fem/arena_resource.hh"
43 #include "fem/arena_vec.hh"
44 
45 template<unsigned int spacedim> class PatchFEValues;
46 
47 
48 
49 typedef typename std::vector< std::array<uint, 3> > DimPointTable; ///< Holds triplet (dim; bulk/side; idx of point in subtable)
50 
51 
52 template <class ValueType>
53 class ElQ {
54 public:
55  /// Forbidden default constructor
56  ElQ() = delete;
57 
58  /// Constructor
59  ElQ(PatchPointValues<3> &patch_point_vals, unsigned int op_idx)
60  : patch_point_vals_(patch_point_vals), op_idx_(op_idx) {}
61 
62  ValueType operator()(FMT_UNUSED const BulkPoint &point);
63 
64  ValueType operator()(FMT_UNUSED const SidePoint &point);
65 
66 private:
67  PatchPointValues<3> &patch_point_vals_; ///< Reference to PatchPointValues
68  unsigned int op_idx_; ///< Index of operation in patch_point_vals_.operations vector
69 };
70 
71 
72 template <class ValueType>
73 class FeQ {
74 public:
75  /// Forbidden default constructor
76  FeQ() = delete;
77 
78  // Class similar to current FeView
79  FeQ(PatchPointValues<3> &patch_point_vals, unsigned int op_idx, unsigned int n_dofs)
80  : patch_point_vals_(patch_point_vals), op_idx_(op_idx), n_dofs_(n_dofs) {}
81 
82 
83  ValueType operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const BulkPoint &point);
84 
85  ValueType operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const SidePoint &point);
86 
87  // Implementation for EdgePoint, SidePoint, and JoinPoint shoud have a common implementation
88  // resolving to side values
89 
90 private:
91  PatchPointValues<3> &patch_point_vals_; ///< Reference to PatchPointValues
92  unsigned int op_idx_; ///< Index of operation in patch_point_vals_.operations vector
93  unsigned int n_dofs_; ///< Number of DOFs
94 };
95 
96 
97 template <class ValueType>
99 public:
100  /// Default constructor
102  : patch_point_vals_bulk_(nullptr), patch_point_vals_side_(nullptr), join_idx_(-1) {}
103 
104  /**
105  * Constructor
106  *
107  * @param patch_point_vals_bulk Pointer to PatchPointValues bulk object.
108  * @param patch_point_vals_side Pointer to PatchPointValues side object.
109  * @param begin Index of the first component of the bulk Quantity.
110  * @param begin_side Index of the first component of the side Quantity.
111  * @param n_dofs_bulk Number of DOFs of bulk (lower-dim) element.
112  * @param n_dofs_side Number of DOFs of side (higher-dim) element.
113  * @param join_idx Index function.
114  */
115  JoinShapeAccessor(PatchPointValues<3> *patch_point_vals_bulk, PatchPointValues<3> *patch_point_vals_side, unsigned int n_dofs_bulk,
116  unsigned int n_dofs_side, unsigned int op_idx_bulk, unsigned int op_idx_side, unsigned int join_idx)
117  : patch_point_vals_bulk_(patch_point_vals_bulk), patch_point_vals_side_(patch_point_vals_side),
118  n_dofs_high_(n_dofs_side), n_dofs_low_(n_dofs_bulk), op_idx_bulk_(op_idx_bulk), op_idx_side_(op_idx_side), join_idx_(join_idx) {
119  //ASSERT( (patch_point_vals_bulk->dim()==2) || (patch_point_vals_bulk->dim()==3) )(patch_point_vals_bulk->dim() ).error("Invalid dimension, must be 2 or 3!");
120  }
121 
122  /// Return global index of DOF
123  inline unsigned int join_idx() const {
124  return join_idx_;
125  }
126 
127  /// Return local index of DOF (on low / high-dim) - should be private method
128  inline unsigned int local_idx() const {
129  if (this->is_high_dim()) return (join_idx_ - n_dofs_low_);
130  else return join_idx_;
131  }
132 
133  inline unsigned int n_dofs_low() const {
134  return n_dofs_low_;
135  }
136 
137  inline unsigned int n_dofs_high() const {
138  return n_dofs_high_;
139  }
140 
141  inline unsigned int n_dofs_both() const {
142  return n_dofs_high_ + n_dofs_low_;
143  }
144 
145  inline bool is_high_dim() const {
146  return (join_idx_ >= n_dofs_low_);
147  }
148 
149  /// Iterates to next item.
150  inline void inc() {
151  join_idx_++;
152  }
153 
154  /// Comparison of accessors.
155  bool operator==(const JoinShapeAccessor<ValueType>& other) const {
156  return (join_idx_ == other.join_idx_);
157  }
158 
159 
160  ValueType operator()(const BulkPoint &point);
161 
162  ValueType operator()(const SidePoint &point);
163 
164 private:
165  // attributes:
166  PatchPointValues<3> *patch_point_vals_bulk_; ///< Pointer to bulk PatchPointValues
167  PatchPointValues<3> *patch_point_vals_side_; ///< Pointer to side PatchPointValues
168  unsigned int n_dofs_high_; ///< Number of DOFs on high-dim element
169  unsigned int n_dofs_low_; ///< Number of DOFs on low-dim element
170  unsigned int op_idx_bulk_; ///< Index of operation in patch_point_vals_bulk_.operations vector
171  unsigned int op_idx_side_; ///< Index of operation in patch_point_vals_side_.operations vector
172  unsigned int join_idx_; ///< Index of processed DOF
173 };
174 
175 
176 template<unsigned int dim>
178 {
179 protected:
180  // Default constructor
182  {}
183 
184  /// Return FiniteElement of \p component_idx for FESystem or \p fe for other types
185  template<unsigned int FE_dim>
186  std::shared_ptr<FiniteElement<FE_dim>> fe_comp(std::shared_ptr< FiniteElement<FE_dim> > fe, uint component_idx) {
187  FESystem<FE_dim> *fe_sys = dynamic_cast<FESystem<FE_dim>*>( fe.get() );
188  if (fe_sys != nullptr) {
189  return fe_sys->fe()[component_idx];
190  } else {
191  ASSERT_EQ(component_idx, 0).warning("Non-zero component_idx can only be used for FESystem.");
192  return fe;
193  }
194  }
195 
196  /**
197  * @brief Precomputed values of basis functions at the bulk quadrature points.
198  *
199  * Dimensions: (no. of quadrature points)
200  * x (no. of dofs)
201  * x (no. of components in ref. cell)
202  */
203  template<unsigned int FE_dim>
205  std::vector<std::vector<arma::vec> > ref_shape_vals( q->size(), vector<arma::vec>(fe->n_dofs()) );
206 
207  arma::mat shape_values(fe->n_dofs(), fe->n_components());
208  for (unsigned int i=0; i<q->size(); i++)
209  {
210  for (unsigned int j=0; j<fe->n_dofs(); j++)
211  {
212  for (unsigned int c=0; c<fe->n_components(); c++)
213  shape_values(j,c) = fe->shape_value(j, q->point<FE_dim>(i), c);
214 
215  ref_shape_vals[i][j] = trans(shape_values.row(j));
216  }
217  }
218 
219  return ref_shape_vals;
220  }
221 
222  /**
223  * @brief Precomputed values of basis functions at the side quadrature points.
224  *
225  * Dimensions: (sides)
226  * x (no. of quadrature points)
227  * x (no. of dofs)
228  * x (no. of components in ref. cell)
229  */
230  template<unsigned int FE_dim>
233 
234  arma::mat shape_values(fe->n_dofs(), fe->n_components());
235 
236  for (unsigned int sid=0; sid<FE_dim+1; sid++) {
237  auto quad = q->make_from_side<dim>(sid);
238  for (unsigned int i=0; i<quad.size(); i++)
239  {
240  for (unsigned int j=0; j<fe->n_dofs(); j++)
241  {
242  for (unsigned int c=0; c<fe->n_components(); c++) {
243  shape_values(j,c) = fe->shape_value(j, quad.template point<FE_dim>(i), c);
244  }
245 
246  ref_shape_vals[sid][i][j] = trans(shape_values.row(j));
247  }
248  }
249  }
250 
251  return ref_shape_vals;
252  }
253 
254 };
255 
256 template<unsigned int dim>
257 class BulkValues : public BaseValues<dim>
258 {
259 public:
260  /// Constructor
262  : BaseValues<dim>(), patch_point_vals_(patch_point_vals) {
263  ASSERT_EQ(patch_point_vals.dim(), dim);
264  fe_ = fe[Dim<dim>{}];
265  }
266 
267  /**
268  * @brief Register the product of Jacobian determinant and the quadrature
269  * weight at bulk quadrature points.
270  *
271  * @param quad Quadrature.
272  */
273  inline ElQ<Scalar> JxW()
274  {
276  }
277 
278  /// Create bulk accessor of coords entity
280  {
282  }
283 
284 // inline ElQ<Tensor> jacobian(std::initializer_list<Quadrature *> quad_list)
285 // {}
286 
287  /// Create bulk accessor of jac determinant entity
289  {
291  }
292 
293  /**
294  * @brief Return the value of the @p function_no-th shape function at
295  * the @p p bulk quadrature point.
296  *
297  * @param component_idx Number of the shape function.
298  */
299  inline FeQ<Scalar> scalar_shape(uint component_idx = 0)
300  {
301  auto fe_component = this->fe_comp(fe_, component_idx);
302  ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
303 
304  // use lambda reinit function
305  std::vector< std::vector<double> > shape_values( fe_component->n_dofs(), vector<double>(patch_point_vals_.get_quadrature()->size()) );
306  auto ref_shape_vals = this->ref_shape_values_bulk(patch_point_vals_.get_quadrature(), fe_component);
307  for (unsigned int i = 0; i < patch_point_vals_.get_quadrature()->size(); i++)
308  for (unsigned int j = 0; j < fe_component->n_dofs(); j++) {
309  shape_values[j][i] = ref_shape_vals[i][j][0];
310  }
311  uint scalar_shape_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
312  auto lambda_scalar_shape = [shape_values, scalar_shape_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
313  bulk_reinit::ptop_scalar_shape(operations, shape_values, scalar_shape_op_idx);
314  };
315  patch_point_vals_.make_fe_op({1}, lambda_scalar_shape, {}, fe_component->n_dofs());
316  uint op_idx = patch_point_vals_.operations_.size()-1;
317 
318  return FeQ<Scalar>(patch_point_vals_, op_idx, fe_component->n_dofs());
319  }
320 
321 // inline FeQ<Vector> vector_shape(uint component_idx = 0)
322 // {}
323 
324 // inline FeQ<Tensor> tensor_shape(uint component_idx = 0)
325 // {}
326 
327  /**
328  * @brief Return the value of the @p function_no-th gradient shape function at
329  * the @p p bulk quadrature point.
330  *
331  * @param component_idx Number of the shape function.
332  */
333  inline FeQ<Vector> grad_scalar_shape(uint component_idx=0)
334  {
335  auto fe_component = this->fe_comp(fe_, component_idx);
336  ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of grad_scalar_shape accessor must be FEScalar!\n");
337 
338  // use lambda reinit function
339  auto ref_shape_grads = this->ref_shape_gradients(fe_component);
340  uint scalar_shape_grads_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
341  auto lambda_scalar_shape_grad = [ref_shape_grads, scalar_shape_grads_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
342  bulk_reinit::ptop_scalar_shape_grads<dim>(operations, ref_shape_grads, scalar_shape_grads_op_idx);
343  };
344  patch_point_vals_.make_fe_op({3}, lambda_scalar_shape_grad, {FeBulk::BulkOps::opInvJac}, fe_component->n_dofs());
345  uint op_idx = patch_point_vals_.operations_.size()-1;
346 
347  return FeQ<Vector>(patch_point_vals_, op_idx, fe_component->n_dofs());
348  }
349 
350 // inline FeQ<Tensor> grad_vector_shape(std::initializer_list<Quadrature *> quad_list, unsigned int i_comp=0)
351 // {}
352 
353 private:
354  /**
355  * @brief Precomputed gradients of basis functions at the quadrature points.
356  *
357  * Dimensions: (no. of quadrature points)
358  * x (no. of dofs)
359  * x ((dim_ of. ref. cell)x(no. of components in ref. cell))
360  */
363  std::vector<std::vector<arma::mat> > ref_shape_grads( q->size(), vector<arma::mat>(fe->n_dofs()) );
364 
365  arma::mat grad(dim, fe->n_components());
366  for (unsigned int i_pt=0; i_pt<q->size(); i_pt++)
367  {
368  for (unsigned int i_dof=0; i_dof<fe->n_dofs(); i_dof++)
369  {
370  grad.zeros();
371  for (unsigned int c=0; c<fe->n_components(); c++)
372  grad.col(c) += fe->shape_grad(i_dof, q->point<dim>(i_pt), c);
373 
374  ref_shape_grads[i_pt][i_dof] = grad;
375  }
376  }
377 
378  return ref_shape_grads;
379  }
380 
382  std::shared_ptr< FiniteElement<dim> > fe_;
383 };
384 
385 
386 template<unsigned int dim>
387 class SideValues : public BaseValues<dim>
388 {
389 public:
390  /// Constructor
392  : BaseValues<dim>(), patch_point_vals_(patch_point_vals) {
393  ASSERT_EQ(patch_point_vals.dim(), dim);
394  fe_ = fe[Dim<dim>{}];
395  }
396 
397  /// Same as BulkValues::JxW but register at side quadrature points.
398  inline ElQ<Scalar> JxW()
399  {
401  }
402 
403  /**
404  * @brief Register the normal vector to a side at side quadrature points.
405  *
406  * @param quad Quadrature.
407  */
409  {
411  }
412 
413  /// Create side accessor of coords entity
415  {
417  }
418 
419  /// Create bulk accessor of jac determinant entity
421  {
423  }
424 
425  /// Same as BulkValues::scalar_shape but register at side quadrature points.
426  inline FeQ<Scalar> scalar_shape(uint component_idx = 0)
427  {
428  auto fe_component = this->fe_comp(fe_, component_idx);
429  ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
430 
431  // use lambda reinit function
433  dim+1,
435  );
436  auto ref_shape_vals = this->ref_shape_values_side(patch_point_vals_.get_quadrature(), fe_component);
437  for (unsigned int s=0; s<dim+1; ++s)
438  for (unsigned int i = 0; i < patch_point_vals_.get_quadrature()->size(); i++)
439  for (unsigned int j = 0; j < fe_component->n_dofs(); j++) {
440  shape_values[s][i][j] = ref_shape_vals[s][i][j][0];
441  }
442  uint scalar_shape_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
443  auto lambda_scalar_shape = [shape_values, scalar_shape_op_idx](std::vector<ElOp<3>> &operations, IntTableArena &el_table) {
444  side_reinit::ptop_scalar_shape(operations, el_table, shape_values, scalar_shape_op_idx);
445  };
446  patch_point_vals_.make_fe_op({1}, lambda_scalar_shape, {}, fe_component->n_dofs());
447  uint op_idx = patch_point_vals_.operations_.size()-1;
448 
449  return FeQ<Scalar>(patch_point_vals_, op_idx, fe_component->n_dofs());
450  }
451 
452  /// Same as BulkValues::grad_scalar_shape but register at side quadrature points.
453  inline FeQ<Vector> grad_scalar_shape(uint component_idx=0)
454  {
455  auto fe_component = this->fe_comp(fe_, component_idx);
456  ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of grad_scalar_shape accessor must be FEScalar!\n");
457 
458  // use lambda reinit function
459  auto ref_shape_grads = this->ref_shape_gradients(fe_component);
460  uint scalar_shape_grads_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
461  auto lambda_scalar_shape_grad = [ref_shape_grads, scalar_shape_grads_op_idx](std::vector<ElOp<3>> &operations, IntTableArena &el_table) {
462  side_reinit::ptop_scalar_shape_grads<dim>(operations, el_table, ref_shape_grads, scalar_shape_grads_op_idx);
463  };
464  patch_point_vals_.make_fe_op({3}, lambda_scalar_shape_grad, {FeSide::SideOps::opElInvJac}, fe_component->n_dofs());
465  uint op_idx = patch_point_vals_.operations_.size()-1;
466 
467  return FeQ<Vector>(patch_point_vals_, op_idx, fe_component->n_dofs());
468  }
469 
470 private:
471  /**
472  * @brief Precomputed gradients of basis functions at the quadrature points.
473  *
474  * Dimensions: (sides)
475  * x (no. of quadrature points)
476  * x (no. of dofs)
477  * x ((dim_ of. ref. cell)x(no. of components in ref. cell))
478  */
482 
483  arma::mat grad(dim, fe->n_components());
484  for (unsigned int sid=0; sid<dim+1; sid++) {
485  auto quad = q->make_from_side<dim>(sid);
486  for (unsigned int i_pt=0; i_pt<quad.size(); i_pt++)
487  {
488  for (unsigned int i_dof=0; i_dof<fe->n_dofs(); i_dof++)
489  {
490  grad.zeros();
491  for (unsigned int c=0; c<fe->n_components(); c++)
492  grad.col(c) += fe->shape_grad(i_dof, quad.template point<dim>(i_pt), c);
493 
494  ref_shape_grads[sid][i_pt][i_dof] = grad;
495  }
496  }
497  }
498 
499  return ref_shape_grads;
500  }
501 
503  std::shared_ptr< FiniteElement<dim> > fe_;
504 };
505 
506 
507 template<unsigned int dim>
508 class JoinValues : public BaseValues<dim>
509 {
510 public:
511  /// Constructor
512  JoinValues(PatchPointValues<3> *patch_point_vals_bulk, PatchPointValues<3> *patch_point_vals_side, MixedPtr<FiniteElement> fe)
513  : BaseValues<dim>(), patch_point_vals_bulk_(patch_point_vals_bulk), patch_point_vals_side_(patch_point_vals_side) {
514  ASSERT_EQ(patch_point_vals_bulk->dim(), dim-1);
515  ASSERT_EQ(patch_point_vals_side->dim(), dim);
516  fe_high_dim_ = fe[Dim<dim>{}];
517  fe_low_dim_ = fe[Dim<dim-1>{}];
518  }
519 
521  {
522  // element of lower dim (bulk points)
523  auto fe_component_low = this->fe_comp(fe_low_dim_, component_idx);
524  ASSERT_EQ(fe_component_low->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
525  // use lambda reinit function
526  std::vector< std::vector<double> > shape_values_bulk( patch_point_vals_bulk_->get_quadrature()->size(), vector<double>(fe_component_low->n_dofs()) );
527  auto ref_shape_vals_bulk = this->ref_shape_values_bulk(patch_point_vals_bulk_->get_quadrature(), fe_component_low);
528  for (unsigned int i = 0; i < patch_point_vals_bulk_->get_quadrature()->size(); i++)
529  for (unsigned int j = 0; j < fe_component_low->n_dofs(); j++) {
530  shape_values_bulk[i][j] = ref_shape_vals_bulk[i][j][0];
531  }
532  uint scalar_shape_op_idx_bulk = patch_point_vals_bulk_->operations_.size(); // index in operations_ vector
533  auto lambda_scalar_shape_bulk = [shape_values_bulk, scalar_shape_op_idx_bulk](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
534  bulk_reinit::ptop_scalar_shape(operations, shape_values_bulk, scalar_shape_op_idx_bulk);
535  };
536  patch_point_vals_bulk_->make_fe_op({1}, lambda_scalar_shape_bulk, {}, fe_component_low->n_dofs());
537  uint op_idx_bulk = patch_point_vals_bulk_->operations_.size()-1;
538 
539  // element of higher dim (side points)
540  auto fe_component_high = this->fe_comp(fe_high_dim_, component_idx);
541  ASSERT_EQ(fe_component_high->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
542  // use lambda reinit function
543  std::vector< std::vector< std::vector<double> > > shape_values_side(
544  dim+1,
546  );
547  auto ref_shape_vals_side = this->ref_shape_values_side(patch_point_vals_side_->get_quadrature(), fe_component_high);
548  for (unsigned int s=0; s<dim+1; ++s)
549  for (unsigned int i = 0; i < patch_point_vals_side_->get_quadrature()->size(); i++)
550  for (unsigned int j = 0; j < fe_component_high->n_dofs(); j++) {
551  shape_values_side[s][i][j] = ref_shape_vals_side[s][i][j][0];
552  }
553  uint scalar_shape_op_idx_side = patch_point_vals_side_->operations_.size(); // index in operations_ vector
554  auto lambda_scalar_shape_side = [shape_values_side, scalar_shape_op_idx_side](std::vector<ElOp<3>> &operations, IntTableArena &el_table) {
555  side_reinit::ptop_scalar_shape(operations, el_table, shape_values_side, scalar_shape_op_idx_side);
556  };
557  patch_point_vals_side_->make_fe_op({1}, lambda_scalar_shape_side, {}, fe_component_high->n_dofs());
558  uint op_idx_side = patch_point_vals_side_->operations_.size()-1;
559 
560  auto bgn_it = make_iter<JoinShapeAccessor<Scalar>>( JoinShapeAccessor<Scalar>(patch_point_vals_bulk_, patch_point_vals_side_,
561  fe_component_low->n_dofs(), fe_component_high->n_dofs(), op_idx_bulk, op_idx_side, 0) );
562  unsigned int end_idx = fe_component_low->n_dofs() + fe_component_high->n_dofs();
563  auto end_it = make_iter<JoinShapeAccessor<Scalar>>( JoinShapeAccessor<Scalar>(patch_point_vals_bulk_, patch_point_vals_side_,
564  fe_component_low->n_dofs(), fe_component_high->n_dofs(), op_idx_bulk, op_idx_side, end_idx) );
565  return Range<JoinShapeAccessor<Scalar>>(bgn_it, end_it);
566  }
567 
568 private:
571  std::shared_ptr< FiniteElement<dim> > fe_high_dim_;
572  std::shared_ptr< FiniteElement<dim-1> > fe_low_dim_;
573 };
574 
575 /// Template specialization of dim = 1
576 template <>
577 class JoinValues<1> : public BaseValues<1>
578 {
579 public:
580  /// Constructor
582  : BaseValues<1>() {}
583 
585  {
587  make_iter<JoinShapeAccessor<Scalar>>(JoinShapeAccessor<Scalar>()),
589  }
590 };
591 
592 
593 template<unsigned int spacedim = 3>
595 public:
596  /// Struct for pre-computing number of elements, sides, bulk points and side points on each dimension.
597  struct TableSizes {
598  public:
599  /// Constructor
603  }
604 
605  /// Set all values to zero
606  void reset() {
607  std::fill(elem_sizes_[0].begin(), elem_sizes_[0].end(), 0);
608  std::fill(elem_sizes_[1].begin(), elem_sizes_[1].end(), 0);
609  std::fill(point_sizes_[0].begin(), point_sizes_[0].end(), 0);
610  std::fill(point_sizes_[1].begin(), point_sizes_[1].end(), 0);
611  }
612 
613  /**
614  * Holds number of elements and sides on each dimension
615  * Format:
616  * { {n_elements_1D, n_elements_2D, n_elements_3D },
617  * {n_sides_1D, n_sides_2D, n_sides_3D } }
618  */
620 
621  /**
622  * Holds number of bulk and side points on each dimension
623  * Format:
624  * { {n_bulk_points_1D, n_bulk_points_2D, n_bulk_points_3D },
625  * {n_side_points_1D, n_side_points_2D, n_side_points_3D } }
626  */
628  };
629 
631  : asm_arena_(1024 * 1024, 256),
632  patch_arena_(nullptr),
636  patch_point_vals_side_{ {FeSide::PatchPointValues(1, 0, asm_arena_),
637  FeSide::PatchPointValues(2, 0, asm_arena_),
638  FeSide::PatchPointValues(3, 0, asm_arena_)} }
639  {
640  used_quads_[0] = false; used_quads_[1] = false;
641  }
642 
643  PatchFEValues(unsigned int quad_order, MixedPtr<FiniteElement> fe)
644  : asm_arena_(1024 * 1024, 256),
645  patch_arena_(nullptr),
646  patch_point_vals_bulk_{ {FeBulk::PatchPointValues(1, quad_order, asm_arena_),
647  FeBulk::PatchPointValues(2, quad_order, asm_arena_),
648  FeBulk::PatchPointValues(3, quad_order, asm_arena_)} },
649  patch_point_vals_side_{ {FeSide::PatchPointValues(1, quad_order, asm_arena_),
650  FeSide::PatchPointValues(2, quad_order, asm_arena_),
651  FeSide::PatchPointValues(3, quad_order, asm_arena_)} },
652  fe_(fe)
653  {
654  used_quads_[0] = false; used_quads_[1] = false;
655  }
656 
657 
658  /// Destructor
660  {
661  if (patch_arena_!=nullptr)
662  delete patch_arena_;
663  }
664 
665  /// Return bulk or side quadrature of given dimension
666  Quadrature *get_quadrature(uint dim, bool is_bulk) const {
667  if (is_bulk) return patch_point_vals_bulk_[dim-1].get_quadrature();
668  else return patch_point_vals_side_[dim-1].get_quadrature();
669  }
670 
671  /**
672  * @brief Initialize structures and calculates cell-independent data.
673  *
674  * @param _quadrature The quadrature rule for the cell associated
675  * to given finite element or for the cell side.
676  * @param _flags The update flags.
677  */
678  template<unsigned int DIM>
679  void initialize(Quadrature &_quadrature)
680  {
681  if ( _quadrature.dim() == DIM ) {
682  used_quads_[0] = true;
683  patch_point_vals_bulk_[DIM-1].initialize(); // bulk
684  } else {
685  used_quads_[1] = true;
686  patch_point_vals_side_[DIM-1].initialize(); // side
687  }
688  }
689 
690  /// Finalize initialization, creates child (patch) arena and passes it to PatchPointValue objects
691  void init_finalize() {
692  patch_arena_ = asm_arena_.get_child_arena();
693  for (unsigned int i=0; i<3; ++i) {
694  if (used_quads_[0]) patch_point_vals_bulk_[i].init_finalize(patch_arena_);
695  if (used_quads_[1]) patch_point_vals_side_[i].init_finalize(patch_arena_);
696  }
697  }
698 
699  /// Reset PatchpointValues structures
700  void reset()
701  {
702  for (unsigned int i=0; i<3; ++i) {
703  if (used_quads_[0]) patch_point_vals_bulk_[i].reset();
704  if (used_quads_[1]) patch_point_vals_side_[i].reset();
705  }
706  patch_arena_->reset();
707  }
708 
709  /// Reinit data.
711  {
712  for (unsigned int i=0; i<3; ++i) {
713  if (used_quads_[0]) patch_point_vals_bulk_[i].reinit_patch();
714  if (used_quads_[1]) patch_point_vals_side_[i].reinit_patch();
715  }
716  }
717 
718  /**
719  * @brief Returns the number of shape functions.
720  */
721  template<unsigned int dim>
722  inline unsigned int n_dofs() const {
723  ASSERT((dim>=0) && (dim<=3))(dim).error("Dimension must be 0, 1, 2 or 3.");
724  return fe_[Dim<dim>{}]->n_dofs();
725  }
726 
727  /// Return BulkValue object of dimension given by template parameter
728  template<unsigned int dim>
730  ASSERT((dim>0) && (dim<=3))(dim).error("Dimension must be 1, 2 or 3.");
731  return BulkValues<dim>(patch_point_vals_bulk_[dim-1], fe_);
732  }
733 
734  /// Return SideValue object of dimension given by template parameter
735  template<unsigned int dim>
737  ASSERT((dim>0) && (dim<=3))(dim).error("Dimension must be 1, 2 or 3.");
738  return SideValues<dim>(patch_point_vals_side_[dim-1], fe_);
739  }
740 
741  /// Return JoinValue object of dimension given by template parameter
742  template<unsigned int dim>
744  //ASSERT((dim>1) && (dim<=3))(dim).error("Dimension must be 2 or 3.");
745  return JoinValues<dim>(&patch_point_vals_bulk_[dim-2], &patch_point_vals_side_[dim-1], fe_);
746  }
747 
748  /** Following methods are used during update of patch. **/
749 
750  /// Resize tables of patch_point_vals_
751  void resize_tables(TableSizes table_sizes) {
752  for (uint i=0; i<3; ++i) {
753  if (used_quads_[0]) patch_point_vals_bulk_[i].resize_tables(table_sizes.elem_sizes_[0][i], table_sizes.point_sizes_[0][i]);
754  if (used_quads_[1]) patch_point_vals_side_[i].resize_tables(table_sizes.elem_sizes_[1][i], table_sizes.point_sizes_[1][i]);
755  }
756  }
757 
758  /// Register element to patch_point_vals_ table by dimension of element
759  uint register_element(DHCellAccessor cell, uint element_patch_idx) {
760  arma::mat coords;
761  switch (cell.dim()) {
762  case 1:
763  coords = MappingP1<1,spacedim>::element_map(cell.elm());
764  return patch_point_vals_bulk_[0].register_element(coords, element_patch_idx);
765  break;
766  case 2:
767  coords = MappingP1<2,spacedim>::element_map(cell.elm());
768  return patch_point_vals_bulk_[1].register_element(coords, element_patch_idx);
769  break;
770  case 3:
771  coords = MappingP1<3,spacedim>::element_map(cell.elm());
772  return patch_point_vals_bulk_[2].register_element(coords, element_patch_idx);
773  break;
774  default:
775  ASSERT(false);
776  return 0;
777  break;
778  }
779  }
780 
781  /// Register side to patch_point_vals_ table by dimension of side
783  arma::mat side_coords(spacedim, cell_side.dim());
784  for (unsigned int n=0; n<cell_side.dim(); n++)
785  for (unsigned int c=0; c<spacedim; c++)
786  side_coords(c,n) = (*cell_side.side().node(n))[c];
787 
788  arma::mat elm_coords;
789  DHCellAccessor cell = cell_side.cell();
790  switch (cell.dim()) {
791  case 1:
792  elm_coords = MappingP1<1,spacedim>::element_map(cell.elm());
793  return patch_point_vals_side_[0].register_side(elm_coords, side_coords, cell_side.side_idx());
794  break;
795  case 2:
796  elm_coords = MappingP1<2,spacedim>::element_map(cell.elm());
797  return patch_point_vals_side_[1].register_side(elm_coords, side_coords, cell_side.side_idx());
798  break;
799  case 3:
800  elm_coords = MappingP1<3,spacedim>::element_map(cell.elm());
801  return patch_point_vals_side_[2].register_side(elm_coords, side_coords, cell_side.side_idx());
802  break;
803  default:
804  ASSERT(false);
805  return 0;
806  break;
807  }
808  }
809 
810  /// Register bulk point to patch_point_vals_ table by dimension of element
811  uint register_bulk_point(DHCellAccessor cell, uint elem_table_row, uint value_patch_idx, uint i_point_on_elem) {
812  return patch_point_vals_bulk_[cell.dim()-1].register_bulk_point(elem_table_row, value_patch_idx, cell.elm_idx(), i_point_on_elem);
813  }
814 
815  /// Register side point to patch_point_vals_ table by dimension of side
816  uint register_side_point(DHCellSide cell_side, uint elem_table_row, uint value_patch_idx, uint i_point_on_side) {
817  return patch_point_vals_side_[cell_side.dim()-1].register_side_point(elem_table_row, value_patch_idx, cell_side.elem_idx(),
818  cell_side.side_idx(), i_point_on_side);
819  }
820 
821  /// Temporary development method
822  void print_data_tables(ostream& stream, bool points, bool ints, bool only_bulk=true) const {
823  stream << endl << "Table of patch FE data:" << endl;
824  for (uint i=0; i<3; ++i) {
825  stream << std::setfill('-') << setw(100) << "" << endl;
826  stream << "Bulk, dimension " << (i+1) << endl;
827  patch_point_vals_bulk_[i].print_data_tables(stream, points, ints);
828  }
829  if (!only_bulk)
830  for (uint i=0; i<3; ++i) {
831  stream << std::setfill('-') << setw(100) << "" << endl;
832  stream << "Side, dimension " << (i+1) << endl;
833  patch_point_vals_side_[i].print_data_tables(stream, points, ints);
834  }
835  stream << std::setfill('=') << setw(100) << "" << endl;
836  }
837 
838  /// Temporary development method
839  void print_operations(ostream& stream) const {
840  stream << endl << "Table of patch FE operations:" << endl;
841  for (uint i=0; i<3; ++i) {
842  stream << std::setfill('-') << setw(100) << "" << endl;
843  stream << "Bulk, dimension " << (i+1) << endl;
844  patch_point_vals_bulk_[i].print_operations(stream, 0);
845  }
846  for (uint i=0; i<3; ++i) {
847  stream << std::setfill('-') << setw(100) << "" << endl;
848  stream << "Side, dimension " << (i+1) << endl;
849  patch_point_vals_side_[i].print_operations(stream, 1);
850  }
851  stream << std::setfill('=') << setw(100) << "" << endl;
852  }
853 
854 private:
857  std::array<FeBulk::PatchPointValues<spacedim>, 3> patch_point_vals_bulk_; ///< Sub objects of bulk data of dimensions 1,2,3
858  std::array<FeSide::PatchPointValues<spacedim>, 3> patch_point_vals_side_; ///< Sub objects of side data of dimensions 1,2,3
859 
860  MixedPtr<FiniteElement> fe_; ///< Mixed of shared pointers of FiniteElement object
861  bool used_quads_[2]; ///< Pair of flags signs holds info if bulk and side quadratures are used
862 
863  template <class ValueType>
864  friend class ElQ;
865  template <class ValueType>
866  friend class FeQ;
867  template <class ValueType>
868  friend class JoinShapeAccessor;
869 };
870 
871 
872 template <class ValueType>
873 ValueType ElQ<ValueType>::operator()(const BulkPoint &point) {
874  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
875  return patch_point_vals_.scalar_value(op_idx_, value_cache_idx);
876 }
877 
878 template <>
880  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
881  return patch_point_vals_.vector_value(op_idx_, value_cache_idx);
882 }
883 
884 template <>
886  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
887  return patch_point_vals_.tensor_value(op_idx_, value_cache_idx);
888 }
889 
890 template <class ValueType>
891 ValueType ElQ<ValueType>::operator()(const SidePoint &point) {
892  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
893  return patch_point_vals_.scalar_value(op_idx_, value_cache_idx);
894 }
895 
896 template <>
898  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
899  return patch_point_vals_.vector_value(op_idx_, value_cache_idx);
900 }
901 
902 template <>
904  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
905  return patch_point_vals_.tensor_value(op_idx_, value_cache_idx);
906 }
907 
908 template <class ValueType>
909 ValueType FeQ<ValueType>::operator()(unsigned int shape_idx, const BulkPoint &point) {
910  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
911  return patch_point_vals_.scalar_value(op_idx_, value_cache_idx, shape_idx);
912 }
913 
914 template <>
915 inline Vector FeQ<Vector>::operator()(unsigned int shape_idx, const BulkPoint &point) {
916  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
917  return patch_point_vals_.vector_value(op_idx_, value_cache_idx, shape_idx);
918 }
919 
920 template <>
921 inline Tensor FeQ<Tensor>::operator()(unsigned int shape_idx, const BulkPoint &point) {
922  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
923  return patch_point_vals_.tensor_value(op_idx_, value_cache_idx, shape_idx);
924 }
925 
926 template <class ValueType>
927 ValueType FeQ<ValueType>::operator()(unsigned int shape_idx, const SidePoint &point) {
928  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
929  return patch_point_vals_.scalar_value(op_idx_, value_cache_idx, shape_idx);
930 }
931 
932 template <>
933 inline Vector FeQ<Vector>::operator()(unsigned int shape_idx, const SidePoint &point) {
934  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
935  return patch_point_vals_.vector_value(op_idx_, value_cache_idx, shape_idx);
936 }
937 
938 template <>
939 inline Tensor FeQ<Tensor>::operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const SidePoint &point) {
940  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
941  return patch_point_vals_.tensor_value(op_idx_, value_cache_idx, shape_idx);
942 }
943 
944 
945 template <class ValueType>
947  if (this->is_high_dim()) {
948  return 0.0;
949  } else {
950  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
951  return patch_point_vals_bulk_->scalar_value(op_idx_bulk_, value_cache_idx, this->local_idx());
952  }
953 }
954 
955 template <>
957  Vector vect; vect.zeros();
958  return vect;
959 }
960 
961 template <>
963  Tensor tens; tens.zeros();
964  return tens;
965 }
966 
967 template <class ValueType>
969  if (this->is_high_dim()) {
970  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
971  return patch_point_vals_side_->scalar_value(op_idx_side_, value_cache_idx, this->local_idx());
972  } else {
973  return 0.0;
974  }
975 }
976 
977 template <>
979  Vector vect; vect.zeros();
980  return vect;
981 }
982 
983 template <>
985  Tensor tens; tens.zeros();
986  return tens;
987 }
988 
989 
990 #endif /* PATCH_FE_VALUES_HH_ */
#define ASSERT(expr)
Definition: asserts.hh:351
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:333
std::vector< std::vector< std::vector< arma::vec > > > ref_shape_values_side(Quadrature *q, std::shared_ptr< FiniteElement< FE_dim >> fe)
Precomputed values of basis functions at the side quadrature points.
std::shared_ptr< FiniteElement< FE_dim > > fe_comp(std::shared_ptr< FiniteElement< FE_dim > > fe, uint component_idx)
Return FiniteElement of component_idx for FESystem or fe for other types.
std::vector< std::vector< arma::vec > > ref_shape_values_bulk(Quadrature *q, std::shared_ptr< FiniteElement< FE_dim >> fe)
Precomputed values of basis functions at the bulk quadrature points.
Base point accessor class.
Definition: eval_subset.hh:55
const ElementCacheMap * elm_cache_map() const
Definition: eval_subset.hh:73
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
ElQ< Scalar > JxW()
Register the product of Jacobian determinant and the quadrature weight at bulk quadrature points.
PatchPointValues< 3 > & patch_point_vals_
ElQ< Vector > coords()
Create bulk accessor of coords entity.
std::vector< std::vector< arma::mat > > ref_shape_gradients(std::shared_ptr< FiniteElement< dim >> fe)
Precomputed gradients of basis functions at the quadrature points.
ElQ< Scalar > determinant()
Create bulk accessor of jac determinant entity.
std::shared_ptr< FiniteElement< dim > > fe_
FeQ< Vector > grad_scalar_shape(uint component_idx=0)
Return the value of the function_no-th gradient shape function at the p bulk quadrature point.
FeQ< Scalar > scalar_shape(uint component_idx=0)
Return the value of the function_no-th shape function at the p bulk quadrature point.
BulkValues(PatchPointValues< 3 > &patch_point_vals, MixedPtr< FiniteElement > fe)
Constructor.
Cell accessor allow iterate over DOF handler cells.
unsigned int dim() const
Return dimension of element appropriate to cell.
unsigned int elm_idx() const
Return serial idx to element of loc_ele_idx_.
ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
Side accessor allows to iterate over sides of DOF handler cell.
unsigned int elem_idx() const
Side side() const
Return Side of given cell and side_idx.
const DHCellAccessor & cell() const
Return DHCellAccessor appropriate to the side.
unsigned int dim() const
Return dimension of element appropriate to the side.
unsigned int side_idx() const
ElQ(PatchPointValues< 3 > &patch_point_vals, unsigned int op_idx)
Constructor.
unsigned int op_idx_
Index of operation in patch_point_vals_.operations vector.
ElQ()=delete
Forbidden default constructor.
PatchPointValues< 3 > & patch_point_vals_
Reference to PatchPointValues.
ValueType operator()(FMT_UNUSED const SidePoint &point)
ValueType operator()(FMT_UNUSED const BulkPoint &point)
int element_eval_point(unsigned int i_elem_in_cache, unsigned int i_eval_point) const
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:102
const std::vector< std::shared_ptr< FiniteElement< dim > > > & fe() const
Definition: fe_system.hh:147
Bulk data specialization, order of item in operations_ vector corresponds to the BulkOps enum.
PatchPointValues< 3 > & patch_point_vals_
Reference to PatchPointValues.
FeQ(PatchPointValues< 3 > &patch_point_vals, unsigned int op_idx, unsigned int n_dofs)
unsigned int n_dofs_
Number of DOFs.
ValueType operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const BulkPoint &point)
unsigned int op_idx_
Index of operation in patch_point_vals_.operations vector.
FeQ()=delete
Forbidden default constructor.
ValueType operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const SidePoint &point)
Bulk Side specialization, order of item in operations_ vector corresponds to the SideOps enum.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
unsigned int n_dofs_both() const
unsigned int local_idx() const
Return local index of DOF (on low / high-dim) - should be private method.
bool operator==(const JoinShapeAccessor< ValueType > &other) const
Comparison of accessors.
JoinShapeAccessor(PatchPointValues< 3 > *patch_point_vals_bulk, PatchPointValues< 3 > *patch_point_vals_side, unsigned int n_dofs_bulk, unsigned int n_dofs_side, unsigned int op_idx_bulk, unsigned int op_idx_side, unsigned int join_idx)
bool is_high_dim() const
PatchPointValues< 3 > * patch_point_vals_bulk_
Pointer to bulk PatchPointValues.
unsigned int op_idx_side_
Index of operation in patch_point_vals_side_.operations vector.
unsigned int op_idx_bulk_
Index of operation in patch_point_vals_bulk_.operations vector.
unsigned int n_dofs_high() const
JoinShapeAccessor()
Default constructor.
void inc()
Iterates to next item.
unsigned int join_idx_
Index of processed DOF.
unsigned int n_dofs_high_
Number of DOFs on high-dim element.
unsigned int join_idx() const
Return global index of DOF.
unsigned int n_dofs_low() const
ValueType operator()(const BulkPoint &point)
unsigned int n_dofs_low_
Number of DOFs on low-dim element.
PatchPointValues< 3 > * patch_point_vals_side_
Pointer to side PatchPointValues.
JoinValues(FMT_UNUSED PatchPointValues< 3 > *patch_point_vals_bulk, FMT_UNUSED PatchPointValues< 3 > *patch_point_vals_side, FMT_UNUSED MixedPtr< FiniteElement > fe)
Constructor.
Range< JoinShapeAccessor< Scalar > > scalar_join_shape(FMT_UNUSED uint component_idx=0)
std::shared_ptr< FiniteElement< dim > > fe_high_dim_
JoinValues(PatchPointValues< 3 > *patch_point_vals_bulk, PatchPointValues< 3 > *patch_point_vals_side, MixedPtr< FiniteElement > fe)
Constructor.
PatchPointValues< 3 > * patch_point_vals_bulk_
Range< JoinShapeAccessor< Scalar > > scalar_join_shape(uint component_idx=0)
std::shared_ptr< FiniteElement< dim-1 > > fe_low_dim_
PatchPointValues< 3 > * patch_point_vals_side_
static ElementMap element_map(ElementAccessor< 3 > elm)
Definition: mapping_p1.cc:48
uint register_element(DHCellAccessor cell, uint element_patch_idx)
Register element to patch_point_vals_ table by dimension of element.
PatchArena * patch_arena_
uint register_bulk_point(DHCellAccessor cell, uint elem_table_row, uint value_patch_idx, uint i_point_on_elem)
Register bulk point to patch_point_vals_ table by dimension of element.
~PatchFEValues()
Destructor.
std::array< FeBulk::PatchPointValues< spacedim >, 3 > patch_point_vals_bulk_
Sub objects of bulk data of dimensions 1,2,3.
uint register_side(DHCellSide cell_side)
Register side to patch_point_vals_ table by dimension of side.
void init_finalize()
Finalize initialization, creates child (patch) arena and passes it to PatchPointValue objects.
void print_data_tables(ostream &stream, bool points, bool ints, bool only_bulk=true) const
Temporary development method.
JoinValues< dim > join_values()
Return JoinValue object of dimension given by template parameter.
Quadrature * get_quadrature(uint dim, bool is_bulk) const
Return bulk or side quadrature of given dimension.
BulkValues< dim > bulk_values()
Return BulkValue object of dimension given by template parameter.
uint register_side_point(DHCellSide cell_side, uint elem_table_row, uint value_patch_idx, uint i_point_on_side)
Register side point to patch_point_vals_ table by dimension of side.
void initialize(Quadrature &_quadrature)
Initialize structures and calculates cell-independent data.
void print_operations(ostream &stream) const
Temporary development method.
std::array< FeSide::PatchPointValues< spacedim >, 3 > patch_point_vals_side_
Sub objects of side data of dimensions 1,2,3.
MixedPtr< FiniteElement > fe_
Mixed of shared pointers of FiniteElement object.
SideValues< dim > side_values()
Return SideValue object of dimension given by template parameter.
AssemblyArena asm_arena_
void resize_tables(TableSizes table_sizes)
Resize tables of patch_point_vals_.
PatchFEValues(unsigned int quad_order, MixedPtr< FiniteElement > fe)
void reinit_patch()
Reinit data.
unsigned int n_dofs() const
Returns the number of shape functions.
void reset()
Reset PatchpointValues structures.
ElOp< spacedim > & make_fe_op(std::initializer_list< uint > shape, ReinitFunction reinit_f, std::vector< uint > input_ops_vec, uint n_dofs, OpSizeType size_type=pointOp)
std::vector< ElOp< spacedim > > operations_
Vector of all defined operations.
Quadrature * get_quadrature() const
Getter for quadrature.
uint dim() const
Getter for dim_.
Scalar scalar_value(uint op_idx, uint point_idx, uint i_dof=0) const
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
Quadrature make_from_side(unsigned int sid) const
Definition: quadrature.cc:47
unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:86
unsigned int dim() const
Definition: quadrature.hh:72
Armor::ArmaVec< double, point_dim > point(unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:151
Range helper class.
General point a+ side_begin_ + ccessor allow iterate over quadrature points of given side defined in ...
Definition: eval_subset.hh:116
unsigned int eval_point_idx() const
Return index in EvalPoints object.
Definition: eval_subset.hh:143
SideValues(PatchPointValues< 3 > &patch_point_vals, MixedPtr< FiniteElement > fe)
Constructor.
ElQ< Vector > normal_vector()
Register the normal vector to a side at side quadrature points.
std::shared_ptr< FiniteElement< dim > > fe_
ElQ< Scalar > JxW()
Same as BulkValues::JxW but register at side quadrature points.
ElQ< Vector > coords()
Create side accessor of coords entity.
FeQ< Scalar > scalar_shape(uint component_idx=0)
Same as BulkValues::scalar_shape but register at side quadrature points.
FeQ< Vector > grad_scalar_shape(uint component_idx=0)
Same as BulkValues::grad_scalar_shape but register at side quadrature points.
std::vector< std::vector< std::vector< arma::mat > > > ref_shape_gradients(std::shared_ptr< FiniteElement< dim >> fe)
Precomputed gradients of basis functions at the quadrature points.
PatchPointValues< 3 > & patch_point_vals_
ElQ< Scalar > determinant()
Create bulk accessor of jac determinant entity.
NodeAccessor< 3 > node(unsigned int i) const
Returns node for given local index i on the side.
Store finite element data on the actual patch such as shape function values, gradients,...
Eigen::Vector< ArenaVec< uint >, Eigen::Dynamic > IntTableArena
Definition: eigen_tools.hh:39
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
Class FESystem for compound finite elements.
Class FEValues calculates finite element data on the actual cells such as shape function values,...
@ FEScalar
Iter< Object > make_iter(Object obj)
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell.
unsigned int uint
ArmaMat< double, N, M > mat
Definition: armor.hh:936
@ opCoords
operations evaluated on quadrature points
@ opInvJac
inverse Jacobian
@ opJxW
JxW value of quadrature point.
@ opNormalVec
normal vector of quadrature point
std::vector< std::array< uint, 3 > > DimPointTable
Holds triplet (dim; bulk/side; idx of point in subtable)
Store finite element data on the actual patch such as shape function values, gradients,...
#define FMT_UNUSED
Definition: posix.h:75
Definitions of particular quadrature rules on simplices.
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc.
Definition: mixed.hh:25
Struct for pre-computing number of elements, sides, bulk points and side points on each dimension.
std::vector< std::vector< uint > > point_sizes_
void reset()
Set all values to zero.
std::vector< std::vector< uint > > elem_sizes_
static void ptop_scalar_shape(std::vector< ElOp< 3 >> &operations, std::vector< std::vector< double > > shape_values, uint scalar_shape_op_idx)
static void ptop_scalar_shape(std::vector< ElOp< 3 >> &operations, IntTableArena &el_table, std::vector< std::vector< std::vector< double > > > shape_values, uint scalar_shape_op_idx)
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.