Hermes2D  2.0
space.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_SPACE_H
17 #define __H2D_SPACE_H
18 
19 #include "../mesh/mesh.h"
20 #include "../shapeset/shapeset.h"
21 #include "asmlist.h"
22 #include "../mesh/traverse.h"
23 #include "../quadrature/quad_all.h"
24 #include "../boundary_conditions/essential_boundary_conditions.h"
25 
26 using namespace Hermes::Algebra::DenseMatrixOperations;
27 
28 namespace Hermes
29 {
30  namespace Hermes2D
31  {
32  template<typename Scalar> class Adapt;
33  template<typename Scalar> class DiscreteProblem;
34  namespace Views
35  {
36  template<typename Scalar> class BaseView;
37  template<typename Scalar> class VectorBaseView;
38  class Orderizer;
39  class OrderView;
40  };
41  class Shapeset;
42 
43  template<typename Scalar> class L2Space;
44  template<typename Scalar> class H1Space;
45  template<typename Scalar> class HcurlSpace;
46  template<typename Scalar> class HdivSpace;
47 
52 
53 
54 
55 
56 
57 
58 
59 
60 
61 
62 
63 
64 
65 
66 
67 
68 
69 
70 
71 
72 
73 
74 
75 
76 
77 
78 
79 
80 
81 
82 
83 
84 
85 
86 
87 
88 
89 
90 
91 
92 
93 
94 
95 
96 
97 
98 
99 
100 
101 
102 
103 
104  template<typename Scalar>
105  class HERMES_API Space : public Hermes::Mixins::Loggable, public Hermes::Hermes2D::Mixins::StateQueryable, public Hermes::Hermes2D::Mixins::XMLParsing
106  {
107  public:
108  Space();
109  Space(const Mesh* mesh, Shapeset* shapeset, EssentialBCs<Scalar>* essential_bcs);
110 
112  void init();
113 
115  virtual bool isOkay() const;
116  inline std::string getClassName() const { return "Space"; }
117 
119  virtual ~Space();
120 
123  virtual void set_element_order(int id, int order);
124 
126  virtual void set_element_orders(int* elem_orders);
127 
129  int get_element_order(int id) const;
130 
133  void set_uniform_order(int order, std::string marker = HERMES_ANY);
134 
137  void adjust_element_order(int order_change, int min_order);
138 
140  void adjust_element_order(int horizontal_order_change, int vertical_order_change, unsigned int horizontal_min_order, unsigned int vertical_min_order);
141 
144  void unrefine_all_mesh_elements(bool keep_initial_refinements = true);
145 
147  void update_element_orders_after_refinement();
148 
150  virtual void set_shapeset(Shapeset* shapeset) = 0;
151 
153  int get_num_dofs() const;
154 
156  static int get_num_dofs(Hermes::vector<const Space<Scalar>*> spaces);
157  static int get_num_dofs(Hermes::vector<Space<Scalar>*> spaces);
158 
160  static int get_num_dofs(const Space<Scalar>* space);
161  static int get_num_dofs(Space<Scalar>* space);
162 
163  Mesh* get_mesh() const;
164 
166  void set_mesh(Mesh* mesh);
167 
169  void set_mesh_seq(int seq);
170 
172  void set_essential_bcs(EssentialBCs<Scalar>* essential_bcs);
173 
175  EssentialBCs<Scalar>* get_essential_bcs() const;
176 
179  void update_essential_bc_values();
180 
181  Shapeset* get_shapeset() const;
182 
184  bool save(const char *filename) const;
185 
187  static Space<Scalar>* load(const char *filename, Mesh* mesh, bool validate, EssentialBCs<Scalar>* essential_bcs = NULL, Shapeset* shapeset = NULL);
188 
190  virtual void get_element_assembly_list(Element* e, AsmList<Scalar>* al, unsigned int first_dof = 0) const;
191 
193  virtual void copy(const Space<Scalar>* space, Mesh* new_mesh);
194 
196  class HERMES_API ReferenceSpaceCreator
197  {
198  public:
203  ReferenceSpaceCreator(const Space<Scalar>* coarse_space, const Mesh* ref_mesh, unsigned int order_increase = 1);
204 
207  virtual void handle_orders(Space<Scalar>* ref_space);
208 
210  virtual Space<Scalar>* create_ref_space(bool assign_dofs = true);
211 
213  private:
214  L2Space<Scalar>* init_construction_l2();
215  H1Space<Scalar>* init_construction_h1();
216  HcurlSpace<Scalar>* init_construction_hcurl();
217  HdivSpace<Scalar>* init_construction_hdiv();
218 
220  virtual void finish_construction(Space<Scalar>* ref_space);
221 
223  const Space<Scalar>* coarse_space;
224  const Mesh* ref_mesh;
225  unsigned int order_increase;
226  };
227 
230  virtual void set_element_order_internal(int id, int order);
231 
238  virtual int assign_dofs(int first_dof = 0, int stride = 1);
239 
241  static int assign_dofs(Hermes::vector<Space<Scalar>*> spaces);
242 
243  virtual Scalar* get_bc_projection(SurfPos* surf_pos, int order, EssentialBoundaryCondition<Scalar> *bc) = 0;
244 
245  static void update_essential_bc_values(Hermes::vector<Space<Scalar>*> spaces, double time);
246 
247  static void update_essential_bc_values(Space<Scalar>*s, double time);
248 
251  virtual SpaceType get_type() const = 0;
252 
253  static Node* get_mid_edge_vertex_node(Element* e, int i, int j);
254 
256  void distribute_orders(Mesh* mesh, int* parents);
257 
259  virtual int get_edge_order(Element* e, int edge) const;
260 
262  int get_max_dof() const;
263 
265  bool is_up_to_date() const;
266 
268  void get_boundary_assembly_list(Element* e, int surf_num, AsmList<Scalar>* al, unsigned int first_dof = 0) const;
269 
272  void set_uniform_order_internal(int order, int marker);
273 
274  void free();
275 
278  int get_vertex_functions_count();
280  int get_edge_functions_count();
282  int get_bubble_functions_count();
283 
284  protected:
286  int ndof;
287 
288  static const int H2D_UNASSIGNED_DOF = -2;
289  static const int H2D_CONSTRAINED_DOF = -1;
290 
291  Shapeset* shapeset;
292 
294 
297 
299  const Mesh* mesh;
300 
301  int default_tri_order, default_quad_order;
302  int vertex_functions_count, edge_functions_count, bubble_functions_count;
303  int first_dof, next_dof;
304  int stride;
305  int seq, mesh_seq;
306  int was_assigned;
307 
309  {
310  int dof;
311  Scalar coef;
312  };
313 
314  union NodeData
315  {
316  struct // regular node
317  {
318  int dof;
319  union {
320  Scalar* edge_bc_proj;
321  Scalar* vertex_bc_coef;
322  };
323  int n;
324 
325  };
326  struct // constrained vertex node
327  {
328  BaseComponent* baselist;
329  int ncomponents;
330  };
331  struct // constrained edge node
332  {
333  Node* base;
334  int part;
335  };
336  NodeData() : dof(0), edge_bc_proj(NULL) {}
337  };
338 
340  {
341  public:
342  ElementData() : changed_in_last_adaptation(true) {};
343  int order;
344  int bdof, n;
345  bool changed_in_last_adaptation;
346  };
347 
350  int nsize, ndata_allocated;
351  int esize;
352 
353  virtual int get_edge_order_internal(Node* en) const;
354 
362  virtual void resize_tables();
363 
364  void update_orders_recurrent(Element* e, int order);
365 
366  virtual void reset_dof_assignment();
367  virtual void assign_vertex_dofs() = 0;
368  virtual void assign_edge_dofs() = 0;
369  virtual void assign_bubble_dofs() = 0;
370 
371  virtual void get_vertex_assembly_list(Element* e, int iv, AsmList<Scalar>* al) const = 0;
372  virtual void get_boundary_assembly_list_internal(Element* e, int surf_num, AsmList<Scalar>* al) const = 0;
373  virtual void get_bubble_assembly_list(Element* e, AsmList<Scalar>* al) const;
374 
375  double** proj_mat;
376  double* chol_p;
377 
379  Hermes::vector<void*> bc_data;
380 
381  void precalculate_projection_matrix(int nv, double**& mat, double*& p);
382  void update_edge_bc(Element* e, SurfPos* surf_pos);
383 
387  virtual void update_constraints();
388 
391  virtual void post_assign();
392 
393  void free_bc_data();
394 
396  int get_seq() const;
397  template<typename T> friend class OGProjection;
398  template<typename T> friend class NewtonSolver;
399  template<typename T> friend class PicardSolver;
400  template<typename T> friend class LinearSolver;
401  template<typename T> friend class OGProjectionNOX;
402  template<typename T> friend class LocalProjection;
403  template<typename T> friend class Solution;
404  template<typename T> friend class RungeKutta;
405  template<typename T> friend class ExactSolution;
406  template<typename T> friend class NeighborSearch;
407  template<typename T> friend class ExactSolutionScalar;
408  template<typename T> friend class ExactSolutionVector;
409  template<typename T> friend class Views::BaseView;
410  friend class Views::Orderizer;
411  friend class Views::OrderView;
412  template<typename T> friend class Views::VectorBaseView;
413  friend class Adapt<Scalar>;
414  friend class DiscreteProblem<Scalar>;
415  template<typename T> friend class CalculationContinuity;
416  };
417  }
418 }
419 #endif