Flow123d  DF_patch_fe_data_tables-6bd0dc8
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 begin, unsigned int op_idx)
60  : patch_point_vals_(patch_point_vals), begin_(begin), 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 begin_; /// Index of the first component of the bulk Quantity. Size is given by ValueType
69  unsigned int op_idx_; /// Index of operation in patch_point_vals_.operations vector
70 };
71 
72 
73 template <class ValueType>
74 class FeQ {
75 public:
76  /// Forbidden default constructor
77  FeQ() = delete;
78 
79  // Class similar to current FeView
80  FeQ(PatchPointValues<3> &patch_point_vals, unsigned int begin, unsigned int n_dofs)
81  : patch_point_vals_(patch_point_vals), begin_(begin), n_dofs_(n_dofs) {}
82 
83 
84  ValueType operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const BulkPoint &point);
85 
86  ValueType operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const SidePoint &point);
87 
88  // Implementation for EdgePoint, SidePoint, and JoinPoint shoud have a common implementation
89  // resolving to side values
90 
91 private:
92  PatchPointValues<3> &patch_point_vals_; ///< Reference to PatchPointValues
93  unsigned int begin_; ///< Index of the first component of the Quantity. Size is given by ValueType
94  unsigned int n_dofs_; ///< Number of DOFs
95 };
96 
97 
98 template <class ValueType>
100 public:
101  /// Default constructor
103  : patch_point_vals_bulk_(nullptr), patch_point_vals_side_(nullptr), join_idx_(-1) {}
104 
105  /**
106  * Constructor
107  *
108  * @param patch_point_vals_bulk Pointer to PatchPointValues bulk object.
109  * @param patch_point_vals_side Pointer to PatchPointValues side object.
110  * @param begin Index of the first component of the bulk Quantity.
111  * @param begin_side Index of the first component of the side Quantity.
112  * @param n_dofs_bulk Number of DOFs of bulk (lower-dim) element.
113  * @param n_dofs_side Number of DOFs of side (higher-dim) element.
114  * @param join_idx Index function.
115  */
116  JoinShapeAccessor(PatchPointValues<3> *patch_point_vals_bulk, PatchPointValues<3> *patch_point_vals_side,
117  unsigned int begin, unsigned int begin_side, unsigned int n_dofs_bulk, unsigned int n_dofs_side, unsigned int join_idx)
118  : patch_point_vals_bulk_(patch_point_vals_bulk), patch_point_vals_side_(patch_point_vals_side), begin_(begin),
119  begin_side_(begin_side), n_dofs_high_(n_dofs_side), n_dofs_low_(n_dofs_bulk), join_idx_(join_idx) {
120  //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!");
121  }
122 
123  /// Return global index of DOF
124  inline unsigned int join_idx() const {
125  return join_idx_;
126  }
127 
128  /// Return local index of DOF (on low / high-dim) - should be private method
129  inline unsigned int local_idx() const {
130  if (this->is_high_dim()) return (join_idx_ - n_dofs_low_);
131  else return join_idx_;
132  }
133 
134  inline unsigned int n_dofs_low() const {
135  return n_dofs_low_;
136  }
137 
138  inline unsigned int n_dofs_high() const {
139  return n_dofs_high_;
140  }
141 
142  inline unsigned int n_dofs_both() const {
143  return n_dofs_high_ + n_dofs_low_;
144  }
145 
146  inline bool is_high_dim() const {
147  return (join_idx_ >= n_dofs_low_);
148  }
149 
150  /// Iterates to next item.
151  inline void inc() {
152  join_idx_++;
153  }
154 
155  /// Comparison of accessors.
156  bool operator==(const JoinShapeAccessor<ValueType>& other) const {
157  return (join_idx_ == other.join_idx_);
158  }
159 
160 
161  ValueType operator()(const BulkPoint &point);
162 
163  ValueType operator()(const SidePoint &point);
164 
165 private:
166  // attributes:
167  PatchPointValues<3> *patch_point_vals_bulk_; ///< Pointer to bulk PatchPointValues
168  PatchPointValues<3> *patch_point_vals_side_; ///< Pointer to side PatchPointValues
169  unsigned int begin_; ///< Index of the first component of the bulk Quantity. Size is given by ValueType
170  unsigned int begin_side_; ///< Index of the first component of the side Quantity. Size is given by ValueType
171  unsigned int n_dofs_high_; ///< Number of DOFs on high-dim element
172  unsigned int n_dofs_low_; ///< Number of DOFs on low-dim element
173  unsigned int join_idx_; ///< Index of processed DOF
174 };
175 
176 
177 template<unsigned int dim>
179 {
180 protected:
181  // Default constructor
183  {}
184 
185  /// Return FiniteElement of \p component_idx for FESystem or \p fe for other types
186  template<unsigned int FE_dim>
187  std::shared_ptr<FiniteElement<FE_dim>> fe_comp(std::shared_ptr< FiniteElement<FE_dim> > fe, uint component_idx) {
188  FESystem<FE_dim> *fe_sys = dynamic_cast<FESystem<FE_dim>*>( fe.get() );
189  if (fe_sys != nullptr) {
190  return fe_sys->fe()[component_idx];
191  } else {
192  ASSERT_EQ(component_idx, 0).warning("Non-zero component_idx can only be used for FESystem.");
193  return fe;
194  }
195  }
196 
197  /**
198  * @brief Precomputed values of basis functions at the bulk quadrature points.
199  *
200  * Dimensions: (no. of quadrature points)
201  * x (no. of dofs)
202  * x (no. of components in ref. cell)
203  */
204  template<unsigned int FE_dim>
206  std::vector<std::vector<arma::vec> > ref_shape_vals( q->size(), vector<arma::vec>(fe->n_dofs()) );
207 
208  arma::mat shape_values(fe->n_dofs(), fe->n_components());
209  for (unsigned int i=0; i<q->size(); i++)
210  {
211  for (unsigned int j=0; j<fe->n_dofs(); j++)
212  {
213  for (unsigned int c=0; c<fe->n_components(); c++)
214  shape_values(j,c) = fe->shape_value(j, q->point<FE_dim>(i), c);
215 
216  ref_shape_vals[i][j] = trans(shape_values.row(j));
217  }
218  }
219 
220  return ref_shape_vals;
221  }
222 
223  /**
224  * @brief Precomputed values of basis functions at the side quadrature points.
225  *
226  * Dimensions: (sides)
227  * x (no. of quadrature points)
228  * x (no. of dofs)
229  * x (no. of components in ref. cell)
230  */
231  template<unsigned int FE_dim>
234 
235  arma::mat shape_values(fe->n_dofs(), fe->n_components());
236 
237  for (unsigned int sid=0; sid<FE_dim+1; sid++) {
238  auto quad = q->make_from_side<dim>(sid);
239  for (unsigned int i=0; i<quad.size(); i++)
240  {
241  for (unsigned int j=0; j<fe->n_dofs(); j++)
242  {
243  for (unsigned int c=0; c<fe->n_components(); c++) {
244  shape_values(j,c) = fe->shape_value(j, quad.template point<FE_dim>(i), c);
245  }
246 
247  ref_shape_vals[sid][i][j] = trans(shape_values.row(j));
248  }
249  }
250  }
251 
252  return ref_shape_vals;
253  }
254 
255 };
256 
257 template<unsigned int dim>
258 class BulkValues : public BaseValues<dim>
259 {
260 public:
261  /// Constructor
263  : BaseValues<dim>(), patch_point_vals_(patch_point_vals) {
264  ASSERT_EQ(patch_point_vals.dim(), dim);
265  fe_ = fe[Dim<dim>{}];
266  }
267 
268  /**
269  * @brief Register the product of Jacobian determinant and the quadrature
270  * weight at bulk quadrature points.
271  *
272  * @param quad Quadrature.
273  */
274  inline ElQ<Scalar> JxW()
275  {
278  }
279 
280  /// Create bulk accessor of coords entity
282  {
285  }
286 
287 // inline ElQ<Tensor> jacobian(std::initializer_list<Quadrature *> quad_list)
288 // {}
289 
290  /// Create bulk accessor of jac determinant entity
292  {
295  }
296 
297  /**
298  * @brief Return the value of the @p function_no-th shape function at
299  * the @p p bulk quadrature point.
300  *
301  * @param component_idx Number of the shape function.
302  */
303  inline FeQ<Scalar> scalar_shape(uint component_idx = 0)
304  {
305  auto fe_component = this->fe_comp(fe_, component_idx);
306  ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
307 
308  // use lambda reinit function
309  std::vector< std::vector<double> > shape_values( fe_component->n_dofs(), vector<double>(patch_point_vals_.get_quadrature()->size()) );
310  auto ref_shape_vals = this->ref_shape_values_bulk(patch_point_vals_.get_quadrature(), fe_component);
311  for (unsigned int i = 0; i < patch_point_vals_.get_quadrature()->size(); i++)
312  for (unsigned int j = 0; j < fe_component->n_dofs(); j++) {
313  shape_values[j][i] = ref_shape_vals[i][j][0];
314  }
315  uint scalar_shape_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
316  auto lambda_scalar_shape = [shape_values, scalar_shape_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
317  bulk_reinit::ptop_scalar_shape(operations, shape_values, scalar_shape_op_idx);
318  };
319  auto &scalar_shape_bulk_op = patch_point_vals_.make_fe_op({1}, lambda_scalar_shape, {}, fe_component->n_dofs());
320  // lambda_scalar_shape
321  uint begin = scalar_shape_bulk_op.result_row();
322 
323  return FeQ<Scalar>(patch_point_vals_, begin, fe_component->n_dofs());
324  }
325 
326 // inline FeQ<Vector> vector_shape(uint component_idx = 0)
327 // {}
328 
329 // inline FeQ<Tensor> tensor_shape(uint component_idx = 0)
330 // {}
331 
332  /**
333  * @brief Return the value of the @p function_no-th gradient shape function at
334  * the @p p bulk quadrature point.
335  *
336  * @param component_idx Number of the shape function.
337  */
338  inline FeQ<Vector> grad_scalar_shape(uint component_idx=0)
339  {
340  auto fe_component = this->fe_comp(fe_, component_idx);
341  ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of grad_scalar_shape accessor must be FEScalar!\n");
342 
343  // use lambda reinit function
344  auto ref_shape_grads = this->ref_shape_gradients(fe_component);
345  uint scalar_shape_grads_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
346  auto lambda_scalar_shape_grad = [ref_shape_grads, scalar_shape_grads_op_idx](std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
347  bulk_reinit::ptop_scalar_shape_grads<dim>(operations, op_results, ref_shape_grads, scalar_shape_grads_op_idx);
348  };
349  auto &grad_scalar_shape_bulk_op = patch_point_vals_.make_fe_op({3}, lambda_scalar_shape_grad, {FeBulk::BulkOps::opInvJac}, fe_component->n_dofs());
350  uint begin = grad_scalar_shape_bulk_op.result_row();
351 
352  return FeQ<Vector>(patch_point_vals_, begin, fe_component->n_dofs());
353  }
354 
355 // inline FeQ<Tensor> grad_vector_shape(std::initializer_list<Quadrature *> quad_list, unsigned int i_comp=0)
356 // {}
357 
358 private:
359  /**
360  * @brief Precomputed gradients of basis functions at the quadrature points.
361  *
362  * Dimensions: (no. of quadrature points)
363  * x (no. of dofs)
364  * x ((dim_ of. ref. cell)x(no. of components in ref. cell))
365  */
368  std::vector<std::vector<arma::mat> > ref_shape_grads( q->size(), vector<arma::mat>(fe->n_dofs()) );
369 
370  arma::mat grad(dim, fe->n_components());
371  for (unsigned int i_pt=0; i_pt<q->size(); i_pt++)
372  {
373  for (unsigned int i_dof=0; i_dof<fe->n_dofs(); i_dof++)
374  {
375  grad.zeros();
376  for (unsigned int c=0; c<fe->n_components(); c++)
377  grad.col(c) += fe->shape_grad(i_dof, q->point<dim>(i_pt), c);
378 
379  ref_shape_grads[i_pt][i_dof] = grad;
380  }
381  }
382 
383  return ref_shape_grads;
384  }
385 
387  std::shared_ptr< FiniteElement<dim> > fe_;
388 };
389 
390 
391 template<unsigned int dim>
392 class SideValues : public BaseValues<dim>
393 {
394 public:
395  /// Constructor
397  : BaseValues<dim>(), patch_point_vals_(patch_point_vals) {
398  ASSERT_EQ(patch_point_vals.dim(), dim);
399  fe_ = fe[Dim<dim>{}];
400  }
401 
402  /// Same as BulkValues::JxW but register at side quadrature points.
403  inline ElQ<Scalar> JxW()
404  {
407  }
408 
409  /**
410  * @brief Register the normal vector to a side at side quadrature points.
411  *
412  * @param quad Quadrature.
413  */
415  {
418  }
419 
420  /// Create side accessor of coords entity
422  {
425  }
426 
427  /// Create bulk accessor of jac determinant entity
429  {
432  }
433 
434  /// Same as BulkValues::scalar_shape but register at side quadrature points.
435  inline FeQ<Scalar> scalar_shape(uint component_idx = 0)
436  {
437  auto fe_component = this->fe_comp(fe_, component_idx);
438  ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
439 
440  // use lambda reinit function
442  dim+1,
444  );
445  auto ref_shape_vals = this->ref_shape_values_side(patch_point_vals_.get_quadrature(), fe_component);
446  for (unsigned int s=0; s<dim+1; ++s)
447  for (unsigned int i = 0; i < patch_point_vals_.get_quadrature()->size(); i++)
448  for (unsigned int j = 0; j < fe_component->n_dofs(); j++) {
449  shape_values[s][i][j] = ref_shape_vals[s][i][j][0];
450  }
451  uint scalar_shape_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
452  auto lambda_scalar_shape = [shape_values, scalar_shape_op_idx](std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
453  side_reinit::ptop_scalar_shape(operations, op_results, el_table, shape_values, scalar_shape_op_idx);
454  };
455  auto &scalar_shape_bulk_op = patch_point_vals_.make_fe_op({1}, lambda_scalar_shape, {}, fe_component->n_dofs());
456  uint begin = scalar_shape_bulk_op.result_row();
457 
458  return FeQ<Scalar>(patch_point_vals_, begin, fe_component->n_dofs());
459  }
460 
461  /// Same as BulkValues::grad_scalar_shape but register at side quadrature points.
462  inline FeQ<Vector> grad_scalar_shape(uint component_idx=0)
463  {
464  auto fe_component = this->fe_comp(fe_, component_idx);
465  ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of grad_scalar_shape accessor must be FEScalar!\n");
466 
467  // use lambda reinit function
468  auto ref_shape_grads = this->ref_shape_gradients(fe_component);
469  uint scalar_shape_grads_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
470  auto lambda_scalar_shape_grad = [ref_shape_grads, scalar_shape_grads_op_idx](std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
471  side_reinit::ptop_scalar_shape_grads<dim>(operations, op_results, el_table, ref_shape_grads, scalar_shape_grads_op_idx);
472  };
473  auto &grad_scalar_shape_side_op = patch_point_vals_.make_fe_op({3}, lambda_scalar_shape_grad, {FeSide::SideOps::opElInvJac}, fe_component->n_dofs());
474  uint begin = grad_scalar_shape_side_op.result_row();
475 
476  return FeQ<Vector>(patch_point_vals_, begin, fe_component->n_dofs());
477  }
478 
479 private:
480  /**
481  * @brief Precomputed gradients of basis functions at the quadrature points.
482  *
483  * Dimensions: (sides)
484  * x (no. of quadrature points)
485  * x (no. of dofs)
486  * x ((dim_ of. ref. cell)x(no. of components in ref. cell))
487  */
491 
492  arma::mat grad(dim, fe->n_components());
493  for (unsigned int sid=0; sid<dim+1; sid++) {
494  auto quad = q->make_from_side<dim>(sid);
495  for (unsigned int i_pt=0; i_pt<quad.size(); i_pt++)
496  {
497  for (unsigned int i_dof=0; i_dof<fe->n_dofs(); i_dof++)
498  {
499  grad.zeros();
500  for (unsigned int c=0; c<fe->n_components(); c++)
501  grad.col(c) += fe->shape_grad(i_dof, quad.template point<dim>(i_pt), c);
502 
503  ref_shape_grads[sid][i_pt][i_dof] = grad;
504  }
505  }
506  }
507 
508  return ref_shape_grads;
509  }
510 
512  std::shared_ptr< FiniteElement<dim> > fe_;
513 };
514 
515 
516 template<unsigned int dim>
517 class JoinValues : public BaseValues<dim>
518 {
519 public:
520  /// Constructor
521  JoinValues(PatchPointValues<3> *patch_point_vals_bulk, PatchPointValues<3> *patch_point_vals_side, MixedPtr<FiniteElement> fe)
522  : BaseValues<dim>(), patch_point_vals_bulk_(patch_point_vals_bulk), patch_point_vals_side_(patch_point_vals_side) {
523  ASSERT_EQ(patch_point_vals_bulk->dim(), dim-1);
524  ASSERT_EQ(patch_point_vals_side->dim(), dim);
525  fe_high_dim_ = fe[Dim<dim>{}];
526  fe_low_dim_ = fe[Dim<dim-1>{}];
527  }
528 
530  {
531  // element of lower dim (bulk points)
532  auto fe_component_low = this->fe_comp(fe_low_dim_, component_idx);
533  ASSERT_EQ(fe_component_low->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
534  // use lambda reinit function
535  std::vector< std::vector<double> > shape_values_bulk( patch_point_vals_bulk_->get_quadrature()->size(), vector<double>(fe_component_low->n_dofs()) );
536  auto ref_shape_vals_bulk = this->ref_shape_values_bulk(patch_point_vals_bulk_->get_quadrature(), fe_component_low);
537  for (unsigned int i = 0; i < patch_point_vals_bulk_->get_quadrature()->size(); i++)
538  for (unsigned int j = 0; j < fe_component_low->n_dofs(); j++) {
539  shape_values_bulk[i][j] = ref_shape_vals_bulk[i][j][0];
540  }
541  uint scalar_shape_op_idx_bulk = patch_point_vals_bulk_->operations_.size(); // index in operations_ vector
542  auto lambda_scalar_shape_bulk = [shape_values_bulk, scalar_shape_op_idx_bulk](std::vector<ElOp<3>> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
543  bulk_reinit::ptop_scalar_shape(operations, shape_values_bulk, scalar_shape_op_idx_bulk);
544  };
545  auto &grad_scalar_shape_bulk_op = patch_point_vals_bulk_->make_fe_op({1}, lambda_scalar_shape_bulk, {}, fe_component_low->n_dofs());
546  uint begin_bulk = grad_scalar_shape_bulk_op.result_row();
547 
548  // element of higher dim (side points)
549  auto fe_component_high = this->fe_comp(fe_high_dim_, component_idx);
550  ASSERT_EQ(fe_component_high->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
551  // use lambda reinit function
552  std::vector< std::vector< std::vector<double> > > shape_values_side(
553  dim+1,
555  );
556  auto ref_shape_vals_side = this->ref_shape_values_side(patch_point_vals_side_->get_quadrature(), fe_component_high);
557  for (unsigned int s=0; s<dim+1; ++s)
558  for (unsigned int i = 0; i < patch_point_vals_side_->get_quadrature()->size(); i++)
559  for (unsigned int j = 0; j < fe_component_high->n_dofs(); j++) {
560  shape_values_side[s][i][j] = ref_shape_vals_side[s][i][j][0];
561  }
562  uint scalar_shape_op_idx_side = patch_point_vals_side_->operations_.size(); // index in operations_ vector
563  auto lambda_scalar_shape_side = [shape_values_side, scalar_shape_op_idx_side](std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
564  side_reinit::ptop_scalar_shape(operations, op_results, el_table, shape_values_side, scalar_shape_op_idx_side);
565  };
566  auto &grad_scalar_shape_side_op = patch_point_vals_side_->make_fe_op({1}, lambda_scalar_shape_side, {}, fe_component_high->n_dofs());
567  uint begin_side = grad_scalar_shape_side_op.result_row();
568 
569  auto bgn_it = make_iter<JoinShapeAccessor<Scalar>>( JoinShapeAccessor<Scalar>(patch_point_vals_bulk_, patch_point_vals_side_,
570  begin_bulk, begin_side, fe_component_low->n_dofs(), fe_component_high->n_dofs(), 0) );
571  unsigned int end_idx = fe_component_low->n_dofs() + fe_component_high->n_dofs();
572  auto end_it = make_iter<JoinShapeAccessor<Scalar>>( JoinShapeAccessor<Scalar>(patch_point_vals_bulk_, patch_point_vals_side_,
573  begin_bulk, begin_side, fe_component_low->n_dofs(), fe_component_high->n_dofs(), end_idx) );
574  return Range<JoinShapeAccessor<Scalar>>(bgn_it, end_it);
575  }
576 
577 private:
580  std::shared_ptr< FiniteElement<dim> > fe_high_dim_;
581  std::shared_ptr< FiniteElement<dim-1> > fe_low_dim_;
582 };
583 
584 /// Template specialization of dim = 1
585 template <>
586 class JoinValues<1> : public BaseValues<1>
587 {
588 public:
589  /// Constructor
591  : BaseValues<1>() {}
592 
594  {
596  make_iter<JoinShapeAccessor<Scalar>>(JoinShapeAccessor<Scalar>()),
598  }
599 };
600 
601 
602 template<unsigned int spacedim = 3>
604 public:
605  /// Struct for pre-computing number of elements, sides, bulk points and side points on each dimension.
606  struct TableSizes {
607  public:
608  /// Constructor
612  }
613 
614  /// Set all values to zero
615  void reset() {
616  std::fill(elem_sizes_[0].begin(), elem_sizes_[0].end(), 0);
617  std::fill(elem_sizes_[1].begin(), elem_sizes_[1].end(), 0);
618  std::fill(point_sizes_[0].begin(), point_sizes_[0].end(), 0);
619  std::fill(point_sizes_[1].begin(), point_sizes_[1].end(), 0);
620  }
621 
622  /**
623  * Holds number of elements and sides on each dimension
624  * Format:
625  * { {n_elements_1D, n_elements_2D, n_elements_3D },
626  * {n_sides_1D, n_sides_2D, n_sides_3D } }
627  */
629 
630  /**
631  * Holds number of bulk and side points on each dimension
632  * Format:
633  * { {n_bulk_points_1D, n_bulk_points_2D, n_bulk_points_3D },
634  * {n_side_points_1D, n_side_points_2D, n_side_points_3D } }
635  */
637  };
638 
640  : asm_arena_(1024 * 1024, 256),
641  patch_arena_(nullptr),
645  patch_point_vals_side_{ {FeSide::PatchPointValues(1, 0, asm_arena_),
646  FeSide::PatchPointValues(2, 0, asm_arena_),
647  FeSide::PatchPointValues(3, 0, asm_arena_)} }
648  {
649  used_quads_[0] = false; used_quads_[1] = false;
650  }
651 
652  PatchFEValues(unsigned int quad_order, MixedPtr<FiniteElement> fe)
653  : asm_arena_(1024 * 1024, 256),
654  patch_arena_(nullptr),
655  patch_point_vals_bulk_{ {FeBulk::PatchPointValues(1, quad_order, asm_arena_),
656  FeBulk::PatchPointValues(2, quad_order, asm_arena_),
657  FeBulk::PatchPointValues(3, quad_order, asm_arena_)} },
658  patch_point_vals_side_{ {FeSide::PatchPointValues(1, quad_order, asm_arena_),
659  FeSide::PatchPointValues(2, quad_order, asm_arena_),
660  FeSide::PatchPointValues(3, quad_order, asm_arena_)} },
661  fe_(fe)
662  {
663  used_quads_[0] = false; used_quads_[1] = false;
664  }
665 
666 
667  /// Destructor
669  {
670  if (patch_arena_!=nullptr)
671  delete patch_arena_;
672  }
673 
674  /// Return bulk or side quadrature of given dimension
675  Quadrature *get_quadrature(uint dim, bool is_bulk) const {
676  if (is_bulk) return patch_point_vals_bulk_[dim-1].get_quadrature();
677  else return patch_point_vals_side_[dim-1].get_quadrature();
678  }
679 
680  /**
681  * @brief Initialize structures and calculates cell-independent data.
682  *
683  * @param _quadrature The quadrature rule for the cell associated
684  * to given finite element or for the cell side.
685  * @param _flags The update flags.
686  */
687  template<unsigned int DIM>
688  void initialize(Quadrature &_quadrature)
689  {
690  if ( _quadrature.dim() == DIM ) {
691  used_quads_[0] = true;
692  patch_point_vals_bulk_[DIM-1].initialize(3); // bulk
693  } else {
694  used_quads_[1] = true;
695  patch_point_vals_side_[DIM-1].initialize(4); // side
696  }
697  }
698 
699  /// Finalize initialization, creates child (patch) arena and passes it to PatchPointValue objects
700  void init_finalize() {
701  patch_arena_ = asm_arena_.get_child_arena();
702  for (unsigned int i=0; i<3; ++i) {
703  if (used_quads_[0]) patch_point_vals_bulk_[i].init_finalize(patch_arena_);
704  if (used_quads_[1]) patch_point_vals_side_[i].init_finalize(patch_arena_);
705  }
706  }
707 
708  /// Reset PatchpointValues structures
709  void reset()
710  {
711  for (unsigned int i=0; i<3; ++i) {
712  if (used_quads_[0]) patch_point_vals_bulk_[i].reset();
713  if (used_quads_[1]) patch_point_vals_side_[i].reset();
714  }
715  patch_arena_->reset();
716  }
717 
718  /// Reinit data.
720  {
721  for (unsigned int i=0; i<3; ++i) {
722  if (used_quads_[0]) patch_point_vals_bulk_[i].reinit_patch();
723  if (used_quads_[1]) patch_point_vals_side_[i].reinit_patch();
724  }
725  }
726 
727  /**
728  * @brief Returns the number of shape functions.
729  */
730  template<unsigned int dim>
731  inline unsigned int n_dofs() const {
732  ASSERT((dim>=0) && (dim<=3))(dim).error("Dimension must be 0, 1, 2 or 3.");
733  return fe_[Dim<dim>{}]->n_dofs();
734  }
735 
736  /// Return BulkValue object of dimension given by template parameter
737  template<unsigned int dim>
739  ASSERT((dim>0) && (dim<=3))(dim).error("Dimension must be 1, 2 or 3.");
740  return BulkValues<dim>(patch_point_vals_bulk_[dim-1], fe_);
741  }
742 
743  /// Return SideValue object of dimension given by template parameter
744  template<unsigned int dim>
746  ASSERT((dim>0) && (dim<=3))(dim).error("Dimension must be 1, 2 or 3.");
747  return SideValues<dim>(patch_point_vals_side_[dim-1], fe_);
748  }
749 
750  /// Return JoinValue object of dimension given by template parameter
751  template<unsigned int dim>
753  //ASSERT((dim>1) && (dim<=3))(dim).error("Dimension must be 2 or 3.");
754  return JoinValues<dim>(&patch_point_vals_bulk_[dim-2], &patch_point_vals_side_[dim-1], fe_);
755  }
756 
757  /** Following methods are used during update of patch. **/
758 
759  /// Resize tables of patch_point_vals_
760  void resize_tables(TableSizes table_sizes) {
761  for (uint i=0; i<3; ++i) {
762  if (used_quads_[0]) patch_point_vals_bulk_[i].resize_tables(table_sizes.elem_sizes_[0][i], table_sizes.point_sizes_[0][i]);
763  if (used_quads_[1]) patch_point_vals_side_[i].resize_tables(table_sizes.elem_sizes_[1][i], table_sizes.point_sizes_[1][i]);
764  }
765  }
766 
767  /// Register element to patch_point_vals_ table by dimension of element
768  uint register_element(DHCellAccessor cell, uint element_patch_idx) {
769  arma::mat coords;
770  switch (cell.dim()) {
771  case 1:
772  coords = MappingP1<1,spacedim>::element_map(cell.elm());
773  return patch_point_vals_bulk_[0].register_element(coords, element_patch_idx);
774  break;
775  case 2:
776  coords = MappingP1<2,spacedim>::element_map(cell.elm());
777  return patch_point_vals_bulk_[1].register_element(coords, element_patch_idx);
778  break;
779  case 3:
780  coords = MappingP1<3,spacedim>::element_map(cell.elm());
781  return patch_point_vals_bulk_[2].register_element(coords, element_patch_idx);
782  break;
783  default:
784  ASSERT(false);
785  return 0;
786  break;
787  }
788  }
789 
790  /// Register side to patch_point_vals_ table by dimension of side
792  arma::mat side_coords(spacedim, cell_side.dim());
793  for (unsigned int n=0; n<cell_side.dim(); n++)
794  for (unsigned int c=0; c<spacedim; c++)
795  side_coords(c,n) = (*cell_side.side().node(n))[c];
796 
797  arma::mat elm_coords;
798  DHCellAccessor cell = cell_side.cell();
799  switch (cell.dim()) {
800  case 1:
801  elm_coords = MappingP1<1,spacedim>::element_map(cell.elm());
802  return patch_point_vals_side_[0].register_side(elm_coords, side_coords);
803  break;
804  case 2:
805  elm_coords = MappingP1<2,spacedim>::element_map(cell.elm());
806  return patch_point_vals_side_[1].register_side(elm_coords, side_coords);
807  break;
808  case 3:
809  elm_coords = MappingP1<3,spacedim>::element_map(cell.elm());
810  return patch_point_vals_side_[2].register_side(elm_coords, side_coords);
811  break;
812  default:
813  ASSERT(false);
814  return 0;
815  break;
816  }
817  }
818 
819  /// Register bulk point to patch_point_vals_ table by dimension of element
820  uint register_bulk_point(DHCellAccessor cell, uint elem_table_row, uint value_patch_idx, uint i_point_on_elem) {
821  return patch_point_vals_bulk_[cell.dim()-1].register_bulk_point(elem_table_row, value_patch_idx, cell.elm_idx(), i_point_on_elem);
822  }
823 
824  /// Register side point to patch_point_vals_ table by dimension of side
825  uint register_side_point(DHCellSide cell_side, uint elem_table_row, uint value_patch_idx) {
826  return patch_point_vals_side_[cell_side.dim()-1].register_side_point(elem_table_row, value_patch_idx, cell_side.elem_idx(), cell_side.side_idx());
827  }
828 
829  /// Temporary development method
830  void print_data_tables(ostream& stream, bool points, bool ints, bool only_bulk=true) const {
831  stream << endl << "Table of patch FE data:" << endl;
832  for (uint i=0; i<3; ++i) {
833  stream << std::setfill('-') << setw(100) << "" << endl;
834  stream << "Bulk, dimension " << (i+1) << endl;
835  patch_point_vals_bulk_[i].print_data_tables(stream, points, ints);
836  }
837  if (!only_bulk)
838  for (uint i=0; i<3; ++i) {
839  stream << std::setfill('-') << setw(100) << "" << endl;
840  stream << "Side, dimension " << (i+1) << endl;
841  patch_point_vals_side_[i].print_data_tables(stream, points, ints);
842  }
843  stream << std::setfill('=') << setw(100) << "" << endl;
844  }
845 
846  /// Temporary development method
847  void print_operations(ostream& stream) const {
848  stream << endl << "Table of patch FE operations:" << endl;
849  for (uint i=0; i<3; ++i) {
850  stream << std::setfill('-') << setw(100) << "" << endl;
851  stream << "Bulk, dimension " << (i+1) << ", n_rows " << patch_point_vals_bulk_[i].n_rows() << endl;
852  patch_point_vals_bulk_[i].print_operations(stream, 0);
853  }
854  for (uint i=0; i<3; ++i) {
855  stream << std::setfill('-') << setw(100) << "" << endl;
856  stream << "Side, dimension " << (i+1) << ", n_rows " << patch_point_vals_side_[i].n_rows() << endl;
857  patch_point_vals_side_[i].print_operations(stream, 1);
858  }
859  stream << std::setfill('=') << setw(100) << "" << endl;
860  }
861 
862 private:
865  std::array<FeBulk::PatchPointValues<spacedim>, 3> patch_point_vals_bulk_; ///< Sub objects of bulk data of dimensions 1,2,3
866  std::array<FeSide::PatchPointValues<spacedim>, 3> patch_point_vals_side_; ///< Sub objects of side data of dimensions 1,2,3
867 
868  MixedPtr<FiniteElement> fe_; ///< Mixed of shared pointers of FiniteElement object
869  bool used_quads_[2]; ///< Pair of flags signs holds info if bulk and side quadratures are used
870 
871  template <class ValueType>
872  friend class ElQ;
873  template <class ValueType>
874  friend class FeQ;
875  template <class ValueType>
876  friend class JoinShapeAccessor;
877 };
878 
879 
880 template <class ValueType>
881 ValueType ElQ<ValueType>::operator()(const BulkPoint &point) {
882  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
883  return patch_point_vals_.scalar_value(op_idx_, value_cache_idx);
884 }
885 
886 template <>
888  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
889  return patch_point_vals_.vector_value(op_idx_, value_cache_idx);
890 }
891 
892 template <>
894  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
895  return patch_point_vals_.tensor_value(op_idx_, value_cache_idx);
896 }
897 
898 template <class ValueType>
899 ValueType ElQ<ValueType>::operator()(const SidePoint &point) {
900  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
901  return patch_point_vals_.scalar_val(begin_, value_cache_idx);
902 }
903 
904 template <>
906  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
907  return patch_point_vals_.vector_val(begin_, value_cache_idx);
908 }
909 
910 template <>
912  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
913  return patch_point_vals_.tensor_val(begin_, value_cache_idx);
914 }
915 
916 template <class ValueType>
917 ValueType FeQ<ValueType>::operator()(unsigned int shape_idx, const BulkPoint &point) {
918  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
919  return patch_point_vals_.scalar_val(begin_+shape_idx, value_cache_idx);
920 }
921 
922 template <>
923 inline Vector FeQ<Vector>::operator()(unsigned int shape_idx, const BulkPoint &point) {
924  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
925  return patch_point_vals_.vector_val(begin_+3*shape_idx, value_cache_idx);
926 }
927 
928 template <>
929 inline Tensor FeQ<Tensor>::operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const BulkPoint &point) {
930  Tensor tens; tens.zeros();
931  return tens;
932 }
933 
934 template <class ValueType>
935 ValueType FeQ<ValueType>::operator()(unsigned int shape_idx, const SidePoint &point) {
936  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
937  return patch_point_vals_.scalar_val(begin_+shape_idx, value_cache_idx);
938 }
939 
940 template <>
941 inline Vector FeQ<Vector>::operator()(unsigned int shape_idx, const SidePoint &point) {
942  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
943  return patch_point_vals_.vector_val(begin_+3*shape_idx, value_cache_idx);
944 }
945 
946 template <>
947 inline Tensor FeQ<Tensor>::operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const SidePoint &point) {
948  Tensor tens; tens.zeros();
949  return tens;
950 }
951 
952 
953 template <class ValueType>
955  if (this->is_high_dim()) {
956  return 0.0;
957  } else {
958  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
959  return patch_point_vals_bulk_->scalar_val(begin_+this->local_idx(), value_cache_idx);
960  }
961 }
962 
963 template <>
965  Vector vect; vect.zeros();
966  return vect;
967 }
968 
969 template <>
971  Tensor tens; tens.zeros();
972  return tens;
973 }
974 
975 template <class ValueType>
977  if (this->is_high_dim()) {
978  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
979  return patch_point_vals_side_->scalar_val(begin_side_+this->local_idx(), value_cache_idx);
980  } else {
981  return 0.0;
982  }
983 }
984 
985 template <>
987  Vector vect; vect.zeros();
988  return vect;
989 }
990 
991 template <>
993  Tensor tens; tens.zeros();
994  return tens;
995 }
996 
997 
998 #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 begin, unsigned int op_idx)
Constructor.
unsigned int op_idx_
Index of the first component of the bulk Quantity. Size is given by ValueType.
ElQ()=delete
Forbidden default constructor.
PatchPointValues< 3 > & patch_point_vals_
Reference to PatchPointValues.
unsigned int begin_
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.
unsigned int begin_
Index of the first component of the Quantity. Size is given by ValueType.
PatchPointValues< 3 > & patch_point_vals_
Reference to PatchPointValues.
unsigned int n_dofs_
Number of DOFs.
ValueType operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const BulkPoint &point)
FeQ(PatchPointValues< 3 > &patch_point_vals, unsigned int begin, unsigned int n_dofs)
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 begin_side_
Index of the first component of the side Quantity. Size is given by ValueType.
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.
bool is_high_dim() const
unsigned int begin_
Index of the first component of the bulk Quantity. Size is given by ValueType.
PatchPointValues< 3 > * patch_point_vals_bulk_
Pointer to bulk PatchPointValues.
JoinShapeAccessor(PatchPointValues< 3 > *patch_point_vals_bulk, PatchPointValues< 3 > *patch_point_vals_side, unsigned int begin, unsigned int begin_side, unsigned int n_dofs_bulk, unsigned int n_dofs_side, unsigned int join_idx)
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.
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.
AssemblyArena * patch_arena_
BulkValues< dim > bulk_values()
Return BulkValue object of dimension given by template parameter.
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.
uint register_side_point(DHCellSide cell_side, uint elem_table_row, uint value_patch_idx)
Register side point to patch_point_vals_ table by dimension of side.
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_val(uint result_row, uint point_idx) 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< ArrayDbl, Eigen::Dynamic > TableDbl
Definition: eigen_tools.hh:36
Eigen::Vector< ArrayInt, Eigen::Dynamic > TableInt
Definition: eigen_tools.hh:37
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, TableDbl &op_results, TableInt &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.