Flow123d  JS_constraints-e651b99
patch_op.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_op.hh
15  * @brief Base class of FE operations.
16  * @author David Flanderka
17  */
18 
19 #ifndef PATCH_OP_HH_
20 #define PATCH_OP_HH_
21 
22 #include <Eigen/Dense>
23 
24 #include "fem/eigen_tools.hh"
25 #include "fem/arena_resource.hh"
26 #include "fem/arena_vec.hh"
27 #include "mesh/ref_element.hh"
28 
29 
30 template<unsigned int spacedim> class PatchFEValues;
31 
32 
33 /// Distinguishes bulk and side domain
35 {
38 };
39 
40 
41 /**
42  * @brief Class represents element or FE operations.
43  */
44 template<unsigned int spacedim = 3>
45 class PatchOp {
46 public:
47  /**
48  * Constructor
49  *
50  * Set all data members.
51  */
52  PatchOp(uint dim, PatchFEValues<spacedim> &pfev, const Quadrature *quad, std::initializer_list<uint> shape, uint n_dofs = 1)
53  : dim_(dim), shape_(set_shape_vec(shape)), n_dofs_(n_dofs), patch_fe_(&pfev), quad_(quad)
54  {
55  ASSERT_GT(n_dofs, 0);
56  result_ = Eigen::Vector<ArenaVec<double>, Eigen::Dynamic>(shape_[0] * shape_[1] * n_dofs_);
57  }
58 
59  /// Destructor
60  virtual ~PatchOp()
61  {}
62 
63 
64  /// Aligns shape_vec to 2 items (equal to matrix number of dimensions)
65  std::vector<uint> set_shape_vec(std::initializer_list<uint> shape) const {
66  std::vector<uint> shape_vec(shape);
67  if (shape_vec.size() == 1) shape_vec.push_back(1);
68  ASSERT_EQ(shape_vec.size(), 2);
69  return shape_vec;
70  }
71 
72  /**
73  * Return number of operation components
74  *
75  * Value is computed from shape_ vector
76  */
77  inline uint n_comp() const {
78  return shape_[0] * shape_[1];
79  }
80 
81 // /// Getter for dimension
82 // inline uint dimension() const { // dim is in conflict with template of some descendants
83 // return dim_;
84 // }
85 
86  /// Getter for bulk_side flag
87  inline ElemDomain domain() const {
88  return domain_;
89  }
90 
91  /// Getter for n_dofs_
92  inline uint n_dofs() const {
93  return n_dofs_;
94  }
95 
96  /// Return pointer to operation of i_op index in input operation vector.
97  inline PatchOp<spacedim> *input_ops(uint i_op) const {
98  return input_ops_[i_op];
99  }
100 
101  /// Getter for shape_
102  inline const std::vector<uint> &shape() const {
103  return shape_;
104  }
105 
106  /**
107  * Format shape to string
108  *
109  * Method is used in output development method.
110  */
111  std::string format_shape() const {
112  std::stringstream ss;
113  ss << shape_[0] << "x" << shape_[1];
114  return ss.str();
115  }
116 
117  inline void allocate_result(size_t data_size, PatchArena &arena) {
118  for (uint i=0; i<n_comp()*n_dofs_; ++i)
119  result_(i) = ArenaVec<double>(data_size, arena);
120  }
121 
122  inline void allocate_const_result(ArenaVec<double> &value_vec) {
123  for (uint i=0; i<n_comp()*n_dofs_; ++i)
124  result_(i) = value_vec;
125  }
126 
127  /// Return map referenced result as Eigen::Vector
128  inline Eigen::Map<Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic>> result_matrix() {
129  return Eigen::Map<Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic>>(result_.data(), shape_[0], shape_[1] * n_dofs_);
130  }
131 
132  /// Return map referenced result of DOF values as Eigen::Matrix
133  inline Eigen::Map<Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic>> result_sub_matrix(uint i_dof) {
134  ASSERT_LT(i_dof, n_dofs_);
135  uint n_dof_comps = shape_[0] * shape_[1];
136  return Eigen::Map<Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic>>(result_.data()+i_dof*n_dof_comps, shape_[0], shape_[1]);
137  }
138 
139  /// Return map referenced result as Eigen::Vector
140  inline Eigen::Vector<ArenaVec<double>, Eigen::Dynamic> &raw_result() {
141  return result_;
142  }
143 
144  /// Same as previous but return const reference
145  inline const Eigen::Vector<ArenaVec<double>, Eigen::Dynamic> &raw_result() const {
146  return result_;
147  }
148 
149  /// Return reference of PatchPointValues
150  inline PatchPointValues<spacedim> &ppv() const {
151  return patch_fe_->patch_point_vals_[domain_][this->dim_-1];
152  }
153 
154  /// Reinit function of operation. Implementation in descendants.
155  virtual void eval() =0;
156 
157  /**
158  * Returns output value of data stored by elements.
159  *
160  * @param point_idx Index of quadrature point in ElementCacheMap
161  */
162  template <class ValueType>
163  inline ValueType elem_value(uint point_idx) const;
164 
165  /**
166  * Returns output value on quadrature point.
167  *
168  * @param point_idx Index of quadrature point in ElementCacheMap
169  * @param i_dof Index of DOF
170  */
171  template <class ValueType>
172  inline ValueType point_value(uint point_idx, uint i_dof=0) const;
173 
174  /// Return size of quadrature
175  inline unsigned int quad_size() const {
176  return quad_->size();
177  }
178 
179  /// return reference to patch arena
180  inline PatchArena &patch_arena() const {
181  return patch_fe_->patch_arena();
182  }
183 
184 
185 protected:
186  uint dim_; ///< Dimension
187  ElemDomain domain_; ///< Flag: BulkOp = 0, SideOp = 1
188  std::vector<uint> shape_; ///< Shape of stored data (size of vector or number of rows and cols of matrix)
189  Eigen::Vector<ArenaVec<double>, Eigen::Dynamic> result_; ///< Result matrix of operation
190  std::vector<PatchOp<spacedim> *> input_ops_; ///< Indices of operations in PatchPointValues::operations_ vector on which PatchOp is depended
191  uint n_dofs_; ///< Number of DOFs of FE operations (or 1 in case of element operations)
192  PatchFEValues<spacedim> *patch_fe_; ///< Pointer to PatchFEValues object
193  const Quadrature* quad_; ///< Pointer to Quadrature
194 
195  friend class PatchFEValues<spacedim>;
196 };
197 
198 
199 #endif /* PATCH_OP_HH_ */
#define ASSERT_LT(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:301
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:333
#define ASSERT_GT(a, b)
Definition of comparative assert macro (Greater Than) only for debug mode.
Definition: asserts.hh:317
Class represents element or FE operations.
Definition: patch_op.hh:45
Eigen::Vector< ArenaVec< double >, Eigen::Dynamic > result_
Result matrix of operation.
Definition: patch_op.hh:189
Eigen::Map< Eigen::Matrix< ArenaVec< double >, Eigen::Dynamic, Eigen::Dynamic > > result_sub_matrix(uint i_dof)
Return map referenced result of DOF values as Eigen::Matrix.
Definition: patch_op.hh:133
const Quadrature * quad_
Pointer to Quadrature.
Definition: patch_op.hh:193
uint n_dofs() const
Getter for n_dofs_.
Definition: patch_op.hh:92
PatchOp(uint dim, PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::initializer_list< uint > shape, uint n_dofs=1)
Definition: patch_op.hh:52
uint n_comp() const
Definition: patch_op.hh:77
virtual ~PatchOp()
Destructor.
Definition: patch_op.hh:60
std::vector< PatchOp< spacedim > * > input_ops_
Indices of operations in PatchPointValues::operations_ vector on which PatchOp is depended.
Definition: patch_op.hh:190
PatchOp< spacedim > * input_ops(uint i_op) const
Return pointer to operation of i_op index in input operation vector.
Definition: patch_op.hh:97
ElemDomain domain() const
Getter for bulk_side flag.
Definition: patch_op.hh:87
Eigen::Vector< ArenaVec< double >, Eigen::Dynamic > & raw_result()
Return map referenced result as Eigen::Vector.
Definition: patch_op.hh:140
uint dim_
Dimension.
Definition: patch_op.hh:186
const std::vector< uint > & shape() const
Getter for shape_.
Definition: patch_op.hh:102
ValueType point_value(uint point_idx, uint i_dof=0) const
std::string format_shape() const
Definition: patch_op.hh:111
std::vector< uint > set_shape_vec(std::initializer_list< uint > shape) const
Aligns shape_vec to 2 items (equal to matrix number of dimensions)
Definition: patch_op.hh:65
Eigen::Map< Eigen::Matrix< ArenaVec< double >, Eigen::Dynamic, Eigen::Dynamic > > result_matrix()
Return map referenced result as Eigen::Vector.
Definition: patch_op.hh:128
std::vector< uint > shape_
Shape of stored data (size of vector or number of rows and cols of matrix)
Definition: patch_op.hh:188
void allocate_const_result(ArenaVec< double > &value_vec)
Definition: patch_op.hh:122
ElemDomain domain_
Flag: BulkOp = 0, SideOp = 1.
Definition: patch_op.hh:187
void allocate_result(size_t data_size, PatchArena &arena)
Definition: patch_op.hh:117
const Eigen::Vector< ArenaVec< double >, Eigen::Dynamic > & raw_result() const
Same as previous but return const reference.
Definition: patch_op.hh:145
unsigned int quad_size() const
Return size of quadrature.
Definition: patch_op.hh:175
PatchPointValues< spacedim > & ppv() const
Return reference of PatchPointValues.
Definition: patch_op.hh:150
PatchArena & patch_arena() const
return reference to patch arena
Definition: patch_op.hh:180
virtual void eval()=0
Reinit function of operation. Implementation in descendants.
uint n_dofs_
Number of DOFs of FE operations (or 1 in case of element operations)
Definition: patch_op.hh:191
ValueType elem_value(uint point_idx) const
PatchFEValues< spacedim > * patch_fe_
Pointer to PatchFEValues object.
Definition: patch_op.hh:192
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
Store finite element data on the actual patch such as shape function values, gradients,...
unsigned int uint
ElemDomain
Distinguishes bulk and side domain.
Definition: patch_op.hh:35
@ domain_side
Definition: patch_op.hh:37
@ domain_bulk
Definition: patch_op.hh:36
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc.