Flow123d  DF_patch_fe_mechanics-b866f5c
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 
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 op_idx)
58  : patch_point_vals_(patch_point_vals), op_idx_(op_idx) {}
59 
60  ValueType operator()(const BulkPoint &point) const;
61 
62  ValueType operator()(const SidePoint &point) const;
63 
64 private:
65  PatchPointValues<3> &patch_point_vals_; ///< Reference to PatchPointValues
66  unsigned int op_idx_; ///< Index of operation in patch_point_vals_.operations vector
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 op_idx, unsigned int n_dofs)
78  : patch_point_vals_(patch_point_vals), op_idx_(op_idx), n_dofs_(n_dofs) {}
79 
80 
81  ValueType operator()(unsigned int shape_idx, FMT_UNUSED const BulkPoint &point) const;
82 
83  ValueType operator()(unsigned int shape_idx, FMT_UNUSED const SidePoint &point) const;
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 op_idx_; ///< Index of operation in patch_point_vals_.operations vector
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, unsigned int n_dofs_bulk,
114  unsigned int n_dofs_side, unsigned int op_idx_bulk, unsigned int op_idx_side, unsigned int join_idx)
115  : patch_point_vals_bulk_(patch_point_vals_bulk), patch_point_vals_side_(patch_point_vals_side),
116  n_dofs_high_(n_dofs_side), n_dofs_low_(n_dofs_bulk), op_idx_bulk_(op_idx_bulk), op_idx_side_(op_idx_side), join_idx_(join_idx) {
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) const;
159 
160  ValueType operator()(const SidePoint &point) const;
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 n_dofs_high_; ///< Number of DOFs on high-dim element
167  unsigned int n_dofs_low_; ///< Number of DOFs on low-dim element
168  unsigned int op_idx_bulk_; ///< Index of operation in patch_point_vals_bulk_.operations vector
169  unsigned int op_idx_side_; ///< Index of operation in patch_point_vals_side_.operations vector
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  if (fe->fe_type() == FEMixedSystem) {
186  FESystem<FE_dim> *fe_sys = dynamic_cast<FESystem<FE_dim>*>( fe.get() );
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<FE_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  * @brief Precomputed gradients of basis functions at the bulk quadrature points.
254  *
255  * Dimensions: (no. of quadrature points)
256  * x (no. of dofs)
257  * x ((dim_ of. ref. cell)x(no. of components in ref. cell))
258  */
259  template<unsigned int FE_dim>
261  std::vector<std::vector<arma::mat> > ref_shape_grads( q->size(), vector<arma::mat>(fe->n_dofs()) );
262 
263  arma::mat grad(FE_dim, fe->n_components());
264  for (unsigned int i_pt=0; i_pt<q->size(); i_pt++)
265  {
266  for (unsigned int i_dof=0; i_dof<fe->n_dofs(); i_dof++)
267  {
268  grad.zeros();
269  for (unsigned int c=0; c<fe->n_components(); c++)
270  grad.col(c) += fe->shape_grad(i_dof, q->point<FE_dim>(i_pt), c);
271 
272  ref_shape_grads[i_pt][i_dof] = grad;
273  }
274  }
275 
276  return ref_shape_grads;
277  }
278 
279  /**
280  * @brief Precomputed gradients of basis functions at the side quadrature points.
281  *
282  * Dimensions: (sides)
283  * x (no. of quadrature points)
284  * x (no. of dofs)
285  * x ((dim_ of. ref. cell)x(no. of components in ref. cell))
286  */
287  template<unsigned int FE_dim>
289  std::vector<std::vector<std::vector<arma::mat> > > ref_shape_grads( FE_dim+1, std::vector<std::vector<arma::mat> >(q->size(), vector<arma::mat>(fe->n_dofs())) );
290 
291  arma::mat grad(dim, fe->n_components());
292  for (unsigned int sid=0; sid<FE_dim+1; sid++) {
293  auto quad = q->make_from_side<FE_dim>(sid);
294  for (unsigned int i_pt=0; i_pt<quad.size(); i_pt++)
295  {
296  for (unsigned int i_dof=0; i_dof<fe->n_dofs(); i_dof++)
297  {
298  grad.zeros();
299  for (unsigned int c=0; c<fe->n_components(); c++)
300  grad.col(c) += fe->shape_grad(i_dof, quad.template point<FE_dim>(i_pt), c);
301 
302  ref_shape_grads[sid][i_pt][i_dof] = grad;
303  }
304  }
305  }
306 
307  return ref_shape_grads;
308  }
309 
310 };
311 
312 template<unsigned int dim>
313 class BulkValues : public BaseValues<dim>
314 {
315 public:
316  /// Constructor
318  : BaseValues<dim>(), patch_point_vals_(patch_point_vals) {
319  ASSERT_EQ(patch_point_vals.dim(), dim);
320  fe_ = fe[Dim<dim>{}];
321  }
322 
323  /**
324  * @brief Register the product of Jacobian determinant and the quadrature
325  * weight at bulk quadrature points.
326  *
327  * @param quad Quadrature.
328  */
329  inline ElQ<Scalar> JxW()
330  {
332  }
333 
334  /// Create bulk accessor of coords entity
336  {
338  }
339 
340 // inline ElQ<Tensor> jacobian(std::initializer_list<Quadrature *> quad_list)
341 // {}
342 
343  /// Create bulk accessor of jac determinant entity
345  {
347  }
348 
349  /**
350  * @brief Return the value of the @p function_no-th shape function at
351  * the @p p bulk quadrature point.
352  *
353  * @param component_idx Number of the shape function.
354  */
355  inline FeQ<Scalar> scalar_shape(uint component_idx = 0)
356  {
357  auto fe_component = this->fe_comp(fe_, component_idx);
358  ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
359 
360  // use lambda reinit function
361  std::vector< std::vector<double> > shape_values( fe_component->n_dofs(), vector<double>(patch_point_vals_.get_quadrature()->size()) );
362  auto ref_shape_vals = this->ref_shape_values_bulk(patch_point_vals_.get_quadrature(), fe_component);
363  for (unsigned int i = 0; i < patch_point_vals_.get_quadrature()->size(); i++)
364  for (unsigned int j = 0; j < fe_component->n_dofs(); j++) {
365  shape_values[j][i] = ref_shape_vals[i][j][0];
366  }
367  uint scalar_shape_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
368  auto lambda_scalar_shape = [shape_values, scalar_shape_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
369  bulk_reinit::ptop_scalar_shape(operations, shape_values, scalar_shape_op_idx);
370  };
371  patch_point_vals_.make_fe_op({1}, lambda_scalar_shape, {}, fe_component->n_dofs());
372  patch_point_vals_.set_fe_op(FEOps::opScalarShape, scalar_shape_op_idx);
373 
374  return FeQ<Scalar>(patch_point_vals_, scalar_shape_op_idx, fe_component->n_dofs());
375  }
376 
377  inline FeQ<Vector> vector_shape(uint component_idx = 0)
378  {
379  auto fe_component = this->fe_comp(fe_, component_idx);
380 
381  // use lambda reinit function
382  std::vector< std::vector<arma::vec3> > shape_values( fe_component->n_dofs(), vector<arma::vec3>(patch_point_vals_.get_quadrature()->size()) );
383  auto ref_shape_vals = this->ref_shape_values_bulk(patch_point_vals_.get_quadrature(), fe_component);
384  for (unsigned int i = 0; i < patch_point_vals_.get_quadrature()->size(); i++)
385  for (unsigned int j = 0; j < fe_component->n_dofs(); j++) {
386  shape_values[j][i] = ref_shape_vals[i][j];
387  }
388  uint vector_shape_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
389 
390  switch (fe_component->fe_type()) {
391  case FEVector:
392  {
393  auto lambda_vector_shape = [shape_values, vector_shape_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
394  bulk_reinit::ptop_vector_shape(operations, shape_values, vector_shape_op_idx);
395  };
396  patch_point_vals_.make_fe_op({3}, lambda_vector_shape, {}, fe_component->n_dofs());
397  break;
398  }
400  {
401  ASSERT_PERMANENT(false).error("Shape vector for FEVectorContravariant is not implemented yet!\n"); // temporary assert
402  auto lambda_contravariant_shape = [shape_values, vector_shape_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
403  bulk_reinit::ptop_vector_contravariant_shape(operations, shape_values, vector_shape_op_idx);
404  };
405  patch_point_vals_.make_fe_op({3}, lambda_contravariant_shape, {FeBulk::BulkOps::opJac}, fe_component->n_dofs());
406  break;
407  }
408  case FEVectorPiola:
409  {
410  ASSERT_PERMANENT(false).error("Shape vector for FEVectorPiola is not implemented yet!\n"); // temporary assert
411  auto lambda_piola_shape = [shape_values, vector_shape_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
412  bulk_reinit::ptop_vector_piola_shape(operations, shape_values, vector_shape_op_idx);
413  };
414  patch_point_vals_.make_fe_op({3}, lambda_piola_shape, {FeBulk::BulkOps::opJac, FeBulk::BulkOps::opJacDet}, fe_component->n_dofs());
415  break;
416  }
417  default:
418  ASSERT(false).error("Type of FiniteElement of grad_vector_shape accessor must be FEVector, FEVectorPiola or FEVectorContravariant!\n");
419  }
420  patch_point_vals_.set_fe_op(FEOps::opVectorShape, vector_shape_op_idx);
421 
422  return FeQ<Vector>(patch_point_vals_, vector_shape_op_idx, fe_component->n_dofs());
423  }
424 
425 // inline FeQ<Tensor> tensor_shape(uint component_idx = 0)
426 // {}
427 
428  /**
429  * @brief Return the value of the @p function_no-th gradient shape function at
430  * the @p p bulk quadrature point.
431  *
432  * @param component_idx Number of the shape function.
433  */
434  inline FeQ<Vector> grad_scalar_shape(uint component_idx=0)
435  {
436  auto fe_component = this->fe_comp(fe_, component_idx);
437  ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of grad_scalar_shape accessor must be FEScalar!\n");
438 
439  // use lambda reinit function
440  auto ref_shape_grads = this->ref_shape_gradients_bulk(patch_point_vals_.get_quadrature(), fe_component);
441  uint scalar_shape_grads_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
442  auto lambda_scalar_shape_grad = [ref_shape_grads, scalar_shape_grads_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
443  bulk_reinit::ptop_scalar_shape_grads<dim>(operations, ref_shape_grads, scalar_shape_grads_op_idx);
444  };
445  patch_point_vals_.make_fe_op({3}, lambda_scalar_shape_grad, {FeBulk::BulkOps::opInvJac}, fe_component->n_dofs());
446  patch_point_vals_.set_fe_op(FEOps::opGradScalarShape, scalar_shape_grads_op_idx);
447 
448  return FeQ<Vector>(patch_point_vals_, scalar_shape_grads_op_idx, fe_component->n_dofs());
449  }
450 
451  /**
452  * @brief Return the value of the @p function_no-th gradient vector shape function
453  * at the @p p bulk quadrature point.
454  *
455  * @param component_idx Number of the shape function.
456  */
457  inline FeQ<Tensor> grad_vector_shape(uint component_idx=0)
458  {
459  auto fe_component = this->fe_comp(fe_, component_idx);
460 
461  // use lambda reinit function
462  auto ref_shape_grads = this->ref_shape_gradients_bulk(patch_point_vals_.get_quadrature(), fe_component);
463  uint vector_shape_grads_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
464 
465  switch (fe_component->fe_type()) {
466  case FEVector:
467  {
468  auto lambda_vector_shape_grad = [ref_shape_grads, vector_shape_grads_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
469  bulk_reinit::ptop_vector_shape_grads<dim>(operations, ref_shape_grads, vector_shape_grads_op_idx);
470  };
472  lambda_vector_shape_grad,
474  fe_component->n_dofs());
475  break;
476  }
478  {
479  ASSERT_PERMANENT(false).error("Grad vector for FEVectorContravariant is not implemented yet!\n"); // temporary assert
480  auto lambda_contravariant_shape_grad = [ref_shape_grads, vector_shape_grads_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
481  bulk_reinit::ptop_vector_contravariant_shape_grads<dim>(operations, ref_shape_grads, vector_shape_grads_op_idx);
482  };
484  lambda_contravariant_shape_grad,
486  fe_component->n_dofs());
487  break;
488  }
489  case FEVectorPiola:
490  {
491  ASSERT_PERMANENT(false).error("Grad vector for FEVectorPiola is not implemented yet!\n"); // temporary assert
492  auto lambda_piola_shape_grad = [ref_shape_grads, vector_shape_grads_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
493  bulk_reinit::ptop_vector_piola_shape_grads<dim>(operations, ref_shape_grads, vector_shape_grads_op_idx);
494  };
496  lambda_piola_shape_grad,
498  fe_component->n_dofs());
499  break;
500  }
501  default:
502  ASSERT(false).error("Type of FiniteElement of grad_vector_shape accessor must be FEVector, FEVectorPiola or FEVectorContravariant!\n");
503  }
504  patch_point_vals_.set_fe_op(FEOps::opGradVectorShape, vector_shape_grads_op_idx);
505 
506  return FeQ<Tensor>(patch_point_vals_, vector_shape_grads_op_idx, fe_component->n_dofs());
507  }
508 
509  /**
510  * @brief Return the value of the @p function_no-th vector symmetric gradient
511  * at the @p p bulk quadrature point.
512  *
513  * @param component_idx Number of the shape function.
514  */
515  inline FeQ<Tensor> vector_sym_grad(uint component_idx=0)
516  {
517  auto fe_component = this->fe_comp(fe_, component_idx);
518  //ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of grad_scalar_shape accessor must be FEScalar!\n");
519 
520  // use lambda reinit function
521  uint vector_sym_grad_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
523  auto lambda_vector_sym_grad = [vector_sym_grad_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
524  common_reinit::ptop_vector_sym_grad(operations, vector_sym_grad_op_idx);
525  };
526  patch_point_vals_.make_fe_op({3,3}, lambda_vector_sym_grad, {grad_vector_op_idx}, fe_component->n_dofs());
527  patch_point_vals_.set_fe_op(FEOps::opVectorSymGrad, vector_sym_grad_op_idx);
528 
529  return FeQ<Tensor>(patch_point_vals_, vector_sym_grad_op_idx, fe_component->n_dofs());
530  }
531 
532  /**
533  * @brief Return the value of the @p function_no-th vector divergence at
534  * the @p p bulk quadrature point.
535  *
536  * @param component_idx Number of the shape function.
537  */
538  inline FeQ<Scalar> vector_divergence(uint component_idx=0)
539  {
540  auto fe_component = this->fe_comp(fe_, component_idx);
541  //ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of grad_scalar_shape accessor must be FEScalar!\n");
542 
543  // use lambda reinit function
544  uint vector_divergence_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
546  auto lambda_vector_divergence = [vector_divergence_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
547  common_reinit::ptop_vector_divergence(operations, vector_divergence_op_idx);
548  };
549  patch_point_vals_.make_fe_op({1}, lambda_vector_divergence, {grad_vector_op_idx}, fe_component->n_dofs());
550  patch_point_vals_.set_fe_op(FEOps::opVectorDivergence, vector_divergence_op_idx);
551 
552  return FeQ<Scalar>(patch_point_vals_, vector_divergence_op_idx, fe_component->n_dofs());
553  }
554 
555 private:
557  std::shared_ptr< FiniteElement<dim> > fe_;
558 };
559 
560 
561 template<unsigned int dim>
562 class SideValues : public BaseValues<dim>
563 {
564 public:
565  /// Constructor
567  : BaseValues<dim>(), patch_point_vals_(patch_point_vals) {
568  ASSERT_EQ(patch_point_vals.dim(), dim);
569  fe_ = fe[Dim<dim>{}];
570  }
571 
572  /// Same as BulkValues::JxW but register at side quadrature points.
573  inline ElQ<Scalar> JxW()
574  {
576  }
577 
578  /**
579  * @brief Register the normal vector to a side at side quadrature points.
580  *
581  * @param quad Quadrature.
582  */
584  {
586  }
587 
588  /// Create side accessor of coords entity
590  {
592  }
593 
594  /// Create bulk accessor of jac determinant entity
596  {
598  }
599 
600  /// Same as BulkValues::scalar_shape but register at side quadrature points.
601  inline FeQ<Scalar> scalar_shape(uint component_idx = 0)
602  {
603  auto fe_component = this->fe_comp(fe_, component_idx);
604  ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
605 
606  // use lambda reinit function
608  dim+1,
610  );
611  auto ref_shape_vals = this->ref_shape_values_side(patch_point_vals_.get_quadrature(), fe_component);
612  for (unsigned int s=0; s<dim+1; ++s)
613  for (unsigned int i = 0; i < patch_point_vals_.get_quadrature()->size(); i++)
614  for (unsigned int j = 0; j < fe_component->n_dofs(); j++) {
615  shape_values[s][i][j] = ref_shape_vals[s][i][j][0];
616  }
617  uint scalar_shape_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
618  auto lambda_scalar_shape = [shape_values, scalar_shape_op_idx](std::vector<ElOp<3>> &operations, IntTableArena &el_table) {
619  side_reinit::ptop_scalar_shape(operations, el_table, shape_values, scalar_shape_op_idx);
620  };
621  patch_point_vals_.make_fe_op({1}, lambda_scalar_shape, {}, fe_component->n_dofs());
622  patch_point_vals_.set_fe_op(FEOps::opScalarShape, scalar_shape_op_idx);
623 
624  return FeQ<Scalar>(patch_point_vals_, scalar_shape_op_idx, fe_component->n_dofs());
625  }
626 
627  /// Same as BulkValues::vector_shape but register at side quadrature points.
628  inline FeQ<Vector> vector_shape(uint component_idx = 0)
629  {
630  auto fe_component = this->fe_comp(fe_, component_idx);
631  //ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
632 
633  // use lambda reinit function
635  dim+1,
637  );
638  auto ref_shape_vals = this->ref_shape_values_side(patch_point_vals_.get_quadrature(), fe_component);
639  for (unsigned int s=0; s<dim+1; ++s)
640  for (unsigned int i = 0; i < patch_point_vals_.get_quadrature()->size(); i++)
641  for (unsigned int j = 0; j < fe_component->n_dofs(); j++) {
642  shape_values[s][i][j] = ref_shape_vals[s][i][j];
643  }
644  uint vector_shape_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
645 
646  switch (fe_component->fe_type()) {
647  case FEVector:
648  {
649  auto lambda_vector_shape = [shape_values, vector_shape_op_idx](std::vector<ElOp<3>> &operations, IntTableArena &el_table) {
650  side_reinit::ptop_vector_shape(operations, el_table, shape_values, vector_shape_op_idx);
651  };
652  patch_point_vals_.make_fe_op({3}, lambda_vector_shape, {}, fe_component->n_dofs());
653  patch_point_vals_.set_fe_op(FEOps::opVectorShape, vector_shape_op_idx);
654  break;
655  }
657  {
658  ASSERT_PERMANENT(false).error("Shape vector for FEVectorContravariant is not implemented yet!\n"); // temporary assert
659  auto lambda_contravariant_shape = [shape_values, vector_shape_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
660  side_reinit::ptop_vector_contravariant_shape(operations, el_table, shape_values, vector_shape_op_idx);
661  };
662  patch_point_vals_.make_fe_op({3}, lambda_contravariant_shape, {FeBulk::BulkOps::opJac}, fe_component->n_dofs());
663  break;
664  }
665  case FEVectorPiola:
666  {
667  ASSERT_PERMANENT(false).error("Shape vector for FEVectorPiola is not implemented yet!\n"); // temporary assert
668  auto lambda_piola_shape = [shape_values, vector_shape_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
669  side_reinit::ptop_vector_piola_shape(operations, el_table, shape_values, vector_shape_op_idx);
670  };
671  patch_point_vals_.make_fe_op({3}, lambda_piola_shape, {FeBulk::BulkOps::opJac, FeBulk::BulkOps::opJacDet}, fe_component->n_dofs());
672  break;
673  }
674  default:
675  ASSERT(false).error("Type of FiniteElement of grad_vector_shape accessor must be FEVector, FEVectorPiola or FEVectorContravariant!\n");
676  }
677  return FeQ<Vector>(patch_point_vals_, vector_shape_op_idx, fe_component->n_dofs());
678  }
679 
680  /// Same as BulkValues::grad_scalar_shape but register at side quadrature points.
681  inline FeQ<Vector> grad_scalar_shape(uint component_idx=0)
682  {
683  auto fe_component = this->fe_comp(fe_, component_idx);
684  ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of grad_scalar_shape accessor must be FEScalar!\n");
685 
686  // use lambda reinit function
687  auto ref_shape_grads = this->ref_shape_gradients_side(patch_point_vals_.get_quadrature(), fe_component);
688  uint scalar_shape_grads_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
689  auto lambda_scalar_shape_grad = [ref_shape_grads, scalar_shape_grads_op_idx](std::vector<ElOp<3>> &operations, IntTableArena &el_table) {
690  side_reinit::ptop_scalar_shape_grads<dim>(operations, el_table, ref_shape_grads, scalar_shape_grads_op_idx);
691  };
692  patch_point_vals_.make_fe_op({3}, lambda_scalar_shape_grad, {FeSide::SideOps::opElInvJac}, fe_component->n_dofs());
693  patch_point_vals_.set_fe_op(FEOps::opGradScalarShape, scalar_shape_grads_op_idx);
694 
695  return FeQ<Vector>(patch_point_vals_, scalar_shape_grads_op_idx, fe_component->n_dofs());
696  }
697 
698  /**
699  * @brief Return the value of the @p function_no-th gradient vector shape function
700  * at the @p p bulk quadrature point.
701  *
702  * @param component_idx Number of the shape function.
703  */
704  inline FeQ<Tensor> grad_vector_shape(uint component_idx=0)
705  {
706  auto fe_component = this->fe_comp(fe_, component_idx);
707 
708  // use lambda reinit function
709  auto ref_shape_grads = this->ref_shape_gradients_side(patch_point_vals_.get_quadrature(), fe_component);
710  uint vector_shape_grads_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
711 
712  switch (fe_component->fe_type()) {
713  case FEVector:
714  {
715  auto lambda_vector_shape_grad = [ref_shape_grads, vector_shape_grads_op_idx](std::vector<ElOp<3>> &operations, IntTableArena &el_table) {
716  side_reinit::ptop_vector_shape_grads<dim>(operations, el_table, ref_shape_grads, vector_shape_grads_op_idx);
717  };
719  lambda_vector_shape_grad,
721  fe_component->n_dofs());
722  break;
723  }
725  {
726  ASSERT_PERMANENT(false).error("Grad vector for FEVectorContravariant is not implemented yet!\n"); // temporary assert
727  auto lambda_contravariant_shape_grad = [ref_shape_grads, vector_shape_grads_op_idx](std::vector<ElOp<3>> &operations, IntTableArena &el_table) {
728  side_reinit::ptop_vector_contravariant_shape_grads<dim>(operations, el_table, ref_shape_grads, vector_shape_grads_op_idx);
729  };
731  lambda_contravariant_shape_grad,
733  fe_component->n_dofs());
734  break;
735  }
736  case FEVectorPiola:
737  {
738  ASSERT_PERMANENT(false).error("Grad vector for FEVectorPiola is not implemented yet!\n"); // temporary assert
739  auto lambda_piola_shape_grad = [ref_shape_grads, vector_shape_grads_op_idx](std::vector<ElOp<3>> &operations, IntTableArena &el_table) {
740  side_reinit::ptop_vector_piola_shape_grads<dim>(operations, el_table, ref_shape_grads, vector_shape_grads_op_idx);
741  };
743  lambda_piola_shape_grad,
745  fe_component->n_dofs());
746  break;
747  }
748  default:
749  ASSERT(false).error("Type of FiniteElement of grad_vector_shape accessor must be FEVector, FEVectorPiola or FEVectorContravariant!\n");
750  }
751  patch_point_vals_.set_fe_op(FEOps::opGradVectorShape, vector_shape_grads_op_idx);
752 
753  return FeQ<Tensor>(patch_point_vals_, vector_shape_grads_op_idx, fe_component->n_dofs());
754  }
755 
756  /**
757  * @brief Return the value of the @p function_no-th vector symmetric gradient
758  * at the @p p side quadrature point.
759  *
760  * @param component_idx Number of the shape function.
761  */
762  inline FeQ<Tensor> vector_sym_grad(uint component_idx=0)
763  {
764  auto fe_component = this->fe_comp(fe_, component_idx);
765  //ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of grad_scalar_shape accessor must be FEScalar!\n");
766 
767  // use lambda reinit function
768  uint vector_sym_grad_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
770  auto lambda_vector_sym_grad = [vector_sym_grad_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
771  common_reinit::ptop_vector_sym_grad(operations, vector_sym_grad_op_idx);
772  };
773  patch_point_vals_.make_fe_op({3,3}, lambda_vector_sym_grad, {grad_vector_op_idx}, fe_component->n_dofs());
774  patch_point_vals_.set_fe_op(FEOps::opVectorSymGrad, vector_sym_grad_op_idx);
775 
776  return FeQ<Tensor>(patch_point_vals_, vector_sym_grad_op_idx, fe_component->n_dofs());
777  }
778 
779  /**
780  * @brief Return the value of the @p function_no-th vector divergence at
781  * the @p p side quadrature point.
782  *
783  * @param component_idx Number of the shape function.
784  */
785  inline FeQ<Scalar> vector_divergence(uint component_idx=0)
786  {
787  auto fe_component = this->fe_comp(fe_, component_idx);
788  //ASSERT_EQ(fe_component->fe_type(), FEType::FEScalar).error("Type of FiniteElement of grad_scalar_shape accessor must be FEScalar!\n");
789 
790  // use lambda reinit function
791  uint vector_divergence_op_idx = patch_point_vals_.operations_.size(); // index in operations_ vector
793  auto lambda_vector_divergence = [vector_divergence_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
794  common_reinit::ptop_vector_divergence(operations, vector_divergence_op_idx);
795  };
796  patch_point_vals_.make_fe_op({1}, lambda_vector_divergence, {grad_vector_op_idx}, fe_component->n_dofs());
797  patch_point_vals_.set_fe_op(FEOps::opVectorDivergence, vector_divergence_op_idx);
798 
799  return FeQ<Scalar>(patch_point_vals_, vector_divergence_op_idx, fe_component->n_dofs());
800  }
801 
802 private:
804  std::shared_ptr< FiniteElement<dim> > fe_;
805 };
806 
807 
808 template<unsigned int dim>
809 class JoinValues : public BaseValues<dim>
810 {
811 public:
812  /// Constructor
813  JoinValues(PatchPointValues<3> *patch_point_vals_bulk, PatchPointValues<3> *patch_point_vals_side, MixedPtr<FiniteElement> fe)
814  : BaseValues<dim>(), patch_point_vals_bulk_(patch_point_vals_bulk), patch_point_vals_side_(patch_point_vals_side) {
815  ASSERT_EQ(patch_point_vals_bulk->dim(), dim-1);
816  ASSERT_EQ(patch_point_vals_side->dim(), dim);
817  fe_high_dim_ = fe[Dim<dim>{}];
818  fe_low_dim_ = fe[Dim<dim-1>{}];
819  }
820 
822  {
823  // element of lower dim (bulk points)
824  auto fe_component_low = this->fe_comp(fe_low_dim_, component_idx);
825  ASSERT_EQ(fe_component_low->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
826  // use lambda reinit function
827  std::vector< std::vector<double> > shape_values_bulk( fe_component_low->n_dofs(), vector<double>(patch_point_vals_bulk_->get_quadrature()->size()) );
828  auto ref_shape_vals_bulk = this->ref_shape_values_bulk(patch_point_vals_bulk_->get_quadrature(), fe_component_low);
829  for (unsigned int i = 0; i < patch_point_vals_bulk_->get_quadrature()->size(); i++)
830  for (unsigned int j = 0; j < fe_component_low->n_dofs(); j++) {
831  shape_values_bulk[j][i] = ref_shape_vals_bulk[i][j][0];
832  }
833  uint scalar_shape_op_idx_bulk = patch_point_vals_bulk_->operations_.size(); // index in operations_ vector
834  auto lambda_scalar_shape_bulk = [shape_values_bulk, scalar_shape_op_idx_bulk](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
835  bulk_reinit::ptop_scalar_shape(operations, shape_values_bulk, scalar_shape_op_idx_bulk);
836  };
837  patch_point_vals_bulk_->make_fe_op({1}, lambda_scalar_shape_bulk, {}, fe_component_low->n_dofs());
838  uint op_idx_bulk = patch_point_vals_bulk_->operations_.size()-1;
839 
840  // element of higher dim (side points)
841  auto fe_component_high = this->fe_comp(fe_high_dim_, component_idx);
842  ASSERT_EQ(fe_component_high->fe_type(), FEType::FEScalar).error("Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
843  // use lambda reinit function
844  std::vector< std::vector< std::vector<double> > > shape_values_side(
845  dim+1,
847  );
848  auto ref_shape_vals_side = this->ref_shape_values_side(patch_point_vals_side_->get_quadrature(), fe_component_high);
849  for (unsigned int s=0; s<dim+1; ++s)
850  for (unsigned int i = 0; i < patch_point_vals_side_->get_quadrature()->size(); i++)
851  for (unsigned int j = 0; j < fe_component_high->n_dofs(); j++) {
852  shape_values_side[s][i][j] = ref_shape_vals_side[s][i][j][0];
853  }
854  uint scalar_shape_op_idx_side = patch_point_vals_side_->operations_.size(); // index in operations_ vector
855  auto lambda_scalar_shape_side = [shape_values_side, scalar_shape_op_idx_side](std::vector<ElOp<3>> &operations, IntTableArena &el_table) {
856  side_reinit::ptop_scalar_shape(operations, el_table, shape_values_side, scalar_shape_op_idx_side);
857  };
858  patch_point_vals_side_->make_fe_op({1}, lambda_scalar_shape_side, {}, fe_component_high->n_dofs());
859  uint op_idx_side = patch_point_vals_side_->operations_.size()-1;
860 
861  auto bgn_it = make_iter<JoinShapeAccessor<Scalar>>( JoinShapeAccessor<Scalar>(patch_point_vals_bulk_, patch_point_vals_side_,
862  fe_component_low->n_dofs(), fe_component_high->n_dofs(), op_idx_bulk, op_idx_side, 0) );
863  unsigned int end_idx = fe_component_low->n_dofs() + fe_component_high->n_dofs();
864  auto end_it = make_iter<JoinShapeAccessor<Scalar>>( JoinShapeAccessor<Scalar>(patch_point_vals_bulk_, patch_point_vals_side_,
865  fe_component_low->n_dofs(), fe_component_high->n_dofs(), op_idx_bulk, op_idx_side, end_idx) );
866  return Range<JoinShapeAccessor<Scalar>>(bgn_it, end_it);
867  }
868 
870  {
871  // element of lower dim (bulk points)
872  auto fe_component_low = this->fe_comp(fe_low_dim_, component_idx);
873  std::vector< std::vector<arma::vec3> > shape_values_bulk( fe_component_low->n_dofs(), vector<arma::vec3>(patch_point_vals_bulk_->get_quadrature()->size()) );
874  auto ref_shape_vals_bulk = this->ref_shape_values_bulk(patch_point_vals_bulk_->get_quadrature(), fe_component_low);
875  for (unsigned int i = 0; i < patch_point_vals_bulk_->get_quadrature()->size(); i++)
876  for (unsigned int j = 0; j < fe_component_low->n_dofs(); j++) {
877  shape_values_bulk[j][i] = ref_shape_vals_bulk[i][j];
878  }
879  uint op_idx_bulk = patch_point_vals_bulk_->operations_.size(); // index in operations_ vector
880 
881  // element of higher dim (side points)
882  auto fe_component_high = this->fe_comp(fe_high_dim_, component_idx);
884  dim+1,
886  );
887  auto ref_shape_vals_side = this->ref_shape_values_side(patch_point_vals_side_->get_quadrature(), fe_component_high);
888  for (unsigned int s=0; s<dim+1; ++s)
889  for (unsigned int i = 0; i < patch_point_vals_side_->get_quadrature()->size(); i++)
890  for (unsigned int j = 0; j < fe_component_high->n_dofs(); j++) {
891  shape_values_side[s][i][j] = ref_shape_vals_side[s][i][j];
892  }
893  uint op_idx_side = patch_point_vals_side_->operations_.size(); // index in operations_ vector
894 
895  ASSERT_EQ(fe_component_high->fe_type(), fe_component_low->fe_type()).error("Type of FiniteElement of low and high element must be same!\n");
896  switch (fe_component_low->fe_type()) {
897  case FEVector:
898  {
899  // use lambda reinit function (lower dim)
900  auto lambda_vector_shape_bulk = [shape_values_bulk, op_idx_bulk](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
901  bulk_reinit::ptop_vector_shape(operations, shape_values_bulk, op_idx_bulk);
902  };
903  patch_point_vals_bulk_->make_fe_op({3}, lambda_vector_shape_bulk, {}, fe_component_low->n_dofs());
904 
905  // use lambda reinit function (higher dim)
906  auto lambda_vector_shape_side = [shape_values_side, op_idx_side](std::vector<ElOp<3>> &operations, IntTableArena &el_table) {
907  side_reinit::ptop_vector_shape(operations, el_table, shape_values_side, op_idx_side);
908  };
909  patch_point_vals_side_->make_fe_op({3}, lambda_vector_shape_side, {}, fe_component_high->n_dofs());
910  break;
911  }
913  {
914  ASSERT_PERMANENT(false).error("Shape vector for FEVectorContravariant is not implemented yet!\n"); // temporary assert
915  //auto lambda_contravariant_shape = [shape_values, vector_shape_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
916  // bulk_reinit::ptop_vector_contravariant_shape(operations, shape_values, vector_shape_op_idx);
917  // };
918  //patch_point_vals_.make_fe_op({3}, lambda_contravariant_shape, {FeBulk::BulkOps::opJac}, fe_component->n_dofs());
919  break;
920  }
921  case FEVectorPiola:
922  {
923  ASSERT_PERMANENT(false).error("Shape vector for FEVectorPiola is not implemented yet!\n"); // temporary assert
924  //auto lambda_piola_shape = [shape_values, vector_shape_op_idx](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
925  // bulk_reinit::ptop_vector_piola_shape(operations, shape_values, vector_shape_op_idx);
926  // };
927  //patch_point_vals_.make_fe_op({3}, lambda_piola_shape, {FeBulk::BulkOps::opJac, FeBulk::BulkOps::opJacDet}, fe_component->n_dofs());
928  break;
929  }
930  default:
931  ASSERT(false).error("Type of FiniteElement of grad_vector_shape accessor must be FEVector, FEVectorPiola or FEVectorContravariant!\n");
932  }
933 
934  auto bgn_it = make_iter<JoinShapeAccessor<Vector>>( JoinShapeAccessor<Vector>(patch_point_vals_bulk_, patch_point_vals_side_,
935  fe_component_low->n_dofs(), fe_component_high->n_dofs(), op_idx_bulk, op_idx_side, 0) );
936  unsigned int end_idx = fe_component_low->n_dofs() + fe_component_high->n_dofs();
937  auto end_it = make_iter<JoinShapeAccessor<Vector>>( JoinShapeAccessor<Vector>(patch_point_vals_bulk_, patch_point_vals_side_,
938  fe_component_low->n_dofs(), fe_component_high->n_dofs(), op_idx_bulk, op_idx_side, end_idx) );
939  return Range<JoinShapeAccessor<Vector>>(bgn_it, end_it);
940  }
941 
943  {
944  // element of lower dim (bulk points)
945  auto fe_component_low = this->fe_comp(fe_low_dim_, component_idx);
946  auto ref_shape_grads_bulk = this->ref_shape_gradients_bulk(patch_point_vals_bulk_->get_quadrature(), fe_component_low);
947  uint op_idx_bulk = patch_point_vals_bulk_->operations_.size(); // index in operations_ vector
948 
949  // element of higher dim (side points)
950  auto fe_component_high = this->fe_comp(fe_high_dim_, component_idx);
951  auto ref_shape_grads_side = this->ref_shape_gradients_side(patch_point_vals_side_->get_quadrature(), fe_component_high);
952  uint op_idx_side = patch_point_vals_side_->operations_.size(); // index in operations_ vector
953 
954  ASSERT_EQ(fe_component_high->fe_type(), fe_component_low->fe_type()).error("Type of FiniteElement of low and high element must be same!\n");
955  switch (fe_component_low->fe_type()) {
956  case FEVector:
957  {
958  // use lambda reinit function (lower dim)
959  auto lambda_vector_grad_bulk = [ref_shape_grads_bulk, op_idx_bulk](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
960  bulk_reinit::ptop_vector_shape_grads<dim-1>(operations, ref_shape_grads_bulk, op_idx_bulk);
961  };
963  lambda_vector_grad_bulk,
965  fe_component_low->n_dofs());
966 
967  // use lambda reinit function (higher dim)
968  auto lambda_vector_grad_side = [ref_shape_grads_side, op_idx_side](std::vector<ElOp<3>> &operations, IntTableArena &el_table) {
969  side_reinit::ptop_vector_shape_grads<dim>(operations, el_table, ref_shape_grads_side, op_idx_side);
970  };
972  lambda_vector_grad_side,
974  fe_component_high->n_dofs());
975  break;
976  }
978  {
979  ASSERT_PERMANENT(false).error("Shape vector for FEVectorContravariant is not implemented yet!\n"); // temporary assert
980 // // use lambda reinit function (lower dim)
981 // auto lambda_contravariant_grad_bulk = [ref_shape_grads_bulk, op_idx_bulk](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
982 // bulk_reinit::ptop_vector_contravariant_shape_grads<dim-1>(operations, ref_shape_grads_bulk, op_idx_bulk);
983 // };
984 // patch_point_vals_bulk_->make_fe_op({3, 3},
985 // lambda_contravariant_grad_bulk,
986 // {FeBulk::BulkOps::opInvJac, FeBulk::BulkOps::opJac},
987 // fe_component_low->n_dofs());
988 //
989 // // use lambda reinit function (higher dim)
990 // auto lambda_contravariant_grad_side = [ref_shape_grads_side, op_idx_side](std::vector<ElOp<3>> &operations, IntTableArena &el_table) {
991 // side_reinit::ptop_vector_contravariant_shape_grads<dim>(operations, el_table, ref_shape_grads_side, op_idx_side);
992 // };
993 // patch_point_vals_side_->make_fe_op({3, 3},
994 // lambda_contravariant_grad_side,
995 // {FeSide::SideOps::opElInvJac, FeSide::SideOps::opElJac},
996 // fe_component_high->n_dofs());
997  break;
998  }
999  case FEVectorPiola:
1000  {
1001  ASSERT_PERMANENT(false).error("Shape vector for FEVectorPiola is not implemented yet!\n"); // temporary assert
1002 // // use lambda reinit function (lower dim)
1003 // auto lambda_piola_grad_bulk = [ref_shape_grads_bulk, op_idx_bulk](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
1004 // bulk_reinit::ptop_vector_piola_shape_grads<dim-1>(operations, ref_shape_grads_bulk, op_idx_bulk);
1005 // };
1006 // patch_point_vals_bulk_->make_fe_op({3, 3},
1007 // lambda_piola_grad_bulk,
1008 // {FeBulk::BulkOps::opInvJac, FeBulk::BulkOps::opJac, FeBulk::BulkOps::opJacDet},
1009 // fe_component_low->n_dofs());
1010 //
1011 // // use lambda reinit function (higher dim)
1012 // auto lambda_piola_grad_side = [ref_shape_grads_side, op_idx_side](std::vector<ElOp<3>> &operations, IntTableArena &el_table) {
1013 // side_reinit::ptop_vector_piola_shape_grads<dim>(operations, el_table, ref_shape_grads_side, op_idx_side);
1014 // };
1015 // patch_point_vals_side_->make_fe_op({3, 3},
1016 // lambda_piola_grad_side,
1017 // {FeSide::SideOps::opElInvJac, FeSide::SideOps::opElJac, FeSide::SideOps::opSideJacDet},
1018 // fe_component_high->n_dofs());
1019  break;
1020  }
1021  default:
1022  ASSERT(false).error("Type of FiniteElement of grad_vector_shape accessor must be FEVector, FEVectorPiola or FEVectorContravariant!\n");
1023  }
1024 
1025  auto bgn_it = make_iter<JoinShapeAccessor<Tensor>>( JoinShapeAccessor<Tensor>(patch_point_vals_bulk_, patch_point_vals_side_,
1026  fe_component_low->n_dofs(), fe_component_high->n_dofs(), op_idx_bulk, op_idx_side, 0) );
1027  unsigned int end_idx = fe_component_low->n_dofs() + fe_component_high->n_dofs();
1028  auto end_it = make_iter<JoinShapeAccessor<Tensor>>( JoinShapeAccessor<Tensor>(patch_point_vals_bulk_, patch_point_vals_side_,
1029  fe_component_low->n_dofs(), fe_component_high->n_dofs(), op_idx_bulk, op_idx_side, end_idx) );
1030  return Range<JoinShapeAccessor<Tensor>>(bgn_it, end_it);
1031  }
1032 
1033 private:
1036  std::shared_ptr< FiniteElement<dim> > fe_high_dim_;
1037  std::shared_ptr< FiniteElement<dim-1> > fe_low_dim_;
1038 };
1039 
1040 /// Template specialization of dim = 1
1041 template <>
1042 class JoinValues<1> : public BaseValues<1>
1043 {
1044 public:
1045  /// Constructor
1047  : BaseValues<1>() {}
1048 
1050  {
1052  make_iter<JoinShapeAccessor<Scalar>>(JoinShapeAccessor<Scalar>()),
1054  }
1055 
1057  {
1059  make_iter<JoinShapeAccessor<Vector>>(JoinShapeAccessor<Vector>()),
1061  }
1062 
1064  {
1066  make_iter<JoinShapeAccessor<Tensor>>(JoinShapeAccessor<Tensor>()),
1068  }
1069 };
1070 
1071 
1072 template<unsigned int spacedim = 3>
1074 public:
1075  /// Struct for pre-computing number of elements, sides, bulk points and side points on each dimension.
1076  struct TableSizes {
1077  public:
1078  /// Constructor
1082  }
1083 
1084  /// Set all values to zero
1085  void reset() {
1086  std::fill(elem_sizes_[0].begin(), elem_sizes_[0].end(), 0);
1087  std::fill(elem_sizes_[1].begin(), elem_sizes_[1].end(), 0);
1088  std::fill(point_sizes_[0].begin(), point_sizes_[0].end(), 0);
1089  std::fill(point_sizes_[1].begin(), point_sizes_[1].end(), 0);
1090  }
1091 
1092  /// Copy values of other TableSizes instance
1093  void copy(const TableSizes &other) {
1094  elem_sizes_[0] = other.elem_sizes_[0];
1095  elem_sizes_[1] = other.elem_sizes_[1];
1096  point_sizes_[0] = other.point_sizes_[0];
1097  point_sizes_[1] = other.point_sizes_[1];
1098  }
1099 
1100  /**
1101  * Holds number of elements and sides on each dimension
1102  * Format:
1103  * { {n_elements_1D, n_elements_2D, n_elements_3D },
1104  * {n_sides_1D, n_sides_2D, n_sides_3D } }
1105  */
1107 
1108  /**
1109  * Holds number of bulk and side points on each dimension
1110  * Format:
1111  * { {n_bulk_points_1D, n_bulk_points_2D, n_bulk_points_3D },
1112  * {n_side_points_1D, n_side_points_2D, n_side_points_3D } }
1113  */
1115  };
1116 
1118  : asm_arena_(1024 * 1024, 256),
1119  patch_arena_(nullptr),
1123  patch_point_vals_side_{ {FeSide::PatchPointValues(1, 0, asm_arena_),
1124  FeSide::PatchPointValues(2, 0, asm_arena_),
1125  FeSide::PatchPointValues(3, 0, asm_arena_)} }
1126  {
1127  used_quads_[0] = false; used_quads_[1] = false;
1128  }
1129 
1130  PatchFEValues(unsigned int quad_order, MixedPtr<FiniteElement> fe)
1131  : asm_arena_(1024 * 1024, 256),
1132  patch_arena_(nullptr),
1133  patch_point_vals_bulk_{ {FeBulk::PatchPointValues(1, quad_order, asm_arena_),
1134  FeBulk::PatchPointValues(2, quad_order, asm_arena_),
1135  FeBulk::PatchPointValues(3, quad_order, asm_arena_)} },
1136  patch_point_vals_side_{ {FeSide::PatchPointValues(1, quad_order, asm_arena_),
1137  FeSide::PatchPointValues(2, quad_order, asm_arena_),
1138  FeSide::PatchPointValues(3, quad_order, asm_arena_)} },
1139  fe_(fe)
1140  {
1141  used_quads_[0] = false; used_quads_[1] = false;
1142  }
1143 
1144 
1145  /// Destructor
1147  {
1148  if (patch_arena_!=nullptr)
1149  delete patch_arena_;
1150  }
1151 
1152  /// Return bulk or side quadrature of given dimension
1153  Quadrature *get_quadrature(uint dim, bool is_bulk) const {
1154  if (is_bulk) return patch_point_vals_bulk_[dim-1].get_quadrature();
1155  else return patch_point_vals_side_[dim-1].get_quadrature();
1156  }
1157 
1158  /**
1159  * @brief Initialize structures and calculates cell-independent data.
1160  *
1161  * @param _quadrature The quadrature rule for the cell associated
1162  * to given finite element or for the cell side.
1163  * @param _flags The update flags.
1164  */
1165  template<unsigned int DIM>
1166  void initialize(Quadrature &_quadrature)
1167  {
1168  if ( _quadrature.dim() == DIM ) {
1169  used_quads_[0] = true;
1170  patch_point_vals_bulk_[DIM-1].initialize(); // bulk
1171  } else {
1172  used_quads_[1] = true;
1173  patch_point_vals_side_[DIM-1].initialize(); // side
1174  }
1175  }
1176 
1177  /// Finalize initialization, creates child (patch) arena and passes it to PatchPointValue objects
1178  void init_finalize() {
1179  patch_arena_ = asm_arena_.get_child_arena();
1180  for (unsigned int i=0; i<3; ++i) {
1181  if (used_quads_[0]) patch_point_vals_bulk_[i].init_finalize(patch_arena_);
1182  if (used_quads_[1]) patch_point_vals_side_[i].init_finalize(patch_arena_);
1183  }
1184  }
1185 
1186  /// Reset PatchpointValues structures
1187  void reset()
1188  {
1189  for (unsigned int i=0; i<3; ++i) {
1190  if (used_quads_[0]) patch_point_vals_bulk_[i].reset();
1191  if (used_quads_[1]) patch_point_vals_side_[i].reset();
1192  }
1193  patch_arena_->reset();
1194  }
1195 
1196  /// Reinit data.
1198  {
1199  for (unsigned int i=0; i<3; ++i) {
1200  if (used_quads_[0]) patch_point_vals_bulk_[i].reinit_patch();
1201  if (used_quads_[1]) patch_point_vals_side_[i].reinit_patch();
1202  }
1203  }
1204 
1205  /**
1206  * @brief Returns the number of shape functions.
1207  */
1208  template<unsigned int dim>
1209  inline unsigned int n_dofs() const {
1210  ASSERT((dim>=0) && (dim<=3))(dim).error("Dimension must be 0, 1, 2 or 3.");
1211  return fe_[Dim<dim>{}]->n_dofs();
1212  }
1213 
1214  /// Return BulkValue object of dimension given by template parameter
1215  template<unsigned int dim>
1217  ASSERT((dim>0) && (dim<=3))(dim).error("Dimension must be 1, 2 or 3.");
1218  return BulkValues<dim>(patch_point_vals_bulk_[dim-1], fe_);
1219  }
1220 
1221  /// Return SideValue object of dimension given by template parameter
1222  template<unsigned int dim>
1224  ASSERT((dim>0) && (dim<=3))(dim).error("Dimension must be 1, 2 or 3.");
1225  return SideValues<dim>(patch_point_vals_side_[dim-1], fe_);
1226  }
1227 
1228  /// Return JoinValue object of dimension given by template parameter
1229  template<unsigned int dim>
1231  //ASSERT((dim>1) && (dim<=3))(dim).error("Dimension must be 2 or 3.");
1232  return JoinValues<dim>(&patch_point_vals_bulk_[dim-2], &patch_point_vals_side_[dim-1], fe_);
1233  }
1234 
1235  /** Following methods are used during update of patch. **/
1236 
1237  /// Resize tables of patch_point_vals_
1238  void resize_tables(TableSizes table_sizes) {
1239  for (uint i=0; i<3; ++i) {
1240  if (used_quads_[0]) patch_point_vals_bulk_[i].resize_tables(table_sizes.elem_sizes_[0][i], table_sizes.point_sizes_[0][i]);
1241  if (used_quads_[1]) patch_point_vals_side_[i].resize_tables(table_sizes.elem_sizes_[1][i], table_sizes.point_sizes_[1][i]);
1242  }
1243  }
1244 
1245  /// Register element to patch_point_vals_ table by dimension of element
1246  uint register_element(DHCellAccessor cell, uint element_patch_idx) {
1247  arma::mat coords;
1248  switch (cell.dim()) {
1249  case 1:
1250  coords = MappingP1<1,spacedim>::element_map(cell.elm());
1251  return patch_point_vals_bulk_[0].register_element(coords, element_patch_idx);
1252  break;
1253  case 2:
1254  coords = MappingP1<2,spacedim>::element_map(cell.elm());
1255  return patch_point_vals_bulk_[1].register_element(coords, element_patch_idx);
1256  break;
1257  case 3:
1258  coords = MappingP1<3,spacedim>::element_map(cell.elm());
1259  return patch_point_vals_bulk_[2].register_element(coords, element_patch_idx);
1260  break;
1261  default:
1262  ASSERT(false);
1263  return 0;
1264  break;
1265  }
1266  }
1267 
1268  /// Register side to patch_point_vals_ table by dimension of side
1270  arma::mat side_coords(spacedim, cell_side.dim());
1271  for (unsigned int n=0; n<cell_side.dim(); n++)
1272  for (unsigned int c=0; c<spacedim; c++)
1273  side_coords(c,n) = (*cell_side.side().node(n))[c];
1274 
1275  arma::mat elm_coords;
1276  DHCellAccessor cell = cell_side.cell();
1277  switch (cell.dim()) {
1278  case 1:
1279  elm_coords = MappingP1<1,spacedim>::element_map(cell.elm());
1280  return patch_point_vals_side_[0].register_side(elm_coords, side_coords, cell_side.side_idx());
1281  break;
1282  case 2:
1283  elm_coords = MappingP1<2,spacedim>::element_map(cell.elm());
1284  return patch_point_vals_side_[1].register_side(elm_coords, side_coords, cell_side.side_idx());
1285  break;
1286  case 3:
1287  elm_coords = MappingP1<3,spacedim>::element_map(cell.elm());
1288  return patch_point_vals_side_[2].register_side(elm_coords, side_coords, cell_side.side_idx());
1289  break;
1290  default:
1291  ASSERT(false);
1292  return 0;
1293  break;
1294  }
1295  }
1296 
1297  /// Register bulk point to patch_point_vals_ table by dimension of element
1298  uint register_bulk_point(DHCellAccessor cell, uint elem_table_row, uint value_patch_idx, uint i_point_on_elem) {
1299  return patch_point_vals_bulk_[cell.dim()-1].register_bulk_point(elem_table_row, value_patch_idx, cell.elm_idx(), i_point_on_elem);
1300  }
1301 
1302  /// Register side point to patch_point_vals_ table by dimension of side
1303  uint register_side_point(DHCellSide cell_side, uint elem_table_row, uint value_patch_idx, uint i_point_on_side) {
1304  return patch_point_vals_side_[cell_side.dim()-1].register_side_point(elem_table_row, value_patch_idx, cell_side.elem_idx(),
1305  cell_side.side_idx(), i_point_on_side);
1306  }
1307 
1308  /// Temporary development method
1309  void print_data_tables(ostream& stream, bool points, bool ints, bool only_bulk=true) const {
1310  stream << endl << "Table of patch FE data:" << endl;
1311  for (uint i=0; i<3; ++i) {
1312  stream << std::setfill('-') << setw(100) << "" << endl;
1313  stream << "Bulk, dimension " << (i+1) << endl;
1314  patch_point_vals_bulk_[i].print_data_tables(stream, points, ints);
1315  }
1316  if (!only_bulk)
1317  for (uint i=0; i<3; ++i) {
1318  stream << std::setfill('-') << setw(100) << "" << endl;
1319  stream << "Side, dimension " << (i+1) << endl;
1320  patch_point_vals_side_[i].print_data_tables(stream, points, ints);
1321  }
1322  stream << std::setfill('=') << setw(100) << "" << endl;
1323  }
1324 
1325  /// Temporary development method
1326  void print_operations(ostream& stream) const {
1327  stream << endl << "Table of patch FE operations:" << endl;
1328  for (uint i=0; i<3; ++i) {
1329  stream << std::setfill('-') << setw(100) << "" << endl;
1330  stream << "Bulk, dimension " << (i+1) << endl;
1331  patch_point_vals_bulk_[i].print_operations(stream, 0);
1332  }
1333  for (uint i=0; i<3; ++i) {
1334  stream << std::setfill('-') << setw(100) << "" << endl;
1335  stream << "Side, dimension " << (i+1) << endl;
1336  patch_point_vals_side_[i].print_operations(stream, 1);
1337  }
1338  stream << std::setfill('=') << setw(100) << "" << endl;
1339  }
1340 
1341 private:
1344  std::array<FeBulk::PatchPointValues<spacedim>, 3> patch_point_vals_bulk_; ///< Sub objects of bulk data of dimensions 1,2,3
1345  std::array<FeSide::PatchPointValues<spacedim>, 3> patch_point_vals_side_; ///< Sub objects of side data of dimensions 1,2,3
1346 
1347  MixedPtr<FiniteElement> fe_; ///< Mixed of shared pointers of FiniteElement object
1348  bool used_quads_[2]; ///< Pair of flags signs holds info if bulk and side quadratures are used
1349 
1350  template <class ValueType>
1351  friend class ElQ;
1352  template <class ValueType>
1353  friend class FeQ;
1354  template <class ValueType>
1355  friend class JoinShapeAccessor;
1356 };
1357 
1358 
1359 template <class ValueType>
1360 ValueType ElQ<ValueType>::operator()(const BulkPoint &point) const {
1361  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1362  return patch_point_vals_.scalar_value(op_idx_, value_cache_idx);
1363 }
1364 
1365 template <>
1366 inline Vector ElQ<Vector>::operator()(const BulkPoint &point) const {
1367  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1368  return patch_point_vals_.vector_value(op_idx_, value_cache_idx);
1369 }
1370 
1371 template <>
1372 inline Tensor ElQ<Tensor>::operator()(const BulkPoint &point) const {
1373  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1374  return patch_point_vals_.tensor_value(op_idx_, value_cache_idx);
1375 }
1376 
1377 template <class ValueType>
1378 ValueType ElQ<ValueType>::operator()(const SidePoint &point) const {
1379  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1380  return patch_point_vals_.scalar_value(op_idx_, value_cache_idx);
1381 }
1382 
1383 template <>
1384 inline Vector ElQ<Vector>::operator()(const SidePoint &point) const {
1385  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1386  return patch_point_vals_.vector_value(op_idx_, value_cache_idx);
1387 }
1388 
1389 template <>
1390 inline Tensor ElQ<Tensor>::operator()(const SidePoint &point) const {
1391  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1392  return patch_point_vals_.tensor_value(op_idx_, value_cache_idx);
1393 }
1394 
1395 template <class ValueType>
1396 ValueType FeQ<ValueType>::operator()(unsigned int shape_idx, const BulkPoint &point) const {
1397  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1398  return patch_point_vals_.scalar_value(op_idx_, value_cache_idx, shape_idx);
1399 }
1400 
1401 template <>
1402 inline Vector FeQ<Vector>::operator()(unsigned int shape_idx, const BulkPoint &point) const {
1403  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1404  return patch_point_vals_.vector_value(op_idx_, value_cache_idx, shape_idx);
1405 }
1406 
1407 template <>
1408 inline Tensor FeQ<Tensor>::operator()(unsigned int shape_idx, const BulkPoint &point) const {
1409  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1410  return patch_point_vals_.tensor_value(op_idx_, value_cache_idx, shape_idx);
1411 }
1412 
1413 template <class ValueType>
1414 ValueType FeQ<ValueType>::operator()(unsigned int shape_idx, const SidePoint &point) const {
1415  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1416  return patch_point_vals_.scalar_value(op_idx_, value_cache_idx, shape_idx);
1417 }
1418 
1419 template <>
1420 inline Vector FeQ<Vector>::operator()(unsigned int shape_idx, const SidePoint &point) const {
1421  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1422  return patch_point_vals_.vector_value(op_idx_, value_cache_idx, shape_idx);
1423 }
1424 
1425 template <>
1426 inline Tensor FeQ<Tensor>::operator()(unsigned int shape_idx, const SidePoint &point) const {
1427  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1428  return patch_point_vals_.tensor_value(op_idx_, value_cache_idx, shape_idx);
1429 }
1430 
1431 
1432 template <class ValueType>
1434  if (this->is_high_dim()) {
1435  return 0.0;
1436  } else {
1437  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1438  return patch_point_vals_bulk_->scalar_value(op_idx_bulk_, value_cache_idx, this->local_idx());
1439  }
1440 }
1441 
1442 template <>
1444  if (this->is_high_dim()) {
1445  Vector vect; vect.zeros();
1446  return vect;
1447  } else {
1448  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1449  return patch_point_vals_bulk_->vector_value(op_idx_bulk_, value_cache_idx, this->local_idx());
1450  }
1451 }
1452 
1453 template <>
1455  if (this->is_high_dim()) {
1456  Tensor tens; tens.zeros();
1457  return tens;
1458  } else {
1459  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1460  return patch_point_vals_bulk_->tensor_value(op_idx_bulk_, value_cache_idx, this->local_idx());
1461  }
1462 }
1463 
1464 template <class ValueType>
1466  if (this->is_high_dim()) {
1467  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1468  return patch_point_vals_side_->scalar_value(op_idx_side_, value_cache_idx, this->local_idx());
1469  } else {
1470  return 0.0;
1471  }
1472 }
1473 
1474 template <>
1476  if (this->is_high_dim()) {
1477  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1478  return patch_point_vals_side_->vector_value(op_idx_side_, value_cache_idx, this->local_idx());
1479  } else {
1480  Vector vect; vect.zeros();
1481  return vect;
1482  }
1483 }
1484 
1485 template <>
1487  if (this->is_high_dim()) {
1488  unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
1489  return patch_point_vals_side_->tensor_value(op_idx_side_, value_cache_idx, this->local_idx());
1490  } else {
1491  Tensor tens; tens.zeros();
1492  return tens;
1493  }
1494 }
1495 
1496 
1497 #endif /* PATCH_FE_VALUES_HH_ */
#define ASSERT(expr)
Definition: asserts.hh:351
#define ASSERT_PERMANENT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
Definition: asserts.hh:348
#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::mat > > ref_shape_gradients_bulk(Quadrature *q, std::shared_ptr< FiniteElement< FE_dim >> fe)
Precomputed gradients of basis functions at the bulk quadrature points.
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.
std::vector< std::vector< std::vector< arma::mat > > > ref_shape_gradients_side(Quadrature *q, std::shared_ptr< FiniteElement< FE_dim >> fe)
Precomputed gradients of basis functions at the side 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.
FeQ< Tensor > vector_sym_grad(uint component_idx=0)
Return the value of the function_no-th vector symmetric gradient at the p bulk quadrature point.
PatchPointValues< 3 > & patch_point_vals_
ElQ< Vector > coords()
Create bulk accessor of coords entity.
FeQ< Tensor > grad_vector_shape(uint component_idx=0)
Return the value of the function_no-th gradient vector shape function at the p bulk quadrature point.
ElQ< Scalar > determinant()
Create bulk accessor of jac determinant entity.
FeQ< Scalar > vector_divergence(uint component_idx=0)
Return the value of the function_no-th vector divergence at the p bulk quadrature point.
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.
FeQ< Vector > vector_shape(uint component_idx=0)
BulkValues(PatchPointValues< 3 > &patch_point_vals, MixedPtr< FiniteElement > fe)
Constructor.
Cell accessor allow iterate over DOF handler cells.
unsigned int dim() const
Return dimension of element appropriate to cell.
unsigned int elm_idx() const
Return serial idx to element of loc_ele_idx_.
ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
Side accessor allows to iterate over sides of DOF handler cell.
unsigned int elem_idx() const
Side side() const
Return Side of given cell and side_idx.
const DHCellAccessor & cell() const
Return DHCellAccessor appropriate to the side.
unsigned int dim() const
Return dimension of element appropriate to the side.
unsigned int side_idx() const
ElQ(PatchPointValues< 3 > &patch_point_vals, unsigned int op_idx)
Constructor.
ValueType operator()(const SidePoint &point) const
ValueType operator()(const BulkPoint &point) const
unsigned int op_idx_
Index of operation in patch_point_vals_.operations vector.
ElQ()=delete
Forbidden default constructor.
PatchPointValues< 3 > & patch_point_vals_
Reference to PatchPointValues.
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.
ValueType operator()(unsigned int shape_idx, FMT_UNUSED const SidePoint &point) const
ValueType operator()(unsigned int shape_idx, FMT_UNUSED const BulkPoint &point) const
PatchPointValues< 3 > & patch_point_vals_
Reference to PatchPointValues.
FeQ(PatchPointValues< 3 > &patch_point_vals, unsigned int op_idx, unsigned int n_dofs)
unsigned int n_dofs_
Number of DOFs.
unsigned int op_idx_
Index of operation in patch_point_vals_.operations vector.
FeQ()=delete
Forbidden default constructor.
Bulk Side specialization, order of item in operations_ vector corresponds to the SideOps enum.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
unsigned int n_dofs_both() const
unsigned int local_idx() const
Return local index of DOF (on low / high-dim) - should be private method.
bool operator==(const JoinShapeAccessor< ValueType > &other) const
Comparison of accessors.
JoinShapeAccessor(PatchPointValues< 3 > *patch_point_vals_bulk, PatchPointValues< 3 > *patch_point_vals_side, unsigned int n_dofs_bulk, unsigned int n_dofs_side, unsigned int op_idx_bulk, unsigned int op_idx_side, unsigned int join_idx)
bool is_high_dim() const
PatchPointValues< 3 > * patch_point_vals_bulk_
Pointer to bulk PatchPointValues.
unsigned int op_idx_side_
Index of operation in patch_point_vals_side_.operations vector.
unsigned int op_idx_bulk_
Index of operation in patch_point_vals_bulk_.operations vector.
unsigned int n_dofs_high() const
JoinShapeAccessor()
Default constructor.
void inc()
Iterates to next item.
unsigned int join_idx_
Index of processed DOF.
unsigned int n_dofs_high_
Number of DOFs on high-dim element.
unsigned int join_idx() const
Return global index of DOF.
unsigned int n_dofs_low() const
ValueType operator()(const BulkPoint &point) const
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< Vector > > vector_join_shape(FMT_UNUSED uint component_idx=0)
Range< JoinShapeAccessor< Tensor > > gradient_vector_join_shape(FMT_UNUSED uint component_idx=0)
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.
Range< JoinShapeAccessor< Vector > > vector_join_shape(uint component_idx=0)
PatchPointValues< 3 > * patch_point_vals_bulk_
Range< JoinShapeAccessor< Tensor > > gradient_vector_join_shape(uint component_idx=0)
Range< JoinShapeAccessor< Scalar > > scalar_join_shape(uint component_idx=0)
std::shared_ptr< FiniteElement< dim-1 > > fe_low_dim_
PatchPointValues< 3 > * patch_point_vals_side_
static ElementMap element_map(ElementAccessor< 3 > elm)
Definition: mapping_p1.cc:48
uint register_element(DHCellAccessor cell, uint element_patch_idx)
Register element to patch_point_vals_ table by dimension of element.
PatchArena * patch_arena_
uint register_bulk_point(DHCellAccessor cell, uint elem_table_row, uint value_patch_idx, uint i_point_on_elem)
Register bulk point to patch_point_vals_ table by dimension of element.
~PatchFEValues()
Destructor.
std::array< FeBulk::PatchPointValues< spacedim >, 3 > patch_point_vals_bulk_
Sub objects of bulk data of dimensions 1,2,3.
uint register_side(DHCellSide cell_side)
Register side to patch_point_vals_ table by dimension of side.
void init_finalize()
Finalize initialization, creates child (patch) arena and passes it to PatchPointValue objects.
void print_data_tables(ostream &stream, bool points, bool ints, bool only_bulk=true) const
Temporary development method.
JoinValues< dim > join_values()
Return JoinValue object of dimension given by template parameter.
Quadrature * get_quadrature(uint dim, bool is_bulk) const
Return bulk or side quadrature of given dimension.
BulkValues< dim > bulk_values()
Return BulkValue object of dimension given by template parameter.
uint register_side_point(DHCellSide cell_side, uint elem_table_row, uint value_patch_idx, uint i_point_on_side)
Register side point to patch_point_vals_ table by dimension of side.
void initialize(Quadrature &_quadrature)
Initialize structures and calculates cell-independent data.
void print_operations(ostream &stream) const
Temporary development method.
std::array< FeSide::PatchPointValues< spacedim >, 3 > patch_point_vals_side_
Sub objects of side data of dimensions 1,2,3.
MixedPtr< FiniteElement > fe_
Mixed of shared pointers of FiniteElement object.
SideValues< dim > side_values()
Return SideValue object of dimension given by template parameter.
AssemblyArena asm_arena_
void resize_tables(TableSizes table_sizes)
Resize tables of patch_point_vals_.
PatchFEValues(unsigned int quad_order, MixedPtr< FiniteElement > fe)
void reinit_patch()
Reinit data.
unsigned int n_dofs() const
Returns the number of shape functions.
void reset()
Reset PatchpointValues structures.
unsigned int get_fe_op(FEOps fe_op) const
Return index of FE operation in operations_ vector.
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)
void set_fe_op(FEOps fe_op, unsigned int op_vec_idx)
Set index of FE operation in operations_ vector.
std::vector< ElOp< spacedim > > operations_
Vector of all defined operations.
Quadrature * get_quadrature() const
Getter for quadrature.
uint dim() const
Getter for dim_.
Tensor tensor_value(uint op_idx, uint point_idx, uint i_dof=0) const
Scalar scalar_value(uint op_idx, uint point_idx, uint i_dof=0) const
Vector vector_value(uint op_idx, uint point_idx, uint i_dof=0) const
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
Quadrature make_from_side(unsigned int sid) const
Definition: quadrature.cc:47
unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:86
unsigned int dim() const
Definition: quadrature.hh:72
Armor::ArmaVec< double, point_dim > point(unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:151
Range helper class.
General point a+ side_begin_ + ccessor allow iterate over quadrature points of given side defined in ...
Definition: eval_subset.hh:116
unsigned int eval_point_idx() const
Return index in EvalPoints object.
Definition: eval_subset.hh:143
SideValues(PatchPointValues< 3 > &patch_point_vals, MixedPtr< FiniteElement > fe)
Constructor.
FeQ< Vector > vector_shape(uint component_idx=0)
Same as BulkValues::vector_shape but register at side quadrature points.
FeQ< Scalar > vector_divergence(uint component_idx=0)
Return the value of the function_no-th vector divergence at the p side quadrature point.
ElQ< Vector > normal_vector()
Register the normal vector to a side at side quadrature points.
std::shared_ptr< FiniteElement< dim > > fe_
FeQ< Tensor > vector_sym_grad(uint component_idx=0)
Return the value of the function_no-th vector symmetric gradient at the p side quadrature point.
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< Tensor > grad_vector_shape(uint component_idx=0)
Return the value of the function_no-th gradient vector shape function at the p bulk quadrature point.
FeQ< Vector > grad_scalar_shape(uint component_idx=0)
Same as BulkValues::grad_scalar_shape but register at side quadrature points.
PatchPointValues< 3 > & patch_point_vals_
ElQ< Scalar > determinant()
Create bulk accessor of jac determinant entity.
NodeAccessor< 3 > node(unsigned int i) const
Returns node for given local index i on the side.
Store finite element data on the actual patch such as shape function values, gradients,...
Eigen::Vector< ArenaVec< uint >, Eigen::Dynamic > IntTableArena
Definition: eigen_tools.hh:39
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
Class FESystem for compound finite elements.
Class FEValues calculates finite element data on the actual cells such as shape function values,...
@ FEScalar
@ FEMixedSystem
@ FEVectorPiola
@ FEVectorContravariant
@ FEVector
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
@ opJac
Jacobian of element.
@ opJxW
JxW value of quadrature point.
@ opNormalVec
normal vector of quadrature point
@ opElJac
Jacobian of element.
Store finite element data on the actual patch such as shape function values, gradients,...
@ opVectorShape
Vector shape operation.
@ opVectorDivergence
Vector divergence.
@ opVectorSymGrad
Vector symmetric gradient.
@ opGradScalarShape
Scalar shape gradient.
@ opGradVectorShape
Vector shape gradient.
@ opScalarShape
Scalar shape operation.
#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.
void copy(const TableSizes &other)
Copy values of other TableSizes instance.
std::vector< std::vector< uint > > elem_sizes_
static void ptop_vector_shape(std::vector< ElOp< 3 >> &operations, std::vector< std::vector< arma::vec3 > > shape_values, uint vector_shape_op_idx)
static void ptop_vector_contravariant_shape(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED std::vector< std::vector< arma::vec3 > > shape_values, FMT_UNUSED uint vector_shape_op_idx)
static void ptop_vector_piola_shape(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED std::vector< std::vector< arma::vec3 > > shape_values, FMT_UNUSED uint vector_shape_op_idx)
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_vector_shape_grads(std::vector< ElOp< 3 >> &operations, std::vector< std::vector< arma::mat > > ref_shape_grads, uint vector_shape_grads_op_idx)
static void ptop_vector_divergence(std::vector< ElOp< 3 >> &operations, uint vector_divergence_op_idx)
Common reinit function of vector divergence on bulk and side points.
static void ptop_vector_sym_grad(std::vector< ElOp< 3 >> &operations, uint vector_sym_grad_op_idx)
Common reinit function of vector symmetric gradient on bulk and side points.
static void ptop_vector_shape(std::vector< ElOp< 3 >> &operations, IntTableArena &el_table, std::vector< std::vector< std::vector< arma::vec3 > > > shape_values, uint vector_shape_op_idx)
static void ptop_scalar_shape(std::vector< ElOp< 3 >> &operations, IntTableArena &el_table, std::vector< std::vector< std::vector< double > > > shape_values, uint scalar_shape_op_idx)
static void ptop_vector_contravariant_shape(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table, FMT_UNUSED std::vector< std::vector< std::vector< arma::vec3 > > > shape_values, FMT_UNUSED uint vector_shape_op_idx)
static void ptop_vector_piola_shape(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table, FMT_UNUSED std::vector< std::vector< std::vector< arma::vec3 > > > shape_values, FMT_UNUSED uint vector_shape_op_idx)
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.