Hermes2D  3.0
function.h
1 // This file is part of Hermes2D.
2 //
3 // Hermes2D is free software: you can redistribute it and/or modify
4 // it under the terms of the GNU General Public License as published by
5 // the Free Software Foundation, either version 2 of the License, or
6 // (at your option) any later version.
7 //
8 // Hermes2D is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11 // GNU General Public License for more details.
12 //
13 // You should have received a copy of the GNU General Public License
14 // along with Hermes2D. If not, see <http://www.gnu.org/licenses/>.
15 
16 #ifndef __H2D_FUNCTION_H
17 #define __H2D_FUNCTION_H
18 
19 #include "transformable.h"
20 #include "../quadrature/quad.h"
21 #include "exceptions.h"
22 namespace Hermes
23 {
24  namespace Hermes2D
25  {
26  struct SurfPos;
27  class PrecalcShapeset;
28  namespace RefinementSelectors{
29  template<typename Scalar> class Selector;
30  template<typename Scalar> class HOnlySelector;
31  template<typename Scalar> class POnlySelector;
32  template<typename Scalar> class OptimumSelector;
33  template<typename Scalar> class ProjBasedSelector;
34  template<typename Scalar> class L2ProjBasedSelector;
35  template<typename Scalar> class H1ProjBasedSelector;
36  template<typename Scalar> class HcurlProjBasedSelector;
37  };
38 
40  enum
41  {
42  H2D_FN_VAL_0 = 0x0001, H2D_FN_VAL_1 = 0x0040, // Function values
43  H2D_FN_DX_0 = 0x0002, H2D_FN_DX_1 = 0x0080, // First derivative
44  H2D_FN_DY_0 = 0x0004, H2D_FN_DY_1 = 0x0100, // First derivative
45 #ifdef H2D_USE_SECOND_DERIVATIVES
46  H2D_FN_DXX_0 = 0x0008, H2D_FN_DXX_1 = 0x0200, // Second derivative
47  H2D_FN_DYY_0 = 0x0010, H2D_FN_DYY_1 = 0x0400, // Second derivative
48  H2D_FN_DXY_0 = 0x0020, H2D_FN_DXY_1 = 0x0800 // Second mixed derivative
49 #endif
50  };
51 
53  const int H2D_FN_VAL = H2D_FN_VAL_0 | H2D_FN_VAL_1;
54  const int H2D_FN_DX = H2D_FN_DX_0 | H2D_FN_DX_1;
55  const int H2D_FN_DY = H2D_FN_DY_0 | H2D_FN_DY_1;
57  const int H2D_FN_DEFAULT = H2D_FN_VAL | H2D_FN_DX | H2D_FN_DY;
58 
59 #ifdef H2D_USE_SECOND_DERIVATIVES
60  const int H2D_FN_COMPONENT_0 = H2D_FN_VAL_0 | H2D_FN_DX_0 | H2D_FN_DY_0 | H2D_FN_DXX_0 | H2D_FN_DYY_0 | H2D_FN_DXY_0;
61  const int H2D_FN_COMPONENT_1 = H2D_FN_VAL_1 | H2D_FN_DX_1 | H2D_FN_DY_1 | H2D_FN_DXX_1 | H2D_FN_DYY_1 | H2D_FN_DXY_1;
62 #else
63  const int H2D_FN_COMPONENT_0 = H2D_FN_VAL_0 | H2D_FN_DX_0 | H2D_FN_DY_0;
64  const int H2D_FN_COMPONENT_1 = H2D_FN_VAL_1 | H2D_FN_DX_1 | H2D_FN_DY_1;
65 #endif
66 
67  const int H2D_GRAD = H2D_FN_DX_0 | H2D_FN_DY_0;
68  const int H2D_CURL = H2D_FN_DX | H2D_FN_DY;
69 
70 #ifdef H2D_USE_SECOND_DERIVATIVES
71  const int H2D_FN_DXX = H2D_FN_DXX_0 | H2D_FN_DXX_1;
72  const int H2D_FN_DYY = H2D_FN_DYY_0 | H2D_FN_DYY_1;
73  const int H2D_FN_DXY = H2D_FN_DXY_0 | H2D_FN_DXY_1;
74 
75  const int H2D_SECOND = H2D_FN_DXX_0 | H2D_FN_DXY_0 | H2D_FN_DYY_0;
76 
78  const int H2D_FN_ALL = H2D_FN_DEFAULT | H2D_FN_DXX | H2D_FN_DYY | H2D_FN_DXY;
79 #endif
80 
105  template<typename Scalar>
106  class HERMES_API Function : public Transformable
107  {
108  public:
109 
111  Function();
112 
115  virtual ~Function();
116 
118  unsigned char get_num_components() const;
119 
123  virtual const Scalar* get_fn_values(int component = 0) const;
124 
128  virtual const Scalar* get_dx_values(int component = 0) const;
129 
133  virtual const Scalar* get_dy_values(int component = 0) const;
134 
135 #ifdef H2D_USE_SECOND_DERIVATIVES
136  virtual const Scalar* get_dxx_values(int component = 0) const;
140 
144  virtual const Scalar* get_dyy_values(int component = 0) const;
145 
149  virtual const Scalar* get_dxy_values(int component = 0) const;
150 #endif
151 
155  Scalar* deep_copy_array(int component = 0, int item = 0) const;
156 
158  Quad2D* get_quad_2d() const;
159 
164  void set_quad_order(unsigned short order, unsigned short mask = H2D_FN_DEFAULT);
165 
166  virtual const Scalar* get_values(int component, int item) const;
167 
169  virtual int get_fn_order() const;
170 
173  virtual void push_transform(int son);
174 
177  virtual void pop_transform();
178 
180  virtual void set_active_element(Element* e);
181 
184  virtual void set_transform(uint64_t idx);
185 
186  protected:
192  virtual void set_quad_2d(Quad2D* quad_2d);
193 
195  virtual void reset_transform();
197  virtual void force_transform(uint64_t sub_idx, Trf* ctm);
198 
200  virtual void free() = 0;
201 
203  Scalar values[H2D_MAX_SOLUTION_COMPONENTS][H2D_NUM_FUNCTION_VALUES][H2D_MAX_INTEGRATION_POINTS_COUNT];
206 
209  virtual int get_edge_fn_order(unsigned char edge) const;
210 
212  virtual void precalculate(unsigned short order, unsigned short mask);
213 
215  int order;
216 
218  unsigned char num_components;
219 
223  void invalidate_values();
224 
226  Quad2D* quads[H2D_MAX_QUADRATURES];
228  int cur_quad;
229 
231  static int idx2mask[H2D_NUM_FUNCTION_VALUES][2];
232 
233  template<typename T> friend class KellyTypeAdapt;
234  template<typename T> friend class RefinementSelectors::H1ProjBasedSelector;
235  template<typename T> friend class RefinementSelectors::L2ProjBasedSelector;
236  template<typename T> friend class RefinementSelectors::HcurlProjBasedSelector;
237  template<typename T> friend class Adapt;
238 
239  template<typename T> friend class DiscontinuousFunc;
240  template<typename T> friend class DiscreteProblem;
241  friend class CurvMap;
242 
243  template<typename T> friend class Func;
244  template<typename T> friend class Filter;
245  template<typename T> friend class SimpleFilter;
246  template<typename T> friend class DXDYFilter;
247  friend class ComplexFilter;
248  friend class VonMisesFilter;
249  friend HERMES_API Func<double>* init_fn(PrecalcShapeset *fu, RefMap *rm, const int order);
250  template<typename T> friend HERMES_API Func<T>* init_fn(MeshFunction<T>*fu, const int order);
251  };
252  }
253 }
254 #endif
Definition: adapt.h:24
A projection-based selector for Hcurl space.
Definition: function.h:36
Calculates the Von Mises stress.
Definition: filter.h:343
Stores one element of a mesh.
Definition: element.h:107
Caches precalculated shape function values.
Definition: precalc.h:32
This class represents a function with jump discontinuity on an interface of two elements.
Definition: forms.h:335
Represents a function defined on a mesh.
Definition: mesh_function.h:56
A projection-based selector for H1 space.
Definition: function.h:35
A selector that selects H-refinements only.
Definition: function.h:30
HERMES_API Func< double > * init_fn(PrecalcShapeset *fu, RefMap *rm, const int order)
Init the shape function for the evaluation of the volumetric/surface integral (transformation of valu...
Definition: forms.cpp:469
A general projection-based selector.
Definition: function.h:33
const int H2D_FN_DEFAULT
default precalculation mask
Definition: function.h:57
2D transformation.
Definition: transformable.h:29
Represents an arbitrary function defined on an element.
Definition: function.h:106
A selector that increases order (i.e., it selects P-refinements only).
Definition: function.h:31
bool values_valid
Flag that the data are not 'dirty'.
Definition: function.h:205
Calculated function values (from the class Function) on an element for assembling.
Definition: forms.h:214
const int H2D_FN_VAL
Both components are usually requested together...
Definition: function.h:53
A selector that chooses an optimal candidates based on a score.
Definition: function.h:32
Represents the reference mapping.
Definition: refmap.h:40
unsigned char num_components
Number of vector components.
Definition: function.h:218
A parent of all refinement selectors. Abstract class.
Definition: function.h:29
A projection-based selector for L2 space.
Definition: function.h:34
int cur_quad
Active quadrature (index into 'quads')
Definition: function.h:228
int order
Current function polynomial order.
Definition: function.h:215