Hermes2D  2.0
space_hdiv.cpp
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 #include "global.h"
17 #include "space_hdiv.h"
18 #include "matrix.h"
19 #include "quad_all.h"
20 #include "shapeset/shapeset_hd_all.h"
21 #include "boundary_conditions/essential_boundary_conditions.h"
22 namespace Hermes
23 {
24  namespace Hermes2D
25  {
26  template<typename Scalar>
27  HdivSpace<Scalar>::HdivSpace() : Space<Scalar>()
28  {
29  }
30 
31  template<typename Scalar>
32  void HdivSpace<Scalar>::init(Shapeset* shapeset, int p_init)
33  {
34  if(shapeset == NULL)
35  {
36  this->shapeset = new HdivShapeset;
37  this->own_shapeset = true;
38  }
39  if(this->shapeset->get_num_components() < 2) throw Hermes::Exceptions::Exception("HdivSpace requires a vector shapeset.");
40 
41  this->precalculate_projection_matrix(0, this->proj_mat, this->chol_p);
42 
43  // set uniform poly order in elements
44  if(p_init < 0) throw Hermes::Exceptions::Exception("P_INIT must be >= 0 in an Hdiv space.");
45  else this->set_uniform_order_internal(p_init, HERMES_ANY_INT);
46 
47  // enumerate basis functions
48  this->assign_dofs();
49  }
50 
51  template<typename Scalar>
53  : Space<Scalar>(mesh, shapeset, essential_bcs)
54  {
55  init(shapeset, p_init);
56  }
57 
58  template<typename Scalar>
59  HdivSpace<Scalar>::HdivSpace(const Mesh* mesh, int p_init, Shapeset* shapeset)
60  : Space<Scalar>(mesh, shapeset, NULL)
61  {
62  init(shapeset, p_init);
63  }
64 
65  template<typename Scalar>
66  HdivSpace<Scalar>::~HdivSpace()
67  {
68  if(this->own_shapeset)
69  delete this->shapeset;
70  }
71 
72  template<typename Scalar>
73  void HdivSpace<Scalar>::copy(const Space<Scalar>* space, Mesh* new_mesh)
74  {
75  Space<Scalar>::copy(space, new_mesh);
76 
77  this->precalculate_projection_matrix(0, this->proj_mat, this->chol_p);
78  }
79 
80  template<typename Scalar>
82  {
83  if(shapeset->get_id() < 30 && shapeset->get_id() > 19)
84  {
85  this->shapeset = shapeset;
86  this->own_shapeset = false;
87  }
88  else
89  throw Hermes::Exceptions::Exception("Wrong shapeset type in HdivSpace<Scalar>::set_shapeset()");
90  }
91 
92  template<typename Scalar>
94  {
95  Node* en;
96  this->edge_functions_count = 0;
97  for_all_edge_nodes(en, this->mesh)
98  {
99  if(en->ref > 1 || en->bnd || this->mesh->peek_vertex_node(en->p1, en->p2) != NULL)
100  {
101  int ndofs = this->get_edge_order_internal(en) + 1;
102  this->ndata[en->id].n = ndofs;
103 
104  if(en->bnd)
105  if(this->essential_bcs != NULL)
106  if(this->essential_bcs->get_boundary_condition(this->mesh->boundary_markers_conversion.get_user_marker(en->marker).marker) != NULL)
107  this->ndata[en->id].dof = this->H2D_CONSTRAINED_DOF;
108  else
109  {
110  this->ndata[en->id].dof = this->next_dof;
111  this->next_dof += ndofs * this->stride;
112  this->edge_functions_count += ndofs;
113  }
114  else
115  {
116  this->ndata[en->id].dof = this->next_dof;
117  this->next_dof += ndofs * this->stride;
118  this->edge_functions_count += ndofs;
119  }
120  else
121  {
122  this->ndata[en->id].dof = this->next_dof;
123  this->next_dof += ndofs * this->stride;
124  this->edge_functions_count += ndofs;
125  }
126  }
127  else
128  {
129  this->ndata[en->id].n = -1;
130  }
131  }
132  }
133 
134  template<typename Scalar>
135  void HdivSpace<Scalar>::assign_bubble_dofs()
136  {
137  Element* e;
138  this->bubble_functions_count = 0;
139  for_all_active_elements(e, this->mesh)
140  {
141  typename Space<Scalar>::ElementData* ed = &this->edata[e->id];
142  ed->bdof = this->next_dof;
143  ed->n = this->shapeset->get_num_bubbles(ed->order, e->get_mode());
144  this->next_dof += ed->n * this->stride;
145  this->bubble_functions_count += ed->n;
146  }
147  }
148 
149  template<typename Scalar>
150  void HdivSpace<Scalar>::get_boundary_assembly_list_internal(Element* e, int surf_num, AsmList<Scalar>* al) const
151  {
152  Node* en = e->en[surf_num];
153  typename Space<Scalar>::NodeData* nd = &this->ndata[en->id];
154 
155  if(nd->n >= 0) // unconstrained
156  {
157  if(nd->dof >= 0)
158  {
159  int ori = (e->vn[surf_num]->id < e->vn[e->next_vert(surf_num)]->id) ? 0 : 1;
160  for (int j = 0, dof = nd->dof; j < nd->n; j++, dof += this->stride)
161  al->add_triplet(this->shapeset->get_edge_index(surf_num, ori, j, e->get_mode()), dof, 1.0);
162  }
163  else
164  {
165  for (int j = 0; j < nd->n; j++)
166  al->add_triplet(this->shapeset->get_edge_index(surf_num, 0, j, e->get_mode()), -1, nd->edge_bc_proj[j]);
167  }
168  }
169  else // constrained
170  {
171  int part = nd->part;
172  int ori = part < 0 ? 1 : 0;
173  if(part < 0) part ^= ~0;
174 
175  nd = &this->ndata[nd->base->id]; // ccc
176  for (int j = 0, dof = nd->dof; j < nd->n; j++, dof += this->stride)
177  al->add_triplet(this->shapeset->get_constrained_edge_index(surf_num, j, ori, part, e->get_mode()), dof, 1.0);
178  }
179  }
180 
181  template<typename Scalar>
182  void HdivSpace<Scalar>::get_bubble_assembly_list(Element* e, AsmList<Scalar>* al) const
183  {
184  typename Space<Scalar>::ElementData* ed = &this->edata[e->id];
185  if(!ed->n) return;
186 
187  int* indices = this->shapeset->get_bubble_indices(ed->order, e->get_mode());
188  for (int i = 0, dof = ed->bdof; i < ed->n; i++, dof += this->stride)
189  al->add_triplet(*indices++, dof, 1.0);
190  }
191 
193 
194  template<typename Scalar>
195  Scalar* HdivSpace<Scalar>::get_bc_projection(SurfPos* surf_pos, int order, EssentialBoundaryCondition<Scalar> *bc)
196  {
197  assert(order >= 0);
198  Scalar* proj = new Scalar[order + 1];
199 
200  Quad1DStd quad1d;
201  Scalar* rhs = proj;
202  int mo = quad1d.get_max_order();
203  double2* pt = quad1d.get_points(mo);
204 
205  Node* vn1 = this->mesh->get_node(surf_pos->v1);
206  Node* vn2 = this->mesh->get_node(surf_pos->v2);
207  double el = sqrt(sqr(vn1->x - vn2->x) + sqr(vn1->y - vn2->y));
208  el *= 0.5 * (surf_pos->hi - surf_pos->lo);
209 
210  // get boundary values at integration points, construct rhs
211  for (int i = 0; i <= order; i++)
212  {
213  rhs[i] = 0.0;
214  int ii = this->shapeset->get_edge_index(0, 0, i, surf_pos->base->get_mode());
215  for (int j = 0; j < quad1d.get_num_points(mo); j++)
216  {
217  double t = (pt[j][0] + 1) * 0.5, s = 1.0 - t;
218  surf_pos->t = surf_pos->lo * s + surf_pos->hi * t;
219 
220  // If the BC on this part of the boundary is constant.
221  EssentialBoundaryCondition<Scalar> *bc = this->essential_bcs->get_boundary_condition(this->mesh->boundary_markers_conversion.get_user_marker(surf_pos->marker).marker);
222 
223  if(bc->get_value_type() == EssentialBoundaryCondition<Scalar>::BC_CONST)
224  {
225  rhs[i] += pt[j][1] * this->shapeset->get_fn_value(ii, pt[j][0], -1.0, 1, surf_pos->base->get_mode())
226  * bc->value_const * el;
227  }
228  // If the BC is not constant.
229  else if(bc->get_value_type() == EssentialBoundaryCondition<Scalar>::BC_FUNCTION)
230  {
231  // Find out the (x, y) coordinate.
232  double x, y, n_x, n_y, t_x, t_y;
233  Nurbs* nurbs = surf_pos->base->is_curved() ? surf_pos->base->cm->nurbs[surf_pos->surf_num] : NULL;
234  CurvMap::nurbs_edge(surf_pos->base, nurbs, surf_pos->surf_num, 2.0*surf_pos->t - 1.0, x, y, n_x, n_y, t_x, t_y);
235  // Calculate.
236  rhs[i] += pt[j][1] * this->shapeset->get_fn_value(ii, pt[j][0], -1.0, 1, surf_pos->base->get_mode())
237  * bc->value(x, y, n_x, n_y, t_x, t_y) * el;
238  }
239  }
240  }
241 
242  // solve the system using a precalculated Cholesky decomposed projection matrix
243  cholsl(this->proj_mat, order + 1, this->chol_p, rhs, rhs);
244 
245  for (int i = 0; i < order + 1; i++)
246  proj[i] = 0.0;
247 
248  return proj;
249  }
250 
251  template<typename Scalar>
252  void HdivSpace<Scalar>::update_constrained_nodes(Element* e, EdgeInfo* ei0, EdgeInfo* ei1, EdgeInfo* ei2, EdgeInfo* ei3)
253  {
254  int j;
255  EdgeInfo* ei[4] = { ei0, ei1, ei2, ei3 };
256  typename Space<Scalar>::NodeData* nd;
257 
258  // on non-refined elements all we have to do is update edge nodes lying on constrained edges
259  if(e->active)
260  {
261  for (unsigned int i = 0; i < e->get_nvert(); i++)
262  {
263  if(ei[i] != NULL)
264  {
265  nd = &this->ndata[e->en[i]->id];
266  nd->base = ei[i]->node;
267  nd->part = ei[i]->part;
268  if(ei[i]->ori) nd->part ^= ~0;
269  }
270  }
271  }
272  // the element has sons - update mid-edge constrained vertex nodes
273  else
274  {
275  // create new edge infos where we don't have them yet
276  EdgeInfo ei_data[4];
277  for (unsigned int i = 0; i < e->get_nvert(); i++)
278  {
279  if(ei[i] == NULL)
280  {
281  j = e->next_vert(i);
282  Node* mid_vn = this->get_mid_edge_vertex_node(e, i, j);
283  if(mid_vn != NULL && mid_vn->is_constrained_vertex())
284  {
285  Node* mid_en = this->mesh->peek_edge_node(e->vn[i]->id, e->vn[j]->id);
286  if(mid_en != NULL)
287  {
288  ei[i] = ei_data + i;
289  ei[i]->node = mid_en;
290  ei[i]->part = -1;
291  ei[i]->lo = -1.0;
292  ei[i]->hi = 1.0;
293  ei[i]->ori = (e->vn[i]->id < e->vn[j]->id) ? 0 : 1;
294  }
295  }
296  }
297  }
298 
299  // create edge infos for half-edges
300  EdgeInfo half_ei_data[4][2];
301  EdgeInfo* half_ei[4][2];
302  for (unsigned int i = 0; i < e->get_nvert(); i++)
303  {
304  if(ei[i] == NULL)
305  {
306  half_ei[i][0] = half_ei[i][1] = NULL;
307  }
308  else
309  {
310  EdgeInfo* h0 = half_ei[i][0] = half_ei_data[i];
311  EdgeInfo* h1 = half_ei[i][1] = half_ei_data[i] + 1;
312 
313  h0->node = h1->node = ei[i]->node;
314  h0->part = (ei[i]->part + 1) * 2;
315  h1->part = h0->part + 1;
316  h0->hi = h1->lo = (ei[i]->lo + ei[i]->hi) / 2;
317  h0->lo = ei[i]->lo;
318  h1->hi = ei[i]->hi;
319  h1->ori = h0->ori = ei[i]->ori;
320  }
321  }
322 
323  // recur to sons
324  if(e->is_triangle())
325  {
326  update_constrained_nodes(e->sons[0], half_ei[0][0], NULL, half_ei[2][1], NULL);
327  update_constrained_nodes(e->sons[1], half_ei[0][1], half_ei[1][0], NULL, NULL);
328  update_constrained_nodes(e->sons[2], NULL, half_ei[1][1], half_ei[2][0], NULL);
329  update_constrained_nodes(e->sons[3], NULL, NULL, NULL, NULL);
330  }
331  else if(e->sons[2] == NULL) // 'horizontally' split quad
332  {
333  update_constrained_nodes(e->sons[0], ei[0], half_ei[1][0], NULL, half_ei[3][1]);
334  update_constrained_nodes(e->sons[1], NULL, half_ei[1][1], ei[2], half_ei[3][0]);
335  }
336  else if(e->sons[0] == NULL) // 'vertically' split quad
337  {
338  update_constrained_nodes(e->sons[2], half_ei[0][0], NULL, half_ei[2][1], ei[3]);
339  update_constrained_nodes(e->sons[3], half_ei[0][1], ei[1], half_ei[2][0], NULL);
340  }
341  else // fully split quad
342  {
343  update_constrained_nodes(e->sons[0], half_ei[0][0], NULL, NULL, half_ei[3][1]);
344  update_constrained_nodes(e->sons[1], half_ei[0][1], half_ei[1][0], NULL, NULL);
345  update_constrained_nodes(e->sons[2], NULL, half_ei[1][1], half_ei[2][0], NULL);
346  update_constrained_nodes(e->sons[3], NULL, NULL, half_ei[2][1], half_ei[3][0]);
347  }
348  }
349  }
350 
351  template<typename Scalar>
353  {
354  Element* e;
355  for_all_base_elements(e, this->mesh)
356  update_constrained_nodes(e, NULL, NULL, NULL, NULL);
357  }
358 
359  template HERMES_API class HdivSpace<double>;
360  template HERMES_API class HdivSpace<std::complex<double> >;
361  }
362 }