Flow123d  DF_patch_fe_data_tables-522f93f
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 
43 template<unsigned int spacedim> class PatchFEValues;
44 
45 
46 
47 typedef typename std::vector< std::array<uint, 3> > DimPointTable; ///< Holds triplet (dim; bulk/side; idx of point in subtable)
48 
49 
50 template <class ValueType>
51 class ElQ {
52 public:
53  /// Forbidden default constructor
54  ElQ() = delete;
55 
56  /// Constructor
57  ElQ(PatchPointValues<3> &patch_point_vals, unsigned int begin)
58  : patch_point_vals_(patch_point_vals), begin_(begin) {}
59 
60  ValueType operator()(FMT_UNUSED const BulkPoint &point);
61 
62  ValueType operator()(FMT_UNUSED const SidePoint &point);
63 
64 private:
65  PatchPointValues<3> &patch_point_vals_; ///< Reference to PatchPointValues
66  unsigned int begin_; /// Index of the first component of the bulk Quantity. Size is given by ValueType
67 };
68 
69 
70 template <class ValueType>
71 class FeQ {
72 public:
73  /// Forbidden default constructor
74  FeQ() = delete;
75 
76  // Class similar to current FeView
77  FeQ(PatchPointValues<3> &patch_point_vals, unsigned int begin, unsigned int n_dofs)
78  : patch_point_vals_(patch_point_vals), begin_(begin), n_dofs_(n_dofs) {}
79 
80 
81  ValueType operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const BulkPoint &point);
82 
83  ValueType operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const SidePoint &point);
84 
85  // Implementation for EdgePoint, SidePoint, and JoinPoint shoud have a common implementation
86  // resolving to side values
87 
88 private:
89  PatchPointValues<3> &patch_point_vals_; ///< Reference to PatchPointValues
90  unsigned int begin_; ///< Index of the first component of the Quantity. Size is given by ValueType
91  unsigned int n_dofs_; ///< Number of DOFs
92 };
93 
94 
95 template <class ValueType>
97 public:
98  /// Default constructor
100  : patch_point_vals_bulk_(nullptr), patch_point_vals_side_(nullptr), join_idx_(-1) {}
101 
102  /**
103  * Constructor
104  *
105  * @param patch_point_vals_bulk Pointer to PatchPointValues bulk object.
106  * @param patch_point_vals_side Pointer to PatchPointValues side object.
107  * @param begin Index of the first component of the bulk Quantity.
108  * @param begin_side Index of the first component of the side Quantity.
109  * @param n_dofs_bulk Number of DOFs of bulk (lower-dim) element.
110  * @param n_dofs_side Number of DOFs of side (higher-dim) element.
111  * @param join_idx Index function.
112  */
113  JoinShapeAccessor(PatchPointValues<3> *patch_point_vals_bulk, PatchPointValues<3> *patch_point_vals_side,
114  unsigned int begin, unsigned int begin_side, unsigned int n_dofs_bulk, unsigned int n_dofs_side, unsigned int join_idx)
115  : patch_point_vals_bulk_(patch_point_vals_bulk), patch_point_vals_side_(patch_point_vals_side), begin_(begin),
116  begin_side_(begin_side), n_dofs_high_(n_dofs_side), n_dofs_low_(n_dofs_bulk), join_idx_(join_idx) {
117  //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!");
118  }
119 
120  /// Return global index of DOF
121  inline unsigned int join_idx() const {
122  return join_idx_;
123  }
124 
125  /// Return local index of DOF (on low / high-dim) - should be private method
126  inline unsigned int local_idx() const {
127  if (this->is_high_dim()) return (join_idx_ - n_dofs_low_);
128  else return join_idx_;
129  }
130 
131  inline unsigned int n_dofs_low() const {
132  return n_dofs_low_;
133  }
134 
135  inline unsigned int n_dofs_high() const {
136  return n_dofs_high_;
137  }
138 
139  inline unsigned int n_dofs_both() const {
140  return n_dofs_high_ + n_dofs_low_;
141  }
142 
143  inline bool is_high_dim() const {
144  return (join_idx_ >= n_dofs_low_);
145  }
146 
147  /// Iterates to next item.
148  inline void inc() {
149  join_idx_++;
150  }
151 
152  /// Comparison of accessors.
153  bool operator==(const JoinShapeAccessor<ValueType>& other) const {
154  return (join_idx_ == other.join_idx_);
155  }
156 
157 
158  ValueType operator()(const BulkPoint &point);
159 
160  ValueType operator()(const SidePoint &point);
161 
162 private:
163  // attributes:
164  PatchPointValues<3> *patch_point_vals_bulk_; ///< Pointer to bulk PatchPointValues
165  PatchPointValues<3> *patch_point_vals_side_; ///< Pointer to side PatchPointValues
166  unsigned int begin_; ///< Index of the first component of the bulk Quantity. Size is given by ValueType
167  unsigned int begin_side_; ///< Index of the first component of the side Quantity. Size is given by ValueType
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 join_idx_; ///< Index of processed DOF
171 };
172 
173 
174 template<unsigned int dim>
176 {
177 protected:
178  // Default constructor
180  {}
181 
182  /// Return FiniteElement of \p component_idx for FESystem or \p fe for other types
183  template<unsigned int FE_dim>
184  std::shared_ptr<FiniteElement<FE_dim>> fe_comp(std::shared_ptr< FiniteElement<FE_dim> > fe, uint component_idx) {
185  FESystem<FE_dim> *fe_sys = dynamic_cast<FESystem<FE_dim>*>( fe.get() );
186  if (fe_sys != nullptr) {
187  return fe_sys->fe()[component_idx];
188  } else {
189  ASSERT_EQ(component_idx, 0).warning("Non-zero component_idx can only be used for FESystem.");
190  return fe;
191  }
192  }
193 
194  /**
195  * @brief Precomputed values of basis functions at the bulk quadrature points.
196  *
197  * Dimensions: (no. of quadrature points)
198  * x (no. of dofs)
199  * x (no. of components in ref. cell)
200  */
201  template<unsigned int FE_dim>
203  std::vector<std::vector<arma::vec> > ref_shape_vals( q->size(), vector<arma::vec>(fe->n_dofs()) );
204 
205  arma::mat shape_values(fe->n_dofs(), fe->n_components());
206  for (unsigned int i=0; i<q->size(); i++)
207  {
208  for (unsigned int j=0; j<fe->n_dofs(); j++)
209  {
210  for (unsigned int c=0; c<fe->n_components(); c++)
211  shape_values(j,c) = fe->shape_value(j, q->point<FE_dim>(i), c);
212 
213  ref_shape_vals[i][j] = trans(shape_values.row(j));
214  }
215  }
216 
217  return ref_shape_vals;
218  }
219 
220  /**
221  * @brief Precomputed values of basis functions at the side quadrature points.
222  *
223  * Dimensions: (sides)
224  * x (no. of quadrature points)
225  * x (no. of dofs)
226  * x (no. of components in ref. cell)
227  */
228  template<unsigned int FE_dim>
231 
232  arma::mat shape_values(fe->n_dofs(), fe->n_components());
233 
234  for (unsigned int sid=0; sid<FE_dim+1; sid++) {
235  auto quad = q->make_from_side<dim>(sid);
236  for (unsigned int i=0; i<quad.size(); i++)
237  {
238  for (unsigned int j=0; j<fe->n_dofs(); j++)
239  {
240  for (unsigned int c=0; c<fe->n_components(); c++) {
241  shape_values(j,c) = fe->shape_value(j, quad.template point<FE_dim>(i), c);
242  }
243 
244  ref_shape_vals[sid][i][j] = trans(shape_values.row(j));
245  }
246  }
247  }
248 
249  return ref_shape_vals;
250  }
251 
252 };
253 
254 template<unsigned int dim>
255 class BulkValues : public BaseValues<dim>
256 {
257 public:
258  /// Constructor
260  : BaseValues<dim>(), patch_point_vals_(patch_point_vals) {
261  ASSERT_EQ(patch_point_vals.dim(), dim);
262  fe_ = fe[Dim<dim>{}];
263  }
264 
265  /**
266  * @brief Register the product of Jacobian determinant and the quadrature
267  * weight at bulk quadrature points.
268  *
269  * @param quad Quadrature.
270  */
271  inline ElQ<Scalar> JxW()
272  {
274  return ElQ<Scalar>(patch_point_vals_, begin);
275  }
276 
277  /// Create bulk accessor of coords entity
279  {
281  return ElQ<Vector>(patch_point_vals_, begin);
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  return ElQ<Scalar>(patch_point_vals_, begin);
292  }
293 
294  /**
295  * @brief Return the value of the @p function_no-th shape function at
296  * the @p p bulk quadrature point.
297  *
298  * @param component_idx Number of the shape function.
299  */
300  inline FeQ<Scalar> scalar_shape(uint component_idx = 0)
301  {
302  auto fe_component = this->fe_comp(fe_, component_idx);
303  ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
304 
305  // use lambda reinit function
306  std::vector< std::vector<double> > shape_values( patch_point_vals_.get_quadrature()->size(), vector<double>(fe_component->n_dofs()) );
307  auto ref_shape_vals = this->ref_shape_values_bulk(patch_point_vals_.get_quadrature(), fe_component);
308  for (unsigned int i = 0; i < patch_point_vals_.get_quadrature()->size(); i++)
309  for (unsigned int j = 0; j < fe_component->n_dofs(); j++) {
310  shape_values[i][j] = ref_shape_vals[i][j][0];
311  }
312  uint scalar_shape_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
313  auto lambda_scalar_shape = [shape_values, scalar_shape_op_idx](std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
314  bulk_reinit::ptop_scalar_shape(operations, op_results, shape_values, scalar_shape_op_idx);
315  };
316  auto &scalar_shape_bulk_op = patch_point_vals_.make_fe_op({1}, lambda_scalar_shape, {}, fe_component->n_dofs());
317  uint begin = scalar_shape_bulk_op.result_row();
318 
319  return FeQ<Scalar>(patch_point_vals_, begin, fe_component->n_dofs());
320  }
321 
322 // inline FeQ<Vector> vector_shape(uint component_idx = 0)
323 // {}
324 
325 // inline FeQ<Tensor> tensor_shape(uint component_idx = 0)
326 // {}
327 
328  /**
329  * @brief Return the value of the @p function_no-th gradient shape function at
330  * the @p p bulk quadrature point.
331  *
332  * @param component_idx Number of the shape function.
333  */
334  inline FeQ<Vector> grad_scalar_shape(uint component_idx=0)
335  {
336  auto fe_component = this->fe_comp(fe_, component_idx);
337  ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of grad_scalar_shape accessor must be FEScalar!\n");
338 
339  // use lambda reinit function
340  auto ref_shape_grads = this->ref_shape_gradients(fe_component);
341  uint scalar_shape_grads_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
342  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) {
343  bulk_reinit::ptop_scalar_shape_grads<dim>(operations, op_results, ref_shape_grads, scalar_shape_grads_op_idx);
344  };
345  auto &grad_scalar_shape_bulk_op = patch_point_vals_.make_fe_op({3}, lambda_scalar_shape_grad, {FeBulk::BulkOps::opInvJac}, fe_component->n_dofs());
346  uint begin = grad_scalar_shape_bulk_op.result_row();
347 
348  return FeQ<Vector>(patch_point_vals_, begin, fe_component->n_dofs());
349  }
350 
351 // inline FeQ<Tensor> grad_vector_shape(std::initializer_list<Quadrature *> quad_list, unsigned int i_comp=0)
352 // {}
353 
354 private:
355  /**
356  * @brief Precomputed gradients of basis functions at the quadrature points.
357  *
358  * Dimensions: (no. of quadrature points)
359  * x (no. of dofs)
360  * x ((dim_ of. ref. cell)x(no. of components in ref. cell))
361  */
364  std::vector<std::vector<arma::mat> > ref_shape_grads( q->size(), vector<arma::mat>(fe->n_dofs()) );
365 
366  arma::mat grad(dim, fe->n_components());
367  for (unsigned int i_pt=0; i_pt<q->size(); i_pt++)
368  {
369  for (unsigned int i_dof=0; i_dof<fe->n_dofs(); i_dof++)
370  {
371  grad.zeros();
372  for (unsigned int c=0; c<fe->n_components(); c++)
373  grad.col(c) += fe->shape_grad(i_dof, q->point<dim>(i_pt), c);
374 
375  ref_shape_grads[i_pt][i_dof] = grad;
376  }
377  }
378 
379  return ref_shape_grads;
380  }
381 
383  std::shared_ptr< FiniteElement<dim> > fe_;
384 };
385 
386 
387 template<unsigned int dim>
388 class SideValues : public BaseValues<dim>
389 {
390 public:
391  /// Constructor
393  : BaseValues<dim>(), patch_point_vals_(patch_point_vals) {
394  ASSERT_EQ(patch_point_vals.dim(), dim);
395  fe_ = fe[Dim<dim>{}];
396  }
397 
398  /// Same as BulkValues::JxW but register at side quadrature points.
399  inline ElQ<Scalar> JxW()
400  {
402  return ElQ<Scalar>(patch_point_vals_, begin);
403  }
404 
405  /**
406  * @brief Register the normal vector to a side at side quadrature points.
407  *
408  * @param quad Quadrature.
409  */
411  {
413  return ElQ<Vector>(patch_point_vals_, begin);
414  }
415 
416  /// Create side accessor of coords entity
418  {
420  return ElQ<Vector>(patch_point_vals_, begin);
421  }
422 
423  /// Create bulk accessor of jac determinant entity
425  {
427  return ElQ<Scalar>(patch_point_vals_, begin);
428  }
429 
430  /// Same as BulkValues::scalar_shape but register at side quadrature points.
431  inline FeQ<Scalar> scalar_shape(uint component_idx = 0)
432  {
433  auto fe_component = this->fe_comp(fe_, component_idx);
434  ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
435 
436  // use lambda reinit function
438  dim+1,
440  );
441  auto ref_shape_vals = this->ref_shape_values_side(patch_point_vals_.get_quadrature(), fe_component);
442  for (unsigned int s=0; s<dim+1; ++s)
443  for (unsigned int i = 0; i < patch_point_vals_.get_quadrature()->size(); i++)
444  for (unsigned int j = 0; j < fe_component->n_dofs(); j++) {
445  shape_values[s][i][j] = ref_shape_vals[s][i][j][0];
446  }
447  uint scalar_shape_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
448  auto lambda_scalar_shape = [shape_values, scalar_shape_op_idx](std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
449  side_reinit::ptop_scalar_shape(operations, op_results, el_table, shape_values, scalar_shape_op_idx);
450  };
451  auto &scalar_shape_bulk_op = patch_point_vals_.make_fe_op({1}, lambda_scalar_shape, {}, fe_component->n_dofs());
452  uint begin = scalar_shape_bulk_op.result_row();
453 
454  return FeQ<Scalar>(patch_point_vals_, begin, fe_component->n_dofs());
455  }
456 
457  /// Same as BulkValues::grad_scalar_shape but register at side quadrature points.
458  inline FeQ<Vector> grad_scalar_shape(uint component_idx=0)
459  {
460  auto fe_component = this->fe_comp(fe_, component_idx);
461  ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of grad_scalar_shape accessor must be FEScalar!\n");
462 
463  // use lambda reinit function
464  auto ref_shape_grads = this->ref_shape_gradients(fe_component);
465  uint scalar_shape_grads_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
466  auto lambda_scalar_shape_grad = [ref_shape_grads, scalar_shape_grads_op_idx](std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
467  side_reinit::ptop_scalar_shape_grads<dim>(operations, op_results, el_table, ref_shape_grads, scalar_shape_grads_op_idx);
468  };
469  auto &grad_scalar_shape_side_op = patch_point_vals_.make_fe_op({3}, lambda_scalar_shape_grad, {FeSide::SideOps::opElInvJac}, fe_component->n_dofs());
470  uint begin = grad_scalar_shape_side_op.result_row();
471 
472  return FeQ<Vector>(patch_point_vals_, begin, fe_component->n_dofs());
473  }
474 
475 private:
476  /**
477  * @brief Precomputed gradients of basis functions at the quadrature points.
478  *
479  * Dimensions: (sides)
480  * x (no. of quadrature points)
481  * x (no. of dofs)
482  * x ((dim_ of. ref. cell)x(no. of components in ref. cell))
483  */
487 
488  arma::mat grad(dim, fe->n_components());
489  for (unsigned int sid=0; sid<dim+1; sid++) {
490  auto quad = q->make_from_side<dim>(sid);
491  for (unsigned int i_pt=0; i_pt<quad.size(); i_pt++)
492  {
493  for (unsigned int i_dof=0; i_dof<fe->n_dofs(); i_dof++)
494  {
495  grad.zeros();
496  for (unsigned int c=0; c<fe->n_components(); c++)
497  grad.col(c) += fe->shape_grad(i_dof, quad.template point<dim>(i_pt), c);
498 
499  ref_shape_grads[sid][i_pt][i_dof] = grad;
500  }
501  }
502  }
503 
504  return ref_shape_grads;
505  }
506 
508  std::shared_ptr< FiniteElement<dim> > fe_;
509 };
510 
511 
512 template<unsigned int dim>
513 class JoinValues : public BaseValues<dim>
514 {
515 public:
516  /// Constructor
517  JoinValues(PatchPointValues<3> *patch_point_vals_bulk, PatchPointValues<3> *patch_point_vals_side, MixedPtr<FiniteElement> fe)
518  : BaseValues<dim>(), patch_point_vals_bulk_(patch_point_vals_bulk), patch_point_vals_side_(patch_point_vals_side) {
519  ASSERT_EQ(patch_point_vals_bulk->dim(), dim-1);
520  ASSERT_EQ(patch_point_vals_side->dim(), dim);
521  fe_high_dim_ = fe[Dim<dim>{}];
522  fe_low_dim_ = fe[Dim<dim-1>{}];
523  }
524 
526  {
527  // element of lower dim (bulk points)
528  auto fe_component_low = this->fe_comp(fe_low_dim_, component_idx);
529  ASSERT_EQ(fe_component_low->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
530  // use lambda reinit function
531  std::vector< std::vector<double> > shape_values_bulk( patch_point_vals_bulk_->get_quadrature()->size(), vector<double>(fe_component_low->n_dofs()) );
532  auto ref_shape_vals_bulk = this->ref_shape_values_bulk(patch_point_vals_bulk_->get_quadrature(), fe_component_low);
533  for (unsigned int i = 0; i < patch_point_vals_bulk_->get_quadrature()->size(); i++)
534  for (unsigned int j = 0; j < fe_component_low->n_dofs(); j++) {
535  shape_values_bulk[i][j] = ref_shape_vals_bulk[i][j][0];
536  }
537  uint scalar_shape_op_idx_bulk = patch_point_vals_bulk_->operations_.size(); // index in operations_ vector
538  auto lambda_scalar_shape_bulk = [shape_values_bulk, scalar_shape_op_idx_bulk](std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
539  bulk_reinit::ptop_scalar_shape(operations, op_results, shape_values_bulk, scalar_shape_op_idx_bulk);
540  };
541  auto &grad_scalar_shape_bulk_op = patch_point_vals_bulk_->make_fe_op({1}, lambda_scalar_shape_bulk, {}, fe_component_low->n_dofs());
542  uint begin_bulk = grad_scalar_shape_bulk_op.result_row();
543 
544  // element of higher dim (side points)
545  auto fe_component_high = this->fe_comp(fe_high_dim_, component_idx);
546  ASSERT_EQ(fe_component_high->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
547  // use lambda reinit function
548  std::vector< std::vector< std::vector<double> > > shape_values_side(
549  dim+1,
551  );
552  auto ref_shape_vals_side = this->ref_shape_values_side(patch_point_vals_side_->get_quadrature(), fe_component_high);
553  for (unsigned int s=0; s<dim+1; ++s)
554  for (unsigned int i = 0; i < patch_point_vals_side_->get_quadrature()->size(); i++)
555  for (unsigned int j = 0; j < fe_component_high->n_dofs(); j++) {
556  shape_values_side[s][i][j] = ref_shape_vals_side[s][i][j][0];
557  }
558  uint scalar_shape_op_idx_side = patch_point_vals_side_->operations_.size(); // index in operations_ vector
559  auto lambda_scalar_shape_side = [shape_values_side, scalar_shape_op_idx_side](std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
560  side_reinit::ptop_scalar_shape(operations, op_results, el_table, shape_values_side, scalar_shape_op_idx_side);
561  };
562  auto &grad_scalar_shape_side_op = patch_point_vals_side_->make_fe_op({1}, lambda_scalar_shape_side, {}, fe_component_high->n_dofs());
563  uint begin_side = grad_scalar_shape_side_op.result_row();
564 
565  auto bgn_it = make_iter<JoinShapeAccessor<Scalar>>( JoinShapeAccessor<Scalar>(patch_point_vals_bulk_, patch_point_vals_side_,
566  begin_bulk, begin_side, fe_component_low->n_dofs(), fe_component_high->n_dofs(), 0) );
567  unsigned int end_idx = fe_component_low->n_dofs() + fe_component_high->n_dofs();
568  auto end_it = make_iter<JoinShapeAccessor<Scalar>>( JoinShapeAccessor<Scalar>(patch_point_vals_bulk_, patch_point_vals_side_,
569  begin_bulk, begin_side, fe_component_low->n_dofs(), fe_component_high->n_dofs(), end_idx) );
570  return Range<JoinShapeAccessor<Scalar>>(bgn_it, end_it);
571  }
572 
573 private:
576  std::shared_ptr< FiniteElement<dim> > fe_high_dim_;
577  std::shared_ptr< FiniteElement<dim-1> > fe_low_dim_;
578 };
579 
580 /// Template specialization of dim = 1
581 template <>
582 class JoinValues<1> : public BaseValues<1>
583 {
584 public:
585  /// Constructor
587  : BaseValues<1>() {}
588 
590  {
592  make_iter<JoinShapeAccessor<Scalar>>(JoinShapeAccessor<Scalar>()),
594  }
595 };
596 
597 
598 template<unsigned int spacedim = 3>
600 private:
603  SideFE = 1
604  };
605 
607  {
608  public:
610 
611  /// Shape functions evaluated at the quadrature points.
613 
614  /// Gradients of shape functions evaluated at the quadrature points.
615  /// Each row of the matrix contains the gradient of one shape function.
617 
618  /// Auxiliary object for calculation of element-dependent data.
619  std::shared_ptr<ElementValues<spacedim> > elm_values_;
620 
621  };
622 
623  /// Subobject holds FE data of one dimension (0,1,2,3)
625  public:
626  /// Constructor
627  DimPatchFEValues(unsigned int max_size=0)
629 
630 
631  inline unsigned int used_size() const {
632  return used_size_;
633  }
634 
635  inline unsigned int max_size() const {
636  return element_data_.size();
637  }
638 
639  void reinit(PatchElementsList patch_elements);
640 
641  /**
642  * @brief Initialize structures and calculates cell-independent data.
643  *
644  * @param _quadrature The quadrature rule for the cell associated
645  * to given finite element or for the cell side.
646  * @param _fe The finite element.
647  * @param _flags The update flags.
648  */
649  template<unsigned int DIM>
650  void initialize(Quadrature &_quadrature,
651  FiniteElement<DIM> &_fe,
652  UpdateFlags _flags);
653 
654  /**
655  * @brief Allocates space for computed data.
656  *
657  * @param n_points Number of quadrature points.
658  * @param _fe The finite element.
659  * @param flags The update flags.
660  */
661  template<unsigned int DIM>
662  void allocate(Quadrature &_quadrature,
663  FiniteElement<DIM> &_fe,
664  UpdateFlags flags);
665 
666  /// Precompute finite element data on reference element.
667  template<unsigned int DIM>
668  std::shared_ptr<FEInternalData> init_fe_data(const FiniteElement<DIM> &fe, const Quadrature &q);
669 
670  /**
671  * @brief Computes the shape function values and gradients on the actual cell
672  * and fills the FEValues structure.
673  *
674  * @param fe_data Precomputed finite element data.
675  */
676  void fill_data(const ElementValues<spacedim> &elm_values, const FEInternalData &fe_data);
677 
678  /**
679  * @brief Computes the shape function values and gradients on the actual cell
680  * and fills the FEValues structure. Specialized variant of previous method for
681  * different FETypes given by template parameter.
682  */
683  template<class MapType>
684  void fill_data_specialized(const ElementValues<spacedim> &elm_values, const FEInternalData &fe_data);
685 
686  /**
687  * !! Temporary function. Use in fill_data.
688  * Set shape value @p val of the @p i_point and @p i_func_comp.
689  */
690  inline void set_shape_value(unsigned int i_point, unsigned int i_func_comp, double val)
691  {
692  element_data_[patch_data_idx_].shape_values_[i_point][i_func_comp] = val;
693  }
694 
695  /**
696  * !! Temporary function. Use in fill_data.
697  * Set shape gradient @p val of the @p i_point and @p i_func_comp.
698  */
699  inline void set_shape_gradient(unsigned int i_point, unsigned int i_func_comp, arma::vec::fixed<spacedim> val)
700  {
701  element_data_[patch_data_idx_].shape_gradients_[i_point][i_func_comp] = val;
702  }
703 
704  /**
705  * !! Temporary function. Use in fill_data.
706  * @brief Return the value of the @p function_no-th shape function at
707  * the @p point_no-th quadrature point.
708  *
709  * @param function_no Number of the shape function.
710  * @param point_no Number of the quadrature point.
711  */
712  inline double shape_value(const unsigned int function_no, const unsigned int point_no) const
713  {
714  ASSERT_LT(function_no, this->n_dofs_);
715  ASSERT_LT(point_no, this->n_points_);
716  return element_data_[patch_data_idx_].shape_values_[point_no][function_no];
717  }
718 
719  /**
720  * !! Temporary function. Use in fill_data.
721  * @brief Return the gradient of the @p function_no-th shape function at
722  * the @p point_no-th quadrature point.
723  *
724  * @param function_no Number of the shape function.
725  * @param point_no Number of the quadrature point.
726  */
727  inline arma::vec::fixed<spacedim> shape_grad(const unsigned int function_no, const unsigned int point_no) const
728  {
729  ASSERT_LT(function_no, this->n_dofs_);
730  ASSERT_LT(point_no, this->n_points_);
731  return element_data_[patch_data_idx_].shape_gradients_[point_no][function_no];
732  }
733 
734  /**
735  * @brief Return the product of Jacobian determinant and the quadrature
736  * weight at given quadrature point.
737  *
738  * @param p BulkPoint corresponds to the quadrature point.
739  */
740  inline double JxW(const BulkPoint &p)
741  {
742  unsigned int patch_data_idx = element_patch_map_.find(p.elem_patch_idx())->second;
743  return element_data_[patch_data_idx].elm_values_->JxW(p.eval_point_idx());
744  }
745 
746  /**
747  * @brief Return the product of Jacobian determinant and the quadrature
748  * weight at given quadrature point.
749  *
750  * @param p SidePoint corresponds to the quadrature point.
751  */
752  inline double JxW(const SidePoint &p)
753  {
754  unsigned int patch_data_idx = element_patch_map_.find(p.elem_patch_idx())->second + p.side_idx();
755  return element_data_[patch_data_idx].elm_values_->side_JxW(p.local_point_idx());
756  }
757 
758  /**
759  * @brief Return the value of the @p function_no-th shape function at
760  * the @p p quadrature point.
761  *
762  * @param function_no Number of the shape function.
763  * @param p BulkPoint corresponds to the quadrature point.
764  */
765  inline double shape_value(const unsigned int function_no, const BulkPoint &p) const
766  {
767  ASSERT_LT(function_no, this->n_dofs_);
768  unsigned int patch_data_idx = element_patch_map_.find(p.elem_patch_idx())->second;
769  return element_data_[patch_data_idx].shape_values_[p.eval_point_idx()][function_no];
770  }
771 
772  /**
773  * @brief Return the value of the @p function_no-th shape function at
774  * the @p p quadrature point.
775  *
776  * @param function_no Number of the shape function.
777  * @param p SidePoint corresponds to the quadrature point.
778  */
779  inline double shape_value(const unsigned int function_no, const SidePoint &p) const
780  {
781  ASSERT_LT(function_no, this->n_dofs_);
782  unsigned int patch_data_idx = element_patch_map_.find(p.elem_patch_idx())->second + p.side_idx();
783  return element_data_[patch_data_idx].shape_values_[p.local_point_idx()][function_no];
784  }
785 
786  /**
787  * @brief Return the gradient of the @p function_no-th shape function at
788  * the @p p quadrature point.
789  *
790  * @param function_no Number of the shape function.
791  * @param p BulkPoint corresponds to the quadrature point.
792  */
793  inline arma::vec::fixed<spacedim> shape_grad(const unsigned int function_no, const BulkPoint &p) const
794  {
795  ASSERT_LT(function_no, this->n_dofs_);
796  unsigned int patch_data_idx = element_patch_map_.find(p.elem_patch_idx())->second;
797  return element_data_[patch_data_idx].shape_gradients_[p.eval_point_idx()][function_no];;
798  }
799 
800  /**
801  * @brief Return the gradient of the @p function_no-th shape function at
802  * the @p p quadrature point.
803  *
804  * @param function_no Number of the shape function.
805  * @param p SidePoint corresponds to the quadrature point.
806  */
807  inline arma::vec::fixed<spacedim> shape_grad(const unsigned int function_no, const SidePoint &p) const
808  {
809  ASSERT_LT(function_no, this->n_dofs_);
810  unsigned int patch_data_idx = element_patch_map_.find(p.elem_patch_idx())->second + p.side_idx();
811  return element_data_[patch_data_idx].shape_gradients_[p.local_point_idx()][function_no];;
812  }
813 
814  /**
815  * @brief Returns the normal vector to a side at given quadrature point.
816  *
817  * @param p SidePoint corresponds to the quadrature point.
818  */
819  inline arma::vec::fixed<spacedim> normal_vector(const SidePoint &p)
820  {
821  unsigned int patch_data_idx = element_patch_map_.find(p.elem_patch_idx())->second + p.side_idx();
822  return element_data_[patch_data_idx].elm_values_->normal_vector(p.local_point_idx());
823  }
824 
825  /// Set size of ElementFEData. Important: Use only during the initialization of FESystem !
826  void resize(unsigned int max_size) {
827  element_data_.resize(max_size);
828  }
829 
830 
831  /// Dimension of reference space.
832  unsigned int dim_;
833 
834  /// Number of integration points.
835  unsigned int n_points_;
836 
837  /// Number of finite element dofs.
838  unsigned int n_dofs_;
839 
840  /// Type of finite element (scalar, vector, tensor).
842 
843  /// Dof indices of FESystem sub-elements.
845 
846  /// Numbers of components of FESystem sub-elements in reference space.
848 
849  /// Numbers of components of FESystem sub-elements in real space.
851 
852  /// Flags that indicate which finite element quantities are to be computed.
854 
855  /// Vector of FEValues for sub-elements of FESystem.
857 
858  /// Number of components of the FE.
859  unsigned int n_components_;
860 
861 // /// Auxiliary storage of FEValuesViews accessors.
862 // ViewsCache views_cache_;
863 
864  /// Precomputed finite element data.
865  std::shared_ptr<FEInternalData> fe_data_;
866 
867  /// Precomputed FE data (shape functions on reference element) for all side quadrature points.
869 
870  /// Patch index of processed element / side.
871  unsigned int patch_data_idx_;
872 
873  /// Map of element patch indexes to element_data_.
875 
876  /// Data of elements / sides on patch
878 
879  /// Number of elements / sides on patch. Must be less or equal to size of element_data vector
880  unsigned int used_size_;
881 
882  /// Maximal number of elements on patch.
883  unsigned int max_n_elem_;
884 
885  /// Distinguishes using of PatchFEValues for storing data of elements or sides.
887  };
888 
889 public:
890 
893  dim_fe_side_vals_({DimPatchFEValues(0), DimPatchFEValues(0), DimPatchFEValues(0)}),
896  used_quads_[0] = false; used_quads_[1] = false;
897  }
898 
899  PatchFEValues(unsigned int n_quad_points, unsigned int quad_order, MixedPtr<FiniteElement> fe)
900  : dim_fe_vals_({DimPatchFEValues(n_quad_points), DimPatchFEValues(n_quad_points), DimPatchFEValues(n_quad_points)}),
901  dim_fe_side_vals_({DimPatchFEValues(n_quad_points), DimPatchFEValues(n_quad_points), DimPatchFEValues(n_quad_points)}),
903  FeBulk::PatchPointValues(2, quad_order),
904  FeBulk::PatchPointValues(3, quad_order)} },
906  FeSide::PatchPointValues(2, quad_order),
907  FeSide::PatchPointValues(3, quad_order)} },
908  fe_(fe) {
909  used_quads_[0] = false; used_quads_[1] = false;
910  }
911 
912 
913  /// Destructor
915  {}
916 
917  /// Return bulk or side quadrature of given dimension
918  Quadrature *get_quadrature(uint dim, bool is_bulk) const {
919  if (is_bulk) return patch_point_vals_bulk_[dim-1].get_quadrature();
920  else return patch_point_vals_side_[dim-1].get_quadrature();
921  }
922 
923  /**
924  * @brief Initialize structures and calculates cell-independent data.
925  *
926  * @param _quadrature The quadrature rule for the cell associated
927  * to given finite element or for the cell side.
928  * @param _flags The update flags.
929  */
930  template<unsigned int DIM>
931  void initialize(Quadrature &_quadrature,
932  UpdateFlags _flags)
933  {
934  if ( _quadrature.dim() == DIM ) {
935  dim_fe_vals_[DIM-1].initialize(_quadrature, *fe_[Dim<DIM>{}], _flags);
936  used_quads_[0] = true;
937  // new data storing
938  patch_point_vals_bulk_[DIM-1].initialize(3); // bulk
939  } else {
940  dim_fe_side_vals_[DIM-1].initialize(_quadrature, *fe_[Dim<DIM>{}], _flags);
941  used_quads_[1] = true;
942  // new data storing
943  patch_point_vals_side_[DIM-1].initialize(4); // side
944  }
945  }
946 
947  /// Reset PatchpointValues structures
948  void reset()
949  {
950  for (unsigned int i=0; i<3; ++i) {
951  if (used_quads_[0]) patch_point_vals_bulk_[i].reset();
952  if (used_quads_[1]) patch_point_vals_side_[i].reset();
953  }
954  }
955 
956  /// Reinit data - old data storing, temporary
957  void reinit(std::array<PatchElementsList, 4> patch_elements)
958  {
959  for (unsigned int i=0; i<3; ++i) {
960  if (used_quads_[0]) dim_fe_vals_[i].reinit(patch_elements[i+1]);
961  if (used_quads_[1]) dim_fe_side_vals_[i].reinit(patch_elements[i+1]);
962  }
963  }
964 
965  /// Reinit data.
967  {
968  for (unsigned int i=0; i<3; ++i) {
969  if (used_quads_[0]) patch_point_vals_bulk_[i].reinit_patch();
970  if (used_quads_[1]) patch_point_vals_side_[i].reinit_patch();
971  }
972  }
973 
974  /**
975  * @brief Returns the number of shape functions.
976  */
977  template<unsigned int dim>
978  inline unsigned int n_dofs() const {
979  ASSERT((dim>=0) && (dim<=3))(dim).error("Dimension must be 0, 1, 2 or 3.");
980  return fe_[Dim<dim>{}]->n_dofs();
981  }
982 
983  /// Return BulkValue object of dimension given by template parameter
984  template<unsigned int dim>
986  ASSERT((dim>0) && (dim<=3))(dim).error("Dimension must be 1, 2 or 3.");
988  }
989 
990  /// Return SideValue object of dimension given by template parameter
991  template<unsigned int dim>
993  ASSERT((dim>0) && (dim<=3))(dim).error("Dimension must be 1, 2 or 3.");
995  }
996 
997  /// Return JoinValue object of dimension given by template parameter
998  template<unsigned int dim>
1000  ASSERT((dim>1) && (dim<=3))(dim).error("Dimension must be 2 or 3.");
1002  }
1003 
1004  /// Resize \p dim_point_table_ if actual size is less than new_size and return reference
1005  inline DimPointTable &dim_point_table(unsigned int new_size) {
1006  if (dim_point_table_.size() < new_size) dim_point_table_.resize(new_size);
1007  for (uint i=0; i<new_size; ++i) dim_point_table_[i][0] = 10; // set invalid dim
1008  return dim_point_table_;
1009  }
1010 
1011  /** Following methods are used during update of patch. **/
1012 
1013  /// Resize tables of patch_point_vals_
1015  ASSERT_EQ(dim_sizes.size(), 2);
1016  ASSERT_EQ(dim_sizes[0].size(), 3);
1017 
1018  for (uint i=0; i<3; ++i) {
1019  patch_point_vals_bulk_[i].resize_tables(dim_sizes[0][i]);
1020  patch_point_vals_side_[i].resize_tables(dim_sizes[1][i]);
1021  }
1022  }
1023 
1024  /// Register element to patch_point_vals_ table by dimension of element
1025  uint register_element(DHCellAccessor cell, uint element_patch_idx) {
1026  arma::mat coords;
1027  switch (cell.dim()) {
1028  case 1:
1029  coords = MappingP1<1,spacedim>::element_map(cell.elm());
1030  return patch_point_vals_bulk_[0].register_element(coords, element_patch_idx);
1031  break;
1032  case 2:
1033  coords = MappingP1<2,spacedim>::element_map(cell.elm());
1034  return patch_point_vals_bulk_[1].register_element(coords, element_patch_idx);
1035  break;
1036  case 3:
1037  coords = MappingP1<3,spacedim>::element_map(cell.elm());
1038  return patch_point_vals_bulk_[2].register_element(coords, element_patch_idx);
1039  break;
1040  default:
1041  ASSERT(false);
1042  return 0;
1043  break;
1044  }
1045  }
1046 
1047  /// Register side to patch_point_vals_ table by dimension of side
1049  arma::mat side_coords(spacedim, cell_side.dim());
1050  for (unsigned int n=0; n<cell_side.dim(); n++)
1051  for (unsigned int c=0; c<spacedim; c++)
1052  side_coords(c,n) = (*cell_side.side().node(n))[c];
1053 
1054  arma::mat elm_coords;
1055  DHCellAccessor cell = cell_side.cell();
1056  switch (cell.dim()) {
1057  case 1:
1058  elm_coords = MappingP1<1,spacedim>::element_map(cell.elm());
1059  return patch_point_vals_side_[0].register_side(elm_coords, side_coords);
1060  break;
1061  case 2:
1062  elm_coords = MappingP1<2,spacedim>::element_map(cell.elm());
1063  return patch_point_vals_side_[1].register_side(elm_coords, side_coords);
1064  break;
1065  case 3:
1066  elm_coords = MappingP1<3,spacedim>::element_map(cell.elm());
1067  return patch_point_vals_side_[2].register_side(elm_coords, side_coords);
1068  break;
1069  default:
1070  ASSERT(false);
1071  return 0;
1072  break;
1073  }
1074  }
1075 
1076  /// Register bulk point to patch_point_vals_ table by dimension of element
1077  uint register_bulk_point(DHCellAccessor cell, uint elem_table_row, uint value_patch_idx) {
1078  return patch_point_vals_bulk_[cell.dim()-1].register_bulk_point(elem_table_row, value_patch_idx, cell.elm_idx());
1079  }
1080 
1081  /// Register side point to patch_point_vals_ table by dimension of side
1082  uint register_side_point(DHCellSide cell_side, uint elem_table_row, uint value_patch_idx) {
1083  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());
1084  }
1085 
1086  /// Temporary development method
1087  void print(bool points, bool ints, bool only_bulk=true) const {
1088  for (uint i=0; i<3; ++i)
1089  patch_point_vals_bulk_[i].print(points, ints);
1090  if (!only_bulk)
1091  for (uint i=0; i<3; ++i)
1092  patch_point_vals_side_[i].print(points, ints);
1093  }
1094 
1095  /// Temporary development method
1096  void print_operations(ostream& stream) const {
1097  stream << endl << "Table of patch FE operations:" << endl;
1098  for (uint i=0; i<3; ++i) {
1099  stream << std::setfill('-') << setw(100) << "" << endl;
1100  stream << "Bulk, dimension " << (i+1) << ", n_rows " << patch_point_vals_bulk_[i].n_rows() << endl;
1101  patch_point_vals_bulk_[i].print_operations(stream, 0);
1102  }
1103  for (uint i=0; i<3; ++i) {
1104  stream << std::setfill('-') << setw(100) << "" << endl;
1105  stream << "Side, dimension " << (i+1) << ", n_rows " << patch_point_vals_side_[i].n_rows() << endl;
1106  patch_point_vals_side_[i].print_operations(stream, 1);
1107  }
1108  stream << std::setfill('=') << setw(100) << "" << endl;
1109  }
1110 
1111 private:
1112  /// Sub objects of dimensions 1,2,3
1113  std::array<DimPatchFEValues, 3> dim_fe_vals_;
1114  std::array<DimPatchFEValues, 3> dim_fe_side_vals_;
1115  std::array<FeBulk::PatchPointValues<spacedim>, 3> patch_point_vals_bulk_;
1116  std::array<FeSide::PatchPointValues<spacedim>, 3> patch_point_vals_side_;
1118 
1119  MixedPtr<FiniteElement> fe_; ///< Mixed of shared pointers of FiniteElement object
1120 
1121  ///< Temporary helper objects used in step between usage old a new implementation
1122  bool used_quads_[2];
1123 
1124  template <class ValueType>
1125  friend class ElQ;
1126  template <class ValueType>
1127  friend class FeQ;
1128  template <class ValueType>
1129  friend class JoinShapeAccessor;
1130 };
1131 
1132 
1133 template <class ValueType>
1134 ValueType ElQ<ValueType>::operator()(const BulkPoint &point) {
1135  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1136  return patch_point_vals_.scalar_val(begin_, value_cache_idx);
1137 }
1138 
1139 template <>
1141  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1142  return patch_point_vals_.vector_val(begin_, value_cache_idx);
1143 }
1144 
1145 template <>
1147  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1148  return patch_point_vals_.tensor_val(begin_, value_cache_idx);
1149 }
1150 
1151 template <class ValueType>
1152 ValueType ElQ<ValueType>::operator()(const SidePoint &point) {
1153  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1154  return patch_point_vals_.scalar_val(begin_, value_cache_idx);
1155 }
1156 
1157 template <>
1159  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1160  return patch_point_vals_.vector_val(begin_, value_cache_idx);
1161 }
1162 
1163 template <>
1165  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1166  return patch_point_vals_.tensor_val(begin_, value_cache_idx);
1167 }
1168 
1169 template <class ValueType>
1170 ValueType FeQ<ValueType>::operator()(unsigned int shape_idx, const BulkPoint &point) {
1171  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1172  return patch_point_vals_.scalar_val(begin_+shape_idx, value_cache_idx);
1173 }
1174 
1175 template <>
1176 inline Vector FeQ<Vector>::operator()(unsigned int shape_idx, const BulkPoint &point) {
1177  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1178  return patch_point_vals_.vector_val(begin_+3*shape_idx, value_cache_idx);
1179 }
1180 
1181 template <>
1182 inline Tensor FeQ<Tensor>::operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const BulkPoint &point) {
1183  Tensor tens; tens.zeros();
1184  return tens;
1185 }
1186 
1187 template <class ValueType>
1188 ValueType FeQ<ValueType>::operator()(unsigned int shape_idx, const SidePoint &point) {
1189  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1190  return patch_point_vals_.scalar_val(begin_+shape_idx, value_cache_idx);
1191 }
1192 
1193 template <>
1194 inline Vector FeQ<Vector>::operator()(unsigned int shape_idx, const SidePoint &point) {
1195  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1196  return patch_point_vals_.vector_val(begin_+3*shape_idx, value_cache_idx);
1197 }
1198 
1199 template <>
1200 inline Tensor FeQ<Tensor>::operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const SidePoint &point) {
1201  Tensor tens; tens.zeros();
1202  return tens;
1203 }
1204 
1205 
1206 template <class ValueType>
1208  if (this->is_high_dim()) {
1209  return 0.0;
1210  } else {
1211  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1212  return patch_point_vals_bulk_->scalar_val(begin_+this->local_idx(), value_cache_idx);
1213  }
1214 }
1215 
1216 template <>
1218  Vector vect; vect.zeros();
1219  return vect;
1220 }
1221 
1222 template <>
1224  Tensor tens; tens.zeros();
1225  return tens;
1226 }
1227 
1228 template <class ValueType>
1230  if (this->is_high_dim()) {
1231  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1232  return patch_point_vals_side_->scalar_val(begin_side_+this->local_idx(), value_cache_idx);
1233  } else {
1234  return 0.0;
1235  }
1236 }
1237 
1238 template <>
1240  Vector vect; vect.zeros();
1241  return vect;
1242 }
1243 
1244 template <>
1246  Tensor tens; tens.zeros();
1247  return tens;
1248 }
1249 
1250 
1251 #endif /* PATCH_FE_VALUES_HH_ */
#define ASSERT(expr)
Definition: asserts.hh:351
#define ASSERT_LT(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:301
#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()=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)
ElQ(PatchPointValues< 3 > &patch_point_vals, unsigned int begin)
Constructor.
int element_eval_point(unsigned int i_elem_in_cache, unsigned int i_eval_point) const
Class for computation of data on cell and side.
Structure for storing the precomputed finite element data.
Definition: fe_values.hh:62
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
Subobject holds FE data of one dimension (0,1,2,3)
arma::vec::fixed< spacedim > normal_vector(const SidePoint &p)
Returns the normal vector to a side at given quadrature point.
void fill_data_specialized(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure....
unsigned int max_n_elem_
Maximal number of elements on patch.
unsigned int n_dofs_
Number of finite element dofs.
std::shared_ptr< FEInternalData > init_fe_data(const FiniteElement< DIM > &fe, const Quadrature &q)
Precompute finite element data on reference element.
unsigned int used_size_
Number of elements / sides on patch. Must be less or equal to size of element_data vector.
double shape_value(const unsigned int function_no, const BulkPoint &p) const
Return the value of the function_no-th shape function at the p quadrature point.
double shape_value(const unsigned int function_no, const unsigned int point_no) const
Return the value of the function_no-th shape function at the point_no-th quadrature point.
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const BulkPoint &p) const
Return the gradient of the function_no-th shape function at the p quadrature point.
void resize(unsigned int max_size)
Set size of ElementFEData. Important: Use only during the initialization of FESystem !
std::vector< ElementFEData > element_data_
Data of elements / sides on patch.
double JxW(const BulkPoint &p)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
std::map< unsigned int, unsigned int > element_patch_map_
Map of element patch indexes to element_data_.
DimPatchFEValues(unsigned int max_size=0)
Constructor.
std::vector< unsigned int > fe_sys_n_components_
Numbers of components of FESystem sub-elements in reference space.
void allocate(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags flags)
Allocates space for computed data.
void fill_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure.
void set_shape_gradient(unsigned int i_point, unsigned int i_func_comp, arma::vec::fixed< spacedim > val)
void reinit(PatchElementsList patch_elements)
std::vector< unsigned int > fe_sys_n_space_components_
Numbers of components of FESystem sub-elements in real space.
void set_shape_value(unsigned int i_point, unsigned int i_func_comp, double val)
void initialize(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
unsigned int n_components_
Number of components of the FE.
MeshObjectType object_type_
Distinguishes using of PatchFEValues for storing data of elements or sides.
double JxW(const SidePoint &p)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
unsigned int n_points_
Number of integration points.
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const unsigned int point_no) const
Return the gradient of the function_no-th shape function at the point_no-th quadrature point.
std::vector< shared_ptr< FEInternalData > > side_fe_data_
Precomputed FE data (shape functions on reference element) for all side quadrature points.
std::shared_ptr< FEInternalData > fe_data_
Precomputed finite element data.
FEType fe_type_
Type of finite element (scalar, vector, tensor).
unsigned int patch_data_idx_
Patch index of processed element / side.
unsigned int dim_
Dimension of reference space.
std::vector< std::vector< unsigned int > > fe_sys_dofs_
Dof indices of FESystem sub-elements.
std::vector< PatchFEValues< spacedim >::DimPatchFEValues > fe_values_vec
Vector of FEValues for sub-elements of FESystem.
double shape_value(const unsigned int function_no, const SidePoint &p) const
Return the value of the function_no-th shape function at the p quadrature point.
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const SidePoint &p) const
Return the gradient of the function_no-th shape function at the p quadrature point.
std::vector< std::vector< arma::vec::fixed< spacedim > > > shape_gradients_
std::shared_ptr< ElementValues< spacedim > > elm_values_
Auxiliary object for calculation of element-dependent data.
std::vector< std::vector< double > > shape_values_
Shape functions evaluated at the quadrature points.
std::array< DimPatchFEValues, 3 > dim_fe_vals_
Sub objects of dimensions 1,2,3.
uint register_element(DHCellAccessor cell, uint element_patch_idx)
Register element to patch_point_vals_ table by dimension of element.
std::array< DimPatchFEValues, 3 > dim_fe_side_vals_
void print(bool points, bool ints, bool only_bulk=true) const
Temporary development method.
~PatchFEValues()
Destructor.
std::array< FeBulk::PatchPointValues< spacedim >, 3 > patch_point_vals_bulk_
uint register_bulk_point(DHCellAccessor cell, uint elem_table_row, uint value_patch_idx)
Register bulk point to patch_point_vals_ table by dimension of element.
uint register_side(DHCellSide cell_side)
Register side to patch_point_vals_ table by dimension of side.
DimPointTable dim_point_table_
void initialize(Quadrature &_quadrature, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
JoinValues< dim > join_values()
Return JoinValue object of dimension given by template parameter.
void reinit(std::array< PatchElementsList, 4 > patch_elements)
Reinit data - old data storing, temporary.
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.
void print_operations(ostream &stream) const
Temporary development method.
void resize_tables(std::vector< std::vector< uint > > dim_sizes)
Resize tables of patch_point_vals_.
std::array< FeSide::PatchPointValues< spacedim >, 3 > patch_point_vals_side_
DimPointTable & dim_point_table(unsigned int new_size)
Resize dim_point_table_ if actual size is less than new_size and return reference.
MixedPtr< FiniteElement > fe_
Mixed of shared pointers of FiniteElement object.
SideValues< dim > side_values()
Return SideValue object of dimension given by template parameter.
void reinit_patch()
Reinit data.
PatchFEValues(unsigned int n_quad_points, unsigned int quad_order, MixedPtr< FiniteElement > fe)
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.
std::vector< ElOp< spacedim > > operations_
Vector of all defined operations.
Quadrature * get_quadrature() const
Return quadrature.
uint dim() const
Getter of dim_.
ElOp< spacedim > & make_fe_op(std::initializer_list< uint > shape, ReinitFunction reinit_f, std::vector< uint > input_ops_vec, uint n_dofs)
Add accessor to operations_ vector.
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
virtual unsigned int side_idx() const =0
Intermediate step in implementation of PatcFEValues.
unsigned int eval_point_idx() const
Return index in EvalPoints object.
Definition: eval_subset.hh:143
unsigned int local_point_idx() const
Return local index in quadrature. Temporary method - intermediate step in implementation of PatcFEVal...
Definition: eval_subset.hh:135
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:33
Eigen::Vector< ArrayInt, Eigen::Dynamic > TableInt
Definition: eigen_tools.hh:34
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,...
FEType
@ 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
static void ptop_scalar_shape(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, 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.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
Definition: update_flags.hh:68