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