17 #include "space_hcurl.h"
20 #include "shapeset_hc_all.h"
21 #include "essential_boundary_conditions.h"
27 template<
typename Scalar>
28 HcurlSpace<Scalar>::HcurlSpace() : Space<Scalar>()
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)
36 init(shapeset, p_init);
39 template<
typename Scalar>
40 HcurlSpace<Scalar>::HcurlSpace(MeshSharedPtr mesh,
int p_init, Shapeset* shapeset)
41 : Space<Scalar>(mesh, shapeset, nullptr)
43 init(shapeset, p_init);
46 template<
typename Scalar>
49 if (shapeset ==
nullptr)
52 this->own_shapeset =
true;
54 if (this->shapeset->
get_num_components() < 2)
throw Hermes::Exceptions::Exception(
"HcurlSpace requires a vector shapeset.");
56 this->precalculate_projection_matrix(0, this->proj_mat, this->chol_p);
63 throw Hermes::Exceptions::Exception(
"P_INIT must be >= 0 in an Hcurl space.");
65 this->set_uniform_order_internal(p_init, HERMES_ANY_INT);
71 template<
typename Scalar>
76 template<
typename Scalar>
79 this->set_shapeset(space->get_shapeset(),
true);
81 this->precalculate_projection_matrix(0, this->proj_mat, this->chol_p);
86 template<
typename Scalar>
89 if (!(shapeset->
get_id() < 20 && shapeset->
get_id() > 9))
90 throw Hermes::Exceptions::Exception(
"Wrong shapeset type in HcurlSpace<Scalar>::set_shapeset()");
94 if (this->own_shapeset)
95 delete this->shapeset;
97 this->shapeset = shapeset->clone();
101 this->shapeset = shapeset;
102 this->own_shapeset =
false;
106 template<
typename Scalar>
110 this->edge_functions_count = 0;
111 for_all_edge_nodes(en, this->mesh)
113 if (en->
ref > 1 || en->
bnd || this->mesh->peek_vertex_node(en->
p1, en->p2) !=
nullptr)
115 int ndofs = this->get_edge_order_internal(en) + 1;
116 this->ndata[en->
id].n = ndofs;
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;
123 this->ndata[en->
id].dof = this->next_dof;
124 this->next_dof += ndofs;
125 this->edge_functions_count += ndofs;
129 this->ndata[en->
id].dof = this->next_dof;
130 this->next_dof += ndofs;
131 this->edge_functions_count += ndofs;
135 this->ndata[en->
id].dof = this->next_dof;
136 this->next_dof += ndofs;
137 this->edge_functions_count += ndofs;
141 this->ndata[en->
id].n = -1;
145 template<
typename Scalar>
146 void HcurlSpace<Scalar>::assign_bubble_dofs()
149 this->bubble_functions_count = 0;
150 for_all_active_elements(e, this->mesh)
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;
160 template<
typename Scalar>
161 void HcurlSpace<Scalar>::get_boundary_assembly_list_internal(Element* e,
int surf_num, AsmList<Scalar>* al)
const
163 Node* en = e->en[surf_num];
164 typename Space<Scalar>::NodeData* nd = &this->ndata[en->id];
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);
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]);
183 int ori = part < 0 ? 1 : 0;
184 if (part < 0) part ^= ~0;
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);
195 template<
typename Scalar>
196 Scalar* HcurlSpace<Scalar>::get_bc_projection(SurfPos* surf_pos,
int order, EssentialBoundaryCondition<Scalar> *bc)
199 Scalar* proj = malloc_with_check<HcurlSpace<Scalar>, Scalar>(order + 1,
this);
203 int mo = quad1d.get_max_order();
204 double2* pt = quad1d.get_points(mo);
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);
212 for (
int i = 0; i <= order; i++)
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++)
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;
221 if (bc->get_value_type() == BC_CONST)
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;
227 else if (bc->get_value_type() == BC_FUNCTION)
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);
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;
241 cholsl(this->proj_mat, order + 1, this->chol_p, rhs, rhs);
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)
249 typename Space<Scalar>::EdgeInfo* ei[4] = { ei0, ei1, ei2, ei3 };
250 typename Space<Scalar>::NodeData* nd;
255 for (
unsigned char i = 0; i < e->get_nvert(); i++)
257 if (ei[i] !=
nullptr)
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;
270 typename Space<Scalar>::EdgeInfo ei_data[4];
271 for (
unsigned char i = 0; i < e->get_nvert(); i++)
273 if (ei[i] ==
nullptr)
276 Node* mid_vn = this->get_mid_edge_vertex_node(e, i, j);
277 if (mid_vn !=
nullptr && mid_vn->is_constrained_vertex())
279 Node* mid_en = this->mesh->peek_edge_node(e->vn[i]->id, e->vn[j]->id);
280 if (mid_en !=
nullptr)
283 ei[i]->node = mid_en;
287 ei[i]->ori = (e->vn[i]->id < e->vn[j]->id) ? 0 : 1;
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++)
298 if (ei[i] ==
nullptr)
300 half_ei[i][0] = half_ei[i][1] =
nullptr;
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;
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;
313 h1->ori = h0->ori = ei[i]->ori;
318 if (e->is_triangle())
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);
325 else if (e->sons[2] ==
nullptr)
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]);
330 else if (e->sons[0] ==
nullptr)
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);
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]);
345 template<
typename Scalar>
349 for_all_base_elements(e, this->mesh)
350 update_constrained_nodes(e,
nullptr,
nullptr,
nullptr,
nullptr);
Stores one element of a mesh.
Used to pass the instances of Space around.
Common definitions for Hermes2D.
Stores one node of a mesh.
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.
unsigned bnd
1 = boundary node; 0 = inner node
unsigned char get_num_components() const
Returns 2 if this is a vector shapeset, 1 otherwise.
Represents a finite element space over a domain.
unsigned ref
the number of elements using the node