Flow123d  DF_patch_fe_values-dbc06cd
integral_acc.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 integral_acc.hh
15  * @brief
16  * @author David Flanderka
17  */
18 
19 #ifndef INTEGRAL_ACC_HH_
20 #define INTEGRAL_ACC_HH_
21 
22 #include <memory>
23 #include <armadillo>
24 #include "fem/eval_points.hh"
25 #include "fem/integral_points.hh"
26 #include "fem/element_cache_map.hh"
27 #include "mesh/range_wrapper.hh"
28 #include "mesh/accessors.hh"
29 #include "fem/dh_cell_accessor.hh"
30 #include "fem/patch_fe_values.hh"
31 #include "fem/op_function.hh"
32 #include "fem/op_accessors_impl.hh"
33 
34 
35 class BulkIntegral;
36 class EdgeIntegral;
37 class CouplingIntegral;
38 class BoundaryIntegral;
39 template <unsigned int qdim> class BulkIntegralAcc;
40 template <unsigned int qdim> class EdgeIntegralAcc;
41 template <unsigned int qdim> class CouplingIntegralAcc;
42 template <unsigned int qdim> class BoundaryIntegralAcc;
43 
44 
45 namespace internal {
46 
47 /**
48  * Base integral class holds common data members and methods.
49  */
50 class BaseIntegral {
51 public:
52  /// Default constructor
53  BaseIntegral() : eval_points_(nullptr), dim_(0) {}
54 
55  /// Constructor of bulk or side subset
56  BaseIntegral(std::shared_ptr<EvalPoints> eval_points, unsigned int dim)
58 
59  /// Destructor
60  virtual ~BaseIntegral();
61 
62  /// Getter of eval_points
63  std::shared_ptr<EvalPoints> eval_points() const {
64  return eval_points_;
65  }
66 
67  /// Returns dimension.
68  unsigned int dim() const {
69  return dim_;
70  }
71 protected:
72  /// Pointer to EvalPoints
73  std::shared_ptr<EvalPoints> eval_points_;
74  /// Dimension of the cell on which points are placed
75  unsigned int dim_;
76 };
77 
78 /**
79  * Provides methods that allows construction of operations..
80  */
81 template<unsigned int dim>
83 {
84 public:
85  // Default constructor
86  IntegralFactory() : patch_fe_values_(nullptr), element_cache_map_(nullptr), quad_(nullptr)
87  {}
88 
89  // Constructor
90  IntegralFactory(PatchFEValues<3> *pfev, ElementCacheMap *element_cache_map, std::shared_ptr< FiniteElement<dim> > fe, Quadrature *quad)
91  : patch_fe_values_(pfev), element_cache_map_(element_cache_map), fe_(fe), quad_(quad)
92  {}
93 
94  /// Factory method. Creates operation of given OpType.
95  template<class OpType>
97  return patch_fe_values_->get< OpType, dim >(quad_);
98  }
99 
100  /// Factory method. Creates element / side operation of given OpType.
101  template<class OpType>
103  return patch_fe_values_->get_for_elem_quad< OpType, dim >();
104  }
105 
106  /// Factory method. Same as previous but creates FE operation.
107  template<class ValueType, template<unsigned int, class, unsigned int> class OpType, class Domain>
108  FeQArray<ValueType> make_qarray(uint component_idx = 0) {
109  std::shared_ptr<FiniteElement<dim>> fe_component = patch_fe_values_->fe_comp(fe_, component_idx);
110  return FeQArray<ValueType>(patch_fe_values_->template get< OpType<dim, Domain, 3>, dim >(quad_, fe_component));
111  }
112 
115  std::shared_ptr< FiniteElement<dim> > fe_;
117 
118 };
119 
120 } // end of namespace internal
121 
122 
123 namespace internal_integrals {
124 
125 class Base {
126 public:
127  /// Default constructor
128  Base() : dim_(0), quad_(nullptr) {}
129 
130  /// Constructor of bulk or side subset
131  Base(Quadrature *quad, unsigned int dim)
132  : dim_(dim), quad_(quad) {}
133 
134  /// Destructor
135  virtual ~Base()
136  {}
137 
138  /// Returns dimension.
139  unsigned int dim() const {
140  return dim_;
141  }
142 
143  /// Returns quadrature.
144  Quadrature *quad() const {
145  return quad_;
146  }
147 
148  /// Comparison operator
149  bool operator==(const Base& other) const {
150  return (dim_ == other.dim_) && (quad()->size() == other.quad()->size());
151  }
152 protected:
153  /// Dimension of the cell on which points are placed
154  unsigned int dim_;
155  /// Pointer to Quadrature that represents quadrature points of integral.
157 };
158 
159 class Bulk : public Base {
160 public:
162  typedef unsigned int MeshItem;
163 
164  /// Default constructor
165  Bulk() : Base() {}
166 
167  /// Constructor of bulk integral- obsolete constructor
168  Bulk(Quadrature *quad, unsigned int dim, std::shared_ptr<EvalPoints> eval_points, unsigned int subset_idx)
169  : Base(quad, dim) {
170  subset_index_ = subset_idx;
171  begin_idx_ = eval_points->subset_begin(dim_, subset_index_);
172  end_idx_ = eval_points->subset_end(dim_, subset_index_);
173  }
174 
175  /// Destructor
177  {}
178 
179  /// Getter of bulk_begin
180  uint begin_idx() const {
181  return begin_idx_;
182  }
183 
184 private:
185  unsigned int subset_index_; ///< Index of data block according to subset in EvalPoints object.
186  uint begin_idx_; ///< Begin index of quadrature points in EvalPoinnts
187  uint end_idx_; ///< Begin index of quadrature points in EvalPoinnts
188 
189  friend class ::BulkIntegral;
190  friend class ::CouplingIntegral;
191  friend class ::BoundaryIntegral;
192  template <unsigned int qdim>
193  friend class ::BulkIntegralAcc;
194  template <unsigned int qdim>
195  friend class ::CouplingIntegralAcc;
196  template <unsigned int qdim>
197  friend class ::BoundaryIntegralAcc;
198 };
199 
200 class Edge : public Base {
201 public:
204 
205  /// Default constructor
206  Edge() : Base() {}
207 
208  /// Constructor of edge integral
209  Edge(Quadrature *quad, unsigned int dim, std::shared_ptr<EvalPoints> eval_points, unsigned int subset_idx)
210  : Base(quad, dim) {
211  subset_index_ = subset_idx;
212 
213  begin_idx_ = eval_points->subset_begin(dim_, subset_index_);
214  uint end_idx = eval_points->subset_end(dim_, subset_index_);
215  n_sides_ = dim_ + 1;
216  //DebugOut() << "begin: " << begin_idx_ << "end: " << end_idx;
217  n_points_per_side_ = (end_idx - begin_idx_) / n_sides_;
218  //DebugOut() << "points per side: " << n_points_per_side_;
219  }
220 
221  /// Destructor
223  {}
224 
225  inline uint side_begin(const DHCellSide &cell_side) const {
226  return begin_idx_ + cell_side.side_idx() * n_points_per_side_;
227  }
228 
229 private:
230  unsigned int subset_index_;
232 
233  /// Number of sides (value 0 indicates bulk set)
234  unsigned int n_sides_;
235  /// Number of points. TODO: pass this to the constructor, avoid extraction from the eval_points
237 
238  friend class ::EdgeIntegral;
239  friend class ::CouplingIntegral;
240  friend class ::BoundaryIntegral;
241  template <unsigned int qdim>
242  friend class ::EdgeIntegralAcc;
243  template <unsigned int qdim>
244  friend class ::CouplingIntegralAcc;
245  template <unsigned int qdim>
246  friend class ::BoundaryIntegralAcc;
247 };
248 
249 }
250 
251 
252 /**
253  * Integral class of bulk points, allows assemblation of volume integrals.
254  */
256 public:
257  /// Default constructor
259 
260  /// Constructor of bulk integral
261  BulkIntegral(std::shared_ptr<EvalPoints> eval_points, Quadrature *quad, unsigned int dim)
263  {
264  switch (dim) {
265  case 1:
266  internal_bulk_ = eval_points->add_bulk_internal<1>(quad);
267  break;
268  case 2:
269  internal_bulk_ = eval_points->add_bulk_internal<2>(quad);
270  break;
271  case 3:
272  internal_bulk_ = eval_points->add_bulk_internal<3>(quad);
273  break;
274  default:
275  ASSERT_PERMANENT(false);
276  }
277  }
278 
279  /// Destructor
280  ~BulkIntegral();
281 
282  /// Return index of data block according to subset in EvalPoints object
283  inline int get_subset_idx() const {
284  return internal_bulk_->subset_index_;
285  }
286 
287 
288  /// Returns range of bulk local points for appropriate cell accessor - obsolete method
289  inline Range< BulkPoint > points(unsigned int element_patch_idx, const ElementCacheMap *elm_cache_map) const {
290  auto bgn_it = make_iter<BulkPoint>( BulkPoint(elm_cache_map, element_patch_idx, internal_bulk_->begin_idx_));
291  auto end_it = make_iter<BulkPoint>( BulkPoint(elm_cache_map, element_patch_idx, internal_bulk_->end_idx_));
292  return Range<BulkPoint>(bgn_it, end_it);
293  }
294 
295  /// Getter for patch_data_
297  return patch_data_;
298  }
299 
300 protected:
301  /// Internal integral object
302  std::shared_ptr<internal_integrals::Bulk> internal_bulk_;
303  /// Set of integral data of given dimension recomputed for each patch
305 };
306 
307 /**
308  * New Integral accessor class, replace of BulkIntegral, will be merged with BulkIntegral
309  *
310  * IN DEVELOPMENT
311  */
312 template<unsigned int qdim>
314 public:
315  /// Default constructor
317 
318  /// Constructor of bulk integral
319  BulkIntegralAcc(std::shared_ptr<EvalPoints> eval_points, Quadrature *quad, PatchFEValues<3> *pfev, ElementCacheMap *element_cache_map)
320  : BulkIntegral(eval_points, quad, qdim), factory_(pfev, element_cache_map, pfev->fe_dim<qdim>(), quad)
321  {
322  ASSERT_EQ(quad->dim(), qdim);
324  }
325 
326  /// Destructor
328  {}
329 
330  /// Returns range of bulk local points for appropriate cell accessor
331  inline Range< BulkPoint > points(unsigned int element_patch_idx) const {
332  auto bgn_it = make_iter<BulkPoint>( BulkPoint(factory_.element_cache_map_, element_patch_idx, internal_bulk_->begin_idx_));
333  auto end_it = make_iter<BulkPoint>( BulkPoint(factory_.element_cache_map_, element_patch_idx, internal_bulk_->end_idx_));
334  return Range<BulkPoint>(bgn_it, end_it);
335  }
336 
337  /**
338  * @brief Register the product of Jacobian determinant and the quadrature
339  * weight at bulk quadrature points.
340  *
341  * @param quad Quadrature.
342  */
343  inline FeQ<Scalar> JxW()
344  {
345  return FeQ<Scalar>(factory_.template make_patch_op< Op::JxW<qdim, Op::BulkDomain, 3> >());
346  }
347 
348  /// Create bulk accessor of coords entity
350  {
351  return FeQ<Vector>(factory_.template make_patch_op< Op::PtCoords<qdim, Op::BulkDomain, 3> >());
352  }
353 
354 // inline ElQ<Tensor> jacobian(std::initializer_list<Quadrature *> quad_list)
355 // {}
356 
357  /// Create bulk accessor of jac determinant entity
359  {
360  return ElQ<Scalar>( factory_.template make_elem_patch_op< Op::JacDet<qdim, Op::BulkDomain, 3> >() );
361  }
362 
363  /**
364  * @brief Return the value of the @p function_no-th shape function at
365  * the @p p bulk quadrature point.
366  *
367  * @param component_idx Number of the shape function.
368  */
369  inline FeQArray<Scalar> scalar_shape(uint component_idx = 0)
370  {
371  return factory_.template make_qarray<Scalar, Op::ScalarShape, Op::BulkDomain>(component_idx);
372  }
373 
374  /**
375  * @brief Return the value of the @p function_no-th vector shape function at
376  * the @p p bulk quadrature point.
377  *
378  * @param component_idx Number of the shape function.
379  */
380  inline FeQArray<Vector> vector_shape(uint component_idx = 0)
381  {
382  return factory_.template make_qarray<Vector, Op::DispatchVectorShape, Op::BulkDomain>(component_idx);
383  }
384 
385 // inline FeQArray<Tensor> tensor_shape(uint component_idx = 0)
386 // {}
387 
388  /**
389  * @brief Return the value of the @p function_no-th gradient shape function at
390  * the @p p bulk quadrature point.
391  *
392  * @param component_idx Number of the shape function.
393  */
394  inline FeQArray<Vector> grad_scalar_shape(uint component_idx=0)
395  {
396  return factory_.template make_qarray<Vector, Op::GradScalarShape, Op::BulkDomain>(component_idx);
397  }
398 
399  /**
400  * @brief Return the value of the @p function_no-th gradient vector shape function
401  * at the @p p bulk quadrature point.
402  *
403  * @param component_idx Number of the shape function.
404  */
405  inline FeQArray<Tensor> grad_vector_shape(uint component_idx=0)
406  {
407  return factory_.template make_qarray<Tensor, Op::DispatchGradVectorShape, Op::BulkDomain>(component_idx);
408  }
409 
410  /**
411  * @brief Return the value of the @p function_no-th vector symmetric gradient
412  * at the @p p bulk quadrature point.
413  *
414  * @param component_idx Number of the shape function.
415  */
416  inline FeQArray<Tensor> vector_sym_grad(uint component_idx=0)
417  {
418  return factory_.template make_qarray<Tensor, Op::VectorSymGrad, Op::BulkDomain>(component_idx);
419  }
420 
421  /**
422  * @brief Return the value of the @p function_no-th vector divergence at
423  * the @p p bulk quadrature point.
424  *
425  * @param component_idx Number of the shape function.
426  */
427  inline FeQArray<Scalar> vector_divergence(uint component_idx=0)
428  {
429  return factory_.template make_qarray<Scalar, Op::VectorDivergence, Op::BulkDomain>(component_idx);
430  }
431 
432 private:
433  /// Defines interface of operation accessors declaration
435 };
436 
437 /**
438  * Integral class of side points, allows assemblation of element - element fluxes.
439  */
441 public:
442  /// Default constructor
444  {
445  ASSERT_PERMANENT(false);
446  }
447 
448  /// Constructor of edge integral
449  EdgeIntegral(std::shared_ptr<EvalPoints> eval_points, Quadrature *quad, unsigned int dim)
451  {
452  switch (dim) {
453  case 1:
454  internal_edge_ = eval_points->add_edge_internal<1>(quad);
455  break;
456  case 2:
457  internal_edge_ = eval_points->add_edge_internal<2>(quad);
458  break;
459  case 3:
460  internal_edge_ = eval_points->add_edge_internal<3>(quad);
461  break;
462  default:
463  ASSERT(false).error("Should not happen!\n");
464  }
465  }
466 
467  /// Destructor
468  ~EdgeIntegral();
469 
470  /// Getter of n_sides
471  inline unsigned int n_sides() const {
472  return internal_edge_->n_sides_;
473  }
474 
475  /// Return index of data block according to subset in EvalPoints object
476  inline int get_subset_idx() const {
477  return internal_edge_->subset_index_;
478  }
479 
480  inline uint side_begin(const DHCellSide &cell_side) const {
481  return internal_edge_->side_begin(cell_side);
482  }
483 
484  /// Returns range of side local points for appropriate cell side accessor - obsolete method
485  inline Range< EdgePoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const {
486  ASSERT_EQ(cell_side.dim(), dim_);
487 
488  //DebugOut() << "points per side: " << internal_edge_->n_points_per_side_;
489  uint element_patch_idx = elm_cache_map->position_in_cache(cell_side.element().idx());
490  uint begin_idx = internal_edge_->side_begin(cell_side);
491  auto bgn_it = make_iter<EdgePoint>( EdgePoint(
492  BulkPoint(elm_cache_map, element_patch_idx, 0), internal_edge_, begin_idx));
493  auto end_it = make_iter<EdgePoint>( EdgePoint(
494  BulkPoint(elm_cache_map, element_patch_idx, internal_edge_->n_points_per_side_), internal_edge_, begin_idx));
495  return Range<EdgePoint>(bgn_it, end_it);
496  }
497 
498  /// Getter for patch_data_
500  return patch_data_;
501  }
502 
503 
504 protected:
505  /// Internal integral object
506  std::shared_ptr<internal_integrals::Edge> internal_edge_;
507  /// Set of integral data of given dimension recomputed for each patch
509 
510  friend class EvalPoints;
511  friend class EdgePoint;
512  friend class CouplingPoint;
513  friend class BoundaryPoint;
514  friend class CouplingIntegral;
515  friend class BoundaryIntegral;
516  template<unsigned int qdim>
517  friend class CouplingIntegralAcc;
518  template<unsigned int qdim>
519  friend class BoundaryIntegralAcc;
520 };
521 
522 /**
523  * New Integral accessor class, replace of EdgeIntegral, will be merged with EdgeIntegral
524  *
525  * IN DEVELOPMENT
526  */
527 template<unsigned int qdim>
529 public:
530  /// Default constructor
532 
533  /// Constructor of edge integral
534  EdgeIntegralAcc(std::shared_ptr<EvalPoints> eval_points, Quadrature *quad, PatchFEValues<3> *pfev, ElementCacheMap *element_cache_map)
535  : EdgeIntegral(eval_points, quad, qdim), factory_(pfev, element_cache_map, pfev->fe_dim<qdim>(), quad)
536  {
537  ASSERT_EQ(quad->dim()+1, qdim);
539  }
540 
541 
542  /// Destructor
544  {}
545 
546  /// Returns range of side local points for appropriate cell side accessor
547  inline Range< EdgePoint > points(const DHCellSide &cell_side) const {
548  ASSERT_EQ(cell_side.dim(), dim_);
549 
550  //DebugOut() << "points per side: " << internal_edge_->n_points_per_side_;
551  uint element_patch_idx = factory_.element_cache_map_->position_in_cache(cell_side.element().idx());
552  uint begin_idx = internal_edge_->side_begin(cell_side);
553  auto bgn_it = make_iter<EdgePoint>( EdgePoint(
554  BulkPoint(factory_.element_cache_map_, element_patch_idx, 0), this->internal_edge_, begin_idx));
555  auto end_it = make_iter<EdgePoint>( EdgePoint(
556  BulkPoint(factory_.element_cache_map_, element_patch_idx, internal_edge_->n_points_per_side_), this->internal_edge_, begin_idx));
557  return Range<EdgePoint>(bgn_it, end_it);
558  }
559 
560  /// Same as BulkValues::JxW but register at side quadrature points.
561  inline FeQ<Scalar> JxW()
562  {
563  return FeQ<Scalar>(factory_.template make_patch_op< Op::JxW<qdim, Op::SideDomain, 3> >());
564  }
565 
566  /**
567  * @brief Register the normal vector to a side at side quadrature points.
568  *
569  * @param quad Quadrature.
570  */
572  {
573  return ElQ<Vector>(factory_.template make_elem_patch_op< Op::NormalVec<qdim, 3> >());
574  }
575 
576  /// Create side accessor of coords entity
578  {
579  return FeQ<Vector>(factory_.template make_patch_op< Op::PtCoords<qdim, Op::SideDomain, 3> >());
580  }
581 
582  /// Create bulk accessor of jac determinant entity
584  {
585  return ElQ<Scalar>(factory_.template make_patch_op< Op::JacDet<qdim, Op::SideDomain, 3> >());
586  }
587 
588  /// Same as BulkValues::scalar_shape but register at side quadrature points.
589  inline FeQArray<Scalar> scalar_shape(uint component_idx = 0)
590  {
591  return factory_.template make_qarray<Scalar, Op::ScalarShape, Op::SideDomain>(component_idx);
592  }
593 
594  /// Same as BulkValues::vector_shape but register at side quadrature points.
595  inline FeQArray<Vector> vector_shape(uint component_idx = 0)
596  {
597  return factory_.template make_qarray<Vector, Op::DispatchVectorShape, Op::SideDomain>(component_idx);
598  }
599 
600  /// Same as BulkValues::grad_scalar_shape but register at side quadrature points.
601  inline FeQArray<Vector> grad_scalar_shape(uint component_idx=0)
602  {
603  return factory_.template make_qarray<Vector, Op::GradScalarShape, Op::SideDomain>(component_idx);
604  }
605 
606  /**
607  * @brief Return the value of the @p function_no-th gradient vector shape function
608  * at the @p p bulk quadrature point.
609  *
610  * @param component_idx Number of the shape function.
611  */
612  inline FeQArray<Tensor> grad_vector_shape(uint component_idx=0)
613  {
614  return factory_.template make_qarray<Tensor, Op::DispatchGradVectorShape, Op::SideDomain>(component_idx);
615  }
616 
617  /**
618  * @brief Return the value of the @p function_no-th vector symmetric gradient
619  * at the @p p side quadrature point.
620  *
621  * @param component_idx Number of the shape function.
622  */
623  inline FeQArray<Tensor> vector_sym_grad(uint component_idx=0)
624  {
625  return factory_.template make_qarray<Tensor, Op::VectorSymGrad, Op::SideDomain>(component_idx);
626  }
627 
628  /**
629  * @brief Return the value of the @p function_no-th vector divergence at
630  * the @p p side quadrature point.
631  *
632  * @param component_idx Number of the shape function.
633  */
634  inline FeQArray<Scalar> vector_divergence(uint component_idx=0)
635  {
636  return factory_.template make_qarray<Scalar, Op::VectorDivergence, Op::SideDomain>(component_idx);
637  }
638 
639 private:
640  /// Defines interface of operation accessors declaration
642 
643  friend class EvalPoints;
644  friend class EdgePoint;
645  friend class CouplingPoint;
646  friend class BoundaryPoint;
647  template <unsigned int quaddim>
648  friend class CouplingIntegralAcc;
649  template <unsigned int quaddim>
650  friend class BoundaryIntegralAcc;
651 };
652 
653 
654 /**
655  * Integral class of neighbour points, allows assemblation of element - side fluxes.
656  *
657  * Dimension corresponds with element of higher dim.
658  */
660 public:
661  /// Default constructor
663 
664  /// Constructor of ngh integral
665  CouplingIntegral(std::shared_ptr<EvalPoints> eval_points, Quadrature *quad, unsigned int dim)
667  {
668  switch (dim) {
669  case 1:
670  internal_bulk_ = eval_points->add_bulk_internal<1>(quad);
671  internal_edge_ = eval_points->add_edge_internal<2>(quad);
672  break;
673  case 2:
674  internal_bulk_ = eval_points->add_bulk_internal<2>(quad);
675  internal_edge_ = eval_points->add_edge_internal<3>(quad);
676  break;
677  default:
678  ASSERT(false)(dim).error("Should not happen!\n");
679  }
680  }
681 
682  /// Destructor
684 
685  /// Return index of data block according to subset of higher dim in EvalPoints object
686  inline int get_subset_high_idx() const {
687  return internal_edge_->subset_index_;
688  }
689 
690  /// Return index of data block according to subset of lower dim in EvalPoints object
691  inline int get_subset_low_idx() const {
692  return internal_bulk_->subset_index_;
693  }
694 
695  inline uint bulk_begin() const {
696  return eval_points_->subset_begin(dim_-1, internal_bulk_->subset_index_);
697  }
698 
699  /// Returns range of side local points for appropriate cell side accessor - obsolete method
700  inline Range< CouplingPoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const {
701  ASSERT_EQ(cell_side.dim(), dim_+1);
702  uint element_patch_idx = elm_cache_map->position_in_cache(cell_side.element().idx());
703  uint side_begin = internal_edge_->side_begin(cell_side);
704  auto bgn_it = make_iter<CouplingPoint>( CouplingPoint(
705  BulkPoint(elm_cache_map, element_patch_idx, 0), internal_bulk_, side_begin) );
706  auto end_it = make_iter<CouplingPoint>( CouplingPoint(
707  BulkPoint(elm_cache_map, element_patch_idx, internal_edge_->n_points_per_side_), internal_bulk_, side_begin) );;
708  return Range<CouplingPoint>(bgn_it, end_it);
709  }
710 
711  /// Getter for patch_data_
713  return patch_data_;
714  }
715 
716 protected:
717  /// Integral according to bulk subset part (element of lower dim) in EvalPoints object.
718  std::shared_ptr<internal_integrals::Bulk> internal_bulk_;
719  /// Integral according to side subset part (element of higher dim) in EvalPoints object.
720  std::shared_ptr<internal_integrals::Edge> internal_edge_;
721  /// Set of integral data of given dimension recomputed for each patch
723 
724  friend class CouplingPoint;
725 };
726 
727 template<unsigned int qdim>
729 public:
730  /// Default constructor
732 
733  /// Constructor of ngh integral
734  CouplingIntegralAcc(std::shared_ptr<EvalPoints> eval_points, Quadrature *quad, PatchFEValues<3> *pfev, ElementCacheMap *element_cache_map)
735  : CouplingIntegral(eval_points, quad, qdim),
736  factory_(pfev, element_cache_map, pfev->fe_dim<qdim>(), quad),
737  factory_high_(pfev, element_cache_map, pfev->fe_dim<qdim+1>(), quad)
738  {
739  ASSERT_EQ(quad->dim(), qdim);
742  }
743 
744  /// Destructor
746  {
747  internal_bulk_.reset();
748  internal_edge_.reset();
749  }
750 
751  /// Return index of data block according to subset of higher dim in EvalPoints object
752  inline int get_subset_high_idx() const {
753  return internal_edge_->subset_index_;
754  }
755 
756  /// Return index of data block according to subset of lower dim in EvalPoints object
757  inline int get_subset_low_idx() const {
758  return internal_bulk_->subset_index_;
759  }
760 
761  inline uint bulk_begin() const {
762  return internal_bulk_->begin_idx_;
763  }
764 
765  /// Returns range of side local points for appropriate cell side accessor
766  inline Range< CouplingPoint > points(const DHCellSide &cell_side) const {
767  ASSERT_EQ(cell_side.dim(), dim_+1);
768 
769  uint element_patch_idx = factory_.element_cache_map_->position_in_cache(cell_side.element().idx());
770  uint side_begin = internal_edge_->side_begin(cell_side);
771  auto bgn_it = make_iter<CouplingPoint>( CouplingPoint(
772  BulkPoint(factory_.element_cache_map_, element_patch_idx, 0), internal_bulk_, side_begin) );
773  auto end_it = make_iter<CouplingPoint>( CouplingPoint(
774  BulkPoint(factory_.element_cache_map_, element_patch_idx, internal_edge_->n_points_per_side_), internal_bulk_, side_begin) );;
775  return Range<CouplingPoint>(bgn_it, end_it);
776  }
777 
778  /// Factory method. Same as previous but creates FE operation.
779  template<class ValueType, template<unsigned int, class, unsigned int> class OpType>
780  FeQJoin<ValueType> make_qjoin(uint component_idx = 0) {
781  // element of lower dim (bulk points)
782  auto fe_component_low = factory_.patch_fe_values_->fe_comp(factory_.fe_, component_idx);
783  auto *low_dim_op = factory_.patch_fe_values_->template get< OpType<qdim, Op::BulkDomain, 3>, qdim >(factory_.quad_, fe_component_low);
784  auto *low_dim_zero_op = factory_.patch_fe_values_->template get< Op::OpZero<qdim, Op::BulkDomain, 3>, qdim >(factory_.quad_, fe_component_low);
785 
786  // element of higher dim (side points)
787  auto fe_component_high = factory_.patch_fe_values_->fe_comp(factory_high_.fe_, component_idx);
788  auto *high_dim_op = factory_.patch_fe_values_->template get< OpType<qdim+1, Op::SideDomain, 3>, qdim+1 >(factory_.quad_, fe_component_high);
789  auto *high_dim_zero_op = factory_.patch_fe_values_->template get< Op::OpZero<qdim+1, Op::SideDomain, 3>, qdim+1 >(factory_.quad_, fe_component_high);
790 
791  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");
792  return FeQJoin<ValueType>(low_dim_op, high_dim_op, low_dim_zero_op, high_dim_zero_op);
793  }
794 
795  /// Same as BulkValues::JxW but register at side quadrature points.
796  inline FeQ<Scalar> JxW()
797  {
798  return FeQ<Scalar>(factory_high_.template make_patch_op< Op::JxW<qdim+1, Op::SideDomain, 3> >());
799  }
800 
801  /**
802  * @brief Register the normal vector to a side at side quadrature points.
803  *
804  * @param quad Quadrature.
805  */
807  {
808  return ElQ<Vector>(factory_high_.template make_elem_patch_op< Op::NormalVec<qdim+1, 3> >());
809  }
810 
811  inline FeQArray<Vector> vector_shape(uint component_idx = 0)
812  {
813  return factory_high_.template make_qarray<Vector, Op::DispatchVectorShape, Op::SideDomain>(component_idx); // error
814  }
815 
816  inline FeQJoin<Scalar> scalar_join_shape(uint component_idx = 0)
817  {
818  return this->template make_qjoin<Scalar, Op::ScalarShape>(component_idx);
819  }
820 
821  inline FeQJoin<Vector> vector_join_shape(uint component_idx = 0)
822  {
823  return this->template make_qjoin<Vector, Op::DispatchVectorShape>(component_idx);
824  }
825 
827  {
828  return this->template make_qjoin<Tensor, Op::DispatchGradVectorShape>(component_idx);
829  }
830 
831 
832 private:
833  /// Defines interface of operation accessors declaration
835 
836  /// Same as prefious but for element of higher dim
838 
839  friend class CouplingPoint;
840 };
841 
842 /// Template specialization of previous class
843 template<>
845 public:
846  /// Default constructor
848 
849  /// Constructor of ngh integral
850  CouplingIntegralAcc(std::shared_ptr<EvalPoints> eval_points, Quadrature *quad, PatchFEValues<3> *pfev, ElementCacheMap *element_cache_map)
851  : CouplingIntegral(),
852  factory_(pfev, element_cache_map, pfev->fe_dim<3>(), quad)
853  {
854  ASSERT_EQ(quad->dim(), 3);
855  this->eval_points_ = eval_points;
856  this->dim_ = 3;
857  }
858 
859  /// Destructor
861  {}
862 
863  /// Returns empty point range
864  inline Range< CouplingPoint > points(FMT_UNUSED const DHCellSide &cell_side) const {
865  auto iter = make_iter<CouplingPoint>( CouplingPoint() );
866  return Range<CouplingPoint>(iter, iter);
867  }
868 
869  /// Define empty operations
870  inline FeQ<Scalar> JxW()
871  {
872  return FeQ<Scalar>();
873  }
874 
876  {
877  return ElQ<Vector>();
878  }
879 
880  inline FeQArray<Vector> vector_shape(FMT_UNUSED uint component_idx = 0)
881  {
882  return FeQArray<Vector>();
883  }
884 
886  {
887  return FeQJoin<Scalar>();
888  }
889 
891  {
892  return FeQJoin<Vector>();
893  }
894 
896  {
897  return FeQJoin<Tensor>();
898  }
899 
900 private:
901  /// Defines interface of operation accessors declaration
903 };
904 
905 
906 /**
907  * Integral class of boundary points, allows assemblation of fluxes between sides and neighbouring boundary elements.
908  */
910 public:
911  /// Default constructor
912  BoundaryIntegral() : internal::BaseIntegral() {}
913 
914  /// Constructor of bulk subset
915  BoundaryIntegral(std::shared_ptr<EvalPoints> eval_points, Quadrature *quad, unsigned int dim);
916 
917  /// Destructor
918  ~BoundaryIntegral();
919 
920  /// Return index of data block according to subset of higher dim in EvalPoints object
921  inline int get_subset_high_idx() const {
922  return internal_edge_->subset_index_;
923  }
924 
925  /// Return index of data block according to subset of lower dim (boundary) in EvalPoints object
926  inline int get_subset_low_idx() const {
927  return internal_bulk_->subset_index_;
928  }
929 
930  inline uint bulk_begin() const {
931  // DebugOut().fmt("edge_begin: {} bdr_begin: {}",
932  // internal_bulk_->get_begin_idx(),
933  // internal_bulk_->get_begin_idx());
934  return internal_bulk_->begin_idx_;
935  }
936 
937  /// Returns range of bulk local points for appropriate cell accessor - obsolete method
938  inline Range< BoundaryPoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const {
939  ASSERT_EQ(cell_side.dim(), dim_);
940  uint element_patch_idx = elm_cache_map->position_in_cache(cell_side.element().idx());
941  uint side_begin = internal_edge_->side_begin(cell_side);
942  auto bgn_it = make_iter<BoundaryPoint>( BoundaryPoint(
943  BulkPoint(elm_cache_map, element_patch_idx, 0), internal_bulk_, side_begin) );
944  auto end_it = make_iter<BoundaryPoint>( BoundaryPoint(
945  BulkPoint(elm_cache_map, element_patch_idx, internal_edge_->n_points_per_side_), internal_bulk_, side_begin) );;
946  return Range<BoundaryPoint>(bgn_it, end_it);
947  }
948 
949  /// Getter for patch_data_
951  return patch_data_;
952  }
953 
954 protected:
955  /// Integral according to kower dim (boundary) element subset part in EvalPoints object.
956  std::shared_ptr<internal_integrals::Bulk> internal_bulk_;
957  /// Integral according to higher dim (bulk) element subset part in EvalPoints object.
958  std::shared_ptr<internal_integrals::Edge> internal_edge_;
959  /// Set of integral data of given dimension recomputed for each patch
961 
962  friend class BoundaryPoint;
963 };
964 
965 
966 template<unsigned int qdim>
968 public:
969  /// Default constructor
971 
972  BoundaryIntegralAcc(std::shared_ptr<EvalPoints> eval_points, Quadrature *quad, PatchFEValues<3> *pfev, ElementCacheMap *element_cache_map)
973  : BoundaryIntegral(eval_points, quad, qdim), factory_(pfev, element_cache_map, pfev->fe_dim<qdim>(), quad)
974  {
975  ASSERT_EQ(quad->dim()+1, qdim);
977  }
978 
979  /// Destructor
981  {}
982 
983  /// Returns range of bulk local points for appropriate cell accessor
984  inline Range< BoundaryPoint > points(const DHCellSide &cell_side) const {
985  ASSERT_EQ(cell_side.dim(), dim_);
986 
987  uint element_patch_idx = factory_.element_cache_map_->position_in_cache(cell_side.element().idx());
988  uint side_begin = internal_edge_->side_begin(cell_side);
989  auto bgn_it = make_iter<BoundaryPoint>( BoundaryPoint(
990  BulkPoint(factory_.element_cache_map_, element_patch_idx, 0), internal_bulk_, side_begin) );
991  auto end_it = make_iter<BoundaryPoint>( BoundaryPoint(
992  BulkPoint(factory_.element_cache_map_, element_patch_idx, internal_edge_->n_points_per_side_), internal_bulk_, side_begin) );;
993  return Range<BoundaryPoint>(bgn_it, end_it);
994  }
995 
996  /// Same as BulkValues::JxW but register at side quadrature points.
997  inline FeQ<Scalar> JxW()
998  {
999  return FeQ<Scalar>(factory_.template make_patch_op< Op::JxW<qdim, Op::SideDomain, 3> >());
1000  }
1001 
1002  /**
1003  * @brief Register the normal vector to a side at side quadrature points.
1004  *
1005  * @param quad Quadrature.
1006  */
1008  {
1009  return ElQ<Vector>(factory_.template make_elem_patch_op< Op::NormalVec<qdim, 3> >());
1010  }
1011 
1012  /// Create side accessor of coords entity
1014  {
1015  return FeQ<Vector>(factory_.template make_patch_op< Op::PtCoords<qdim, Op::SideDomain, 3> >());
1016  }
1017 
1018  /// Create bulk accessor of jac determinant entity
1020  {
1021  return ElQ<Scalar>(factory_.template make_patch_op< Op::JacDet<qdim, Op::SideDomain, 3> >());
1022  }
1023 
1024  /// Same as BulkValues::scalar_shape but register at side quadrature points.
1025  inline FeQArray<Scalar> scalar_shape(uint component_idx = 0)
1026  {
1027  return factory_.template make_qarray<Scalar, Op::ScalarShape, Op::SideDomain>(component_idx);
1028  }
1029 
1030  /// Same as BulkValues::vector_shape but register at side quadrature points.
1031  inline FeQArray<Vector> vector_shape(uint component_idx = 0)
1032  {
1033  return factory_.template make_qarray<Vector, Op::DispatchVectorShape, Op::SideDomain>(component_idx);
1034  }
1035 
1036  /// Same as BulkValues::grad_scalar_shape but register at side quadrature points.
1037  inline FeQArray<Vector> grad_scalar_shape(uint component_idx=0)
1038  {
1039  return factory_.template make_qarray<Vector, Op::GradScalarShape, Op::SideDomain>(component_idx);
1040  }
1041 
1042  /**
1043  * @brief Return the value of the @p function_no-th gradient vector shape function
1044  * at the @p p bulk quadrature point.
1045  *
1046  * @param component_idx Number of the shape function.
1047  */
1048  inline FeQArray<Tensor> grad_vector_shape(uint component_idx=0)
1049  {
1050  return factory_.template make_qarray<Tensor, Op::DispatchGradVectorShape, Op::SideDomain>(component_idx);
1051  }
1052 
1053  /**
1054  * @brief Return the value of the @p function_no-th vector symmetric gradient
1055  * at the @p p side quadrature point.
1056  *
1057  * @param component_idx Number of the shape function.
1058  */
1059  inline FeQArray<Tensor> vector_sym_grad(uint component_idx=0)
1060  {
1061  return factory_.template make_qarray<Tensor, Op::VectorSymGrad, Op::SideDomain>(component_idx);
1062  }
1063 
1064  /**
1065  * @brief Return the value of the @p function_no-th vector divergence at
1066  * the @p p side quadrature point.
1067  *
1068  * @param component_idx Number of the shape function.
1069  */
1070  inline FeQArray<Scalar> vector_divergence(uint component_idx=0)
1071  {
1072  return factory_.template make_qarray<Scalar, Op::VectorDivergence, Op::SideDomain>(component_idx);
1073  }
1074 
1075 private:
1076  /// Defines interface of operation accessors declaration
1078 
1079  friend class BoundaryPoint;
1080 };
1081 
1082 
1083 #endif /* INTEGRAL_ACC_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
FeQArray< Vector > grad_scalar_shape(uint component_idx=0)
Same as BulkValues::grad_scalar_shape but register at side quadrature points.
internal::IntegralFactory< qdim > factory_
Defines interface of operation accessors declaration.
Range< BoundaryPoint > points(const DHCellSide &cell_side) const
Returns range of bulk local points for appropriate cell accessor.
BoundaryIntegralAcc()
Default constructor.
FeQArray< 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.
~BoundaryIntegralAcc()
Destructor.
FeQ< Vector > coords()
Create side accessor of coords entity.
FeQArray< 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.
FeQArray< Vector > vector_shape(uint component_idx=0)
Same as BulkValues::vector_shape but register at side quadrature points.
FeQArray< Scalar > scalar_shape(uint component_idx=0)
Same as BulkValues::scalar_shape but register at side quadrature points.
FeQ< Scalar > JxW()
Same as BulkValues::JxW but register at side quadrature points.
FeQArray< Scalar > vector_divergence(uint component_idx=0)
Return the value of the function_no-th vector divergence at the p side quadrature point.
friend class BoundaryPoint
ElQ< Scalar > determinant()
Create bulk accessor of jac determinant entity.
ElQ< Vector > normal_vector()
Register the normal vector to a side at side quadrature points.
BoundaryIntegralAcc(std::shared_ptr< EvalPoints > eval_points, Quadrature *quad, PatchFEValues< 3 > *pfev, ElementCacheMap *element_cache_map)
RevertableList< BoundaryIntegralData > & patch_data()
Getter for patch_data_.
Range< BoundaryPoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const
Returns range of bulk local points for appropriate cell accessor - obsolete method.
RevertableList< BoundaryIntegralData > patch_data_
Set of integral data of given dimension recomputed for each patch.
int get_subset_low_idx() const
Return index of data block according to subset of lower dim (boundary) in EvalPoints object.
BoundaryIntegral()
Default constructor.
std::shared_ptr< internal_integrals::Bulk > internal_bulk_
Integral according to kower dim (boundary) element subset part in EvalPoints object.
friend class BoundaryPoint
int get_subset_high_idx() const
Return index of data block according to subset of higher dim in EvalPoints object.
std::shared_ptr< internal_integrals::Edge > internal_edge_
Integral according to higher dim (bulk) element subset part in EvalPoints object.
uint bulk_begin() const
Point accessor allow iterate over quadrature points of given side defined in local element coordinate...
BulkIntegralAcc(std::shared_ptr< EvalPoints > eval_points, Quadrature *quad, PatchFEValues< 3 > *pfev, ElementCacheMap *element_cache_map)
Constructor of bulk integral.
FeQArray< Scalar > scalar_shape(uint component_idx=0)
Return the value of the function_no-th shape function at the p bulk quadrature point.
FeQArray< 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 > JxW()
Register the product of Jacobian determinant and the quadrature weight at bulk quadrature points.
ElQ< Scalar > determinant()
Create bulk accessor of jac determinant entity.
FeQArray< 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.
FeQArray< 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.
FeQ< Vector > coords()
Create bulk accessor of coords entity.
Range< BulkPoint > points(unsigned int element_patch_idx) const
Returns range of bulk local points for appropriate cell accessor.
FeQArray< Vector > vector_shape(uint component_idx=0)
Return the value of the function_no-th vector shape function at the p bulk quadrature point.
BulkIntegralAcc()
Default constructor.
internal::IntegralFactory< qdim > factory_
Defines interface of operation accessors declaration.
~BulkIntegralAcc()
Destructor.
FeQArray< Scalar > vector_divergence(uint component_idx=0)
Return the value of the function_no-th vector divergence at the p bulk quadrature point.
int get_subset_idx() const
Return index of data block according to subset in EvalPoints object.
~BulkIntegral()
Destructor.
Definition: integral_acc.cc:38
Range< BulkPoint > points(unsigned int element_patch_idx, const ElementCacheMap *elm_cache_map) const
Returns range of bulk local points for appropriate cell accessor - obsolete method.
RevertableList< BulkIntegralData > & patch_data()
Getter for patch_data_.
BulkIntegral()
Default constructor.
BulkIntegral(std::shared_ptr< EvalPoints > eval_points, Quadrature *quad, unsigned int dim)
Constructor of bulk integral.
RevertableList< BulkIntegralData > patch_data_
Set of integral data of given dimension recomputed for each patch.
std::shared_ptr< internal_integrals::Bulk > internal_bulk_
Internal integral object.
Base point accessor class.
Template specialization of previous class.
CouplingIntegralAcc()
Default constructor.
ElQ< Vector > normal_vector()
FeQJoin< Scalar > scalar_join_shape(FMT_UNUSED uint component_idx=0)
CouplingIntegralAcc(std::shared_ptr< EvalPoints > eval_points, Quadrature *quad, PatchFEValues< 3 > *pfev, ElementCacheMap *element_cache_map)
Constructor of ngh integral.
FeQArray< Vector > vector_shape(FMT_UNUSED uint component_idx=0)
FeQJoin< Vector > vector_join_shape(FMT_UNUSED uint component_idx=0)
~CouplingIntegralAcc()
Destructor.
FeQ< Scalar > JxW()
Define empty operations.
FeQJoin< Tensor > gradient_vector_join_shape(FMT_UNUSED uint component_idx=0)
Range< CouplingPoint > points(FMT_UNUSED const DHCellSide &cell_side) const
Returns empty point range.
internal::IntegralFactory< 3 > factory_
Defines interface of operation accessors declaration.
FeQJoin< Scalar > scalar_join_shape(uint component_idx=0)
FeQJoin< Tensor > gradient_vector_join_shape(uint component_idx=0)
CouplingIntegralAcc(std::shared_ptr< EvalPoints > eval_points, Quadrature *quad, PatchFEValues< 3 > *pfev, ElementCacheMap *element_cache_map)
Constructor of ngh integral.
~CouplingIntegralAcc()
Destructor.
int get_subset_low_idx() const
Return index of data block according to subset of lower dim in EvalPoints object.
FeQJoin< ValueType > make_qjoin(uint component_idx=0)
Factory method. Same as previous but creates FE operation.
FeQ< Scalar > JxW()
Same as BulkValues::JxW but register at side quadrature points.
ElQ< Vector > normal_vector()
Register the normal vector to a side at side quadrature points.
internal::IntegralFactory< qdim+1 > factory_high_
Same as prefious but for element of higher dim.
uint bulk_begin() const
Range< CouplingPoint > points(const DHCellSide &cell_side) const
Returns range of side local points for appropriate cell side accessor.
FeQJoin< Vector > vector_join_shape(uint component_idx=0)
CouplingIntegralAcc()
Default constructor.
internal::IntegralFactory< qdim > factory_
Defines interface of operation accessors declaration.
int get_subset_high_idx() const
Return index of data block according to subset of higher dim in EvalPoints object.
friend class CouplingPoint
FeQArray< Vector > vector_shape(uint component_idx=0)
RevertableList< CouplingIntegralData > patch_data_
Set of integral data of given dimension recomputed for each patch.
Range< CouplingPoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const
Returns range of side local points for appropriate cell side accessor - obsolete method.
int get_subset_high_idx() const
Return index of data block according to subset of higher dim in EvalPoints object.
uint bulk_begin() const
int get_subset_low_idx() const
Return index of data block according to subset of lower dim in EvalPoints object.
CouplingIntegral()
Default constructor.
std::shared_ptr< internal_integrals::Edge > internal_edge_
Integral according to side subset part (element of higher dim) in EvalPoints object.
std::shared_ptr< internal_integrals::Bulk > internal_bulk_
Integral according to bulk subset part (element of lower dim) in EvalPoints object.
~CouplingIntegral()
Destructor.
Definition: integral_acc.cc:55
CouplingIntegral(std::shared_ptr< EvalPoints > eval_points, Quadrature *quad, unsigned int dim)
Constructor of ngh integral.
RevertableList< CouplingIntegralData > & patch_data()
Getter for patch_data_.
friend class CouplingPoint
Point accessor allow iterate over quadrature points of given side defined in local element coordinate...
Side accessor allows to iterate over sides of DOF handler cell.
unsigned int dim() const
Return dimension of element appropriate to the side.
ElementAccessor< 3 > element() const
unsigned int side_idx() const
FeQArray< Scalar > scalar_shape(uint component_idx=0)
Same as BulkValues::scalar_shape but register at side quadrature points.
internal::IntegralFactory< qdim > factory_
Defines interface of operation accessors declaration.
EdgeIntegralAcc(std::shared_ptr< EvalPoints > eval_points, Quadrature *quad, PatchFEValues< 3 > *pfev, ElementCacheMap *element_cache_map)
Constructor of edge integral.
ElQ< Vector > normal_vector()
Register the normal vector to a side at side quadrature points.
FeQArray< Vector > grad_scalar_shape(uint component_idx=0)
Same as BulkValues::grad_scalar_shape but register at side quadrature points.
~EdgeIntegralAcc()
Destructor.
Range< EdgePoint > points(const DHCellSide &cell_side) const
Returns range of side local points for appropriate cell side accessor.
ElQ< Scalar > determinant()
Create bulk accessor of jac determinant entity.
FeQArray< Scalar > vector_divergence(uint component_idx=0)
Return the value of the function_no-th vector divergence at the p side quadrature point.
EdgeIntegralAcc()
Default constructor.
FeQ< Scalar > JxW()
Same as BulkValues::JxW but register at side quadrature points.
FeQArray< Vector > vector_shape(uint component_idx=0)
Same as BulkValues::vector_shape but register at side quadrature points.
FeQArray< 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 > coords()
Create side accessor of coords entity.
friend class EdgePoint
FeQArray< 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.
RevertableList< EdgeIntegralData > & patch_data()
Getter for patch_data_.
RevertableList< EdgeIntegralData > patch_data_
Set of integral data of given dimension recomputed for each patch.
int get_subset_idx() const
Return index of data block according to subset in EvalPoints object.
EdgeIntegral()
Default constructor.
std::shared_ptr< internal_integrals::Edge > internal_edge_
Internal integral object.
EdgeIntegral(std::shared_ptr< EvalPoints > eval_points, Quadrature *quad, unsigned int dim)
Constructor of edge integral.
unsigned int n_sides() const
Getter of n_sides.
~EdgeIntegral()
Destructor.
Definition: integral_acc.cc:47
Range< EdgePoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const
Returns range of side local points for appropriate cell side accessor - obsolete method.
friend class EdgePoint
uint side_begin(const DHCellSide &cell_side) const
Point accessor allow iterate over quadrature points of given side defined in local element coordinate...
unsigned int idx() const
We need this method after replacing Region by RegionIdx, and movinf RegionDB instance into particular...
Definition: accessors.hh:215
Directing class of FieldValueCache.
unsigned int position_in_cache(unsigned mesh_elm_idx, bool bdr=false) const
Return position of element stored in ElementCacheMap.
Class holds local coordinations of evaluating points (bulk and sides) specified by element dimension.
Definition: eval_points.hh:54
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Evaluates Jacobian determinants on Bulk (Element) / Side.
Definition: op_function.hh:130
Evaluates normal vector on side quadrature points.
Definition: op_function.hh:446
PatchOp< spacedim > * get_for_elem_quad()
Returns operation of given dim and OpType, creates it if doesn't exist.
PatchOp< spacedim > * get(const Quadrature *quad)
Returns operation of given dim and OpType, creates it if doesn't exist.
void set_used_domain(fem_domain domain)
Mark domain (bulk or side) as used in assembly class.
std::shared_ptr< FiniteElement< dim > > fe_comp(std::shared_ptr< FiniteElement< dim > > fe, uint component_idx)
Returnd FiniteElement of component_idx for FESystem or fe for other types.
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:86
unsigned int dim() const
Definition: quadrature.hh:72
Range helper class.
BaseIntegral(std::shared_ptr< EvalPoints > eval_points, unsigned int dim)
Constructor of bulk or side subset.
Definition: integral_acc.hh:56
std::shared_ptr< EvalPoints > eval_points() const
Getter of eval_points.
Definition: integral_acc.hh:63
unsigned int dim_
Dimension of the cell on which points are placed.
Definition: integral_acc.hh:75
BaseIntegral()
Default constructor.
Definition: integral_acc.hh:53
std::shared_ptr< EvalPoints > eval_points_
Pointer to EvalPoints.
Definition: integral_acc.hh:73
unsigned int dim() const
Returns dimension.
Definition: integral_acc.hh:68
virtual ~BaseIntegral()
Destructor.
Definition: integral_acc.cc:29
FeQArray< ValueType > make_qarray(uint component_idx=0)
Factory method. Same as previous but creates FE operation.
ElementCacheMap * element_cache_map_
PatchOp< 3 > * make_patch_op()
Factory method. Creates operation of given OpType.
Definition: integral_acc.hh:96
IntegralFactory(PatchFEValues< 3 > *pfev, ElementCacheMap *element_cache_map, std::shared_ptr< FiniteElement< dim > > fe, Quadrature *quad)
Definition: integral_acc.hh:90
PatchOp< 3 > * make_elem_patch_op()
Factory method. Creates element / side operation of given OpType.
std::shared_ptr< FiniteElement< dim > > fe_
PatchFEValues< 3 > * patch_fe_values_
Quadrature * quad_
Pointer to Quadrature that represents quadrature points of integral.
unsigned int dim_
Dimension of the cell on which points are placed.
unsigned int dim() const
Returns dimension.
Base()
Default constructor.
virtual ~Base()
Destructor.
bool operator==(const Base &other) const
Comparison operator.
Quadrature * quad() const
Returns quadrature.
Base(Quadrature *quad, unsigned int dim)
Constructor of bulk or side subset.
unsigned int subset_index_
Index of data block according to subset in EvalPoints object.
Bulk(Quadrature *quad, unsigned int dim, std::shared_ptr< EvalPoints > eval_points, unsigned int subset_idx)
Constructor of bulk integral- obsolete constructor.
uint begin_idx_
Begin index of quadrature points in EvalPoinnts.
uint end_idx_
Begin index of quadrature points in EvalPoinnts.
Bulk()
Default constructor.
uint begin_idx() const
Getter of bulk_begin.
uint n_points_per_side_
Number of points. TODO: pass this to the constructor, avoid extraction from the eval_points.
Edge()
Default constructor.
uint side_begin(const DHCellSide &cell_side) const
unsigned int n_sides_
Number of sides (value 0 indicates bulk set)
Edge(Quadrature *quad, unsigned int dim, std::shared_ptr< EvalPoints > eval_points, unsigned int subset_idx)
Constructor of edge integral.
unsigned int uint
Store finite element reinit functions.
Class FEValues calculates finite element data on the actual cells such as shape function values,...
@ bulk_domain
@ side_domain
#define FMT_UNUSED
Definition: posix.h:75
Implementation of range helper class.