11 #ifndef BEMBEL_SRC_ANSATZSPACE_GLUE_HPP_
12 #define BEMBEL_SRC_ANSATZSPACE_GLUE_HPP_
38 namespace GlueRoutines {
43 std::vector<int> dofs;
49 template <
typename Derived,
unsigned int DF>
51 static std::vector<dofIdentification> makeIdentification(
52 const std::vector<std::array<int, 4>> &edges_,
54 assert(
false &&
"Needs to be specialized");
75 template <
typename Derived>
77 std::vector<std::array<int, 4>> edges_;
79 Eigen::SparseMatrix<double> glue_mat_;
158 std::vector<Eigen::Triplet<double>> trips;
159 std::vector<GlueRoutines::dofIdentification> dof_id =
163 std::vector<bool> dof_is_slave(pre_dofs,
false);
164 std::vector<bool> dof_is_master(pre_dofs,
false);
167 int number_of_slaves = 0;
168 for (
auto dofset : dof_id) {
173 std::sort(dofset.dofs.begin(), dofset.dofs.end());
174 dof_is_master[dofset.dofs[0]] =
true;
175 for (
int i = 1; i < dofset.dofs.size(); ++i) {
176 dof_is_slave[dofset.dofs[i]] =
true;
182 std::sort(dof_id.begin(), dof_id.end(),
185 return a.dofs[0] < b.dofs[0];
188 const int post_dofs = pre_dofs - number_of_slaves;
194 int master_index = 0;
195 for (
int i = 0; i < post_dofs; ++i) {
196 const int post_index = i;
197 while (dof_is_slave[i + skip] && ((i + skip) < pre_dofs)) ++skip;
198 const int pre_index = i + skip;
199 trips.push_back(Eigen::Triplet<double>(pre_index, post_index, 1));
201 assert(!(dof_is_master[pre_index] && dof_is_slave[pre_index]));
202 if (dof_is_master[pre_index]) {
206 assert(pre_index == dof_id[master_index].dofs[0]);
207 const int number_of_partners = dof_id[master_index].dofs.size();
208 for (
int j = 1; j < number_of_partners; ++j) {
209 trips.push_back(Eigen::Triplet<double>(dof_id[master_index].dofs[j],
211 dof_id[master_index].coef));
217 Eigen::SparseMatrix<double> glue_matrix(pre_dofs, post_dofs);
218 glue_matrix.setFromTriplets(trips.begin(), trips.end());
220 assert(trips.size() == pre_dofs);
225 namespace GlueRoutines {
249 std::vector<int> getEdgeDofIndices(
int edgeCase,
int dimXdir,
int dimYdir,
251 std::vector<int> out;
254 out.reserve(dimXdir);
256 "This should be a given. If not something went wrong in Glue.hpp" &&
258 for (
int i = 0; i < dimXdir; ++i) {
259 out.push_back(i + shift);
264 out.reserve(dimYdir);
266 "This should be a given. If not something went wrong in Glue.hpp" &&
268 for (
int i = 0; i < dimYdir; ++i) {
269 out.push_back(dimXdir * (i + 1) - 1 + shift);
274 out.reserve(dimXdir);
276 "This should be a given. If not something went wrong in Glue.hpp" &&
278 for (
int i = 0; i < dimXdir; ++i) {
279 out.push_back(dimXdir * (dimYdir - 1) + i + shift);
284 out.reserve(dimYdir);
286 "This should be a given. If not something went wrong in Glue.hpp" &&
288 for (
int i = 0; i < dimYdir; ++i) {
289 out.push_back(dimXdir * i + shift);
303 inline bool edgeIsForwardParametrized(
int edgeCase) {
304 return edgeCase == 0 || edgeCase == 1;
306 inline bool edgeIsBackwardsParametrized(
int edgeCase) {
307 return !(edgeIsForwardParametrized(edgeCase));
309 inline bool normalComponentIsInwardDirected(
int edgeCase) {
310 return edgeCase == 0 || edgeCase == 3;
312 inline bool normalComponentIsOutwardDirected(
int edgeCase) {
313 return !(normalComponentIsInwardDirected(edgeCase));
316 inline bool reverseParametrized(
const std::array<int, 4> &edge) {
317 return ((edgeIsForwardParametrized(edge[2]) &&
318 edgeIsForwardParametrized(edge[3])) ||
319 (edgeIsBackwardsParametrized(edge[2]) &&
320 edgeIsBackwardsParametrized(edge[3])))
325 inline int glueCoefficientDivergenceConforming(
const std::array<int, 4> &edge) {
326 return ((normalComponentIsInwardDirected(edge[2]) &&
327 normalComponentIsInwardDirected(edge[3])) ||
328 (normalComponentIsOutwardDirected(edge[2]) &&
329 normalComponentIsOutwardDirected(edge[3])))
333 inline bool edgeToBeGluedInFirstComp(
int edgeCase) {
334 return (edgeCase == 0 || edgeCase == 2) ? false :
true;
342 template <
typename Derived>
344 static std::vector<dofIdentification> makeIdentification(
345 const std::vector<std::array<int, 4>> &edges_,
360 template <
typename Derived>
362 static std::vector<dofIdentification> makeIdentification(
363 const std::vector<std::array<int, 4>> &edges_,
366 assert(superspace.get_polynomial_degree() >= 1);
367 const int one_d_dim = superspace.get_polynomial_degree() + 1 +
369 ((1 << superspace.get_refinement_level()) - 1);
370 const int dofs_per_patch = one_d_dim * one_d_dim;
372 std::vector<GlueRoutines::dofIdentification> out;
373 out.reserve(edges_.size() * one_d_dim);
377 for (
auto edge : edges_) {
380 const int shift_e1 = edge[0] * dofs_per_patch;
381 const int shift_e2 = edge[1] * dofs_per_patch;
382 const bool needReversion = reverseParametrized(edge);
384 std::vector<int> dofs_e1 =
385 getEdgeDofIndices(edge[2], one_d_dim, one_d_dim, shift_e1);
386 std::vector<int> dofs_e2 =
387 getEdgeDofIndices(edge[3], one_d_dim, one_d_dim, shift_e2);
390 if (edge[1] > -1 && edge[0] > -1) {
391 for (
int i = 0; i < one_d_dim; ++i) {
393 const int j = needReversion ? one_d_dim - 1 - i : i;
395 assert(dofs_e1[i] != dofs_e2[j] &&
396 "If this happens something went horribly wrong.");
398 const int small_dof = std::min(dofs_e1[i], dofs_e2[j]);
399 const int large_dof = std::max(dofs_e1[i], dofs_e2[j]);
405 if (already_stored_in[small_dof] > -1) {
409 if (already_stored_in[large_dof] == -1) {
410 out[already_stored_in[small_dof]].dofs.push_back(large_dof);
411 already_stored_in[large_dof] = already_stored_in[small_dof];
413 if (!(already_stored_in[small_dof] ==
414 already_stored_in[large_dof])) {
419 const int small_store = std::min(already_stored_in[large_dof],
420 already_stored_in[small_dof]);
421 const int large_store = std::max(already_stored_in[large_dof],
422 already_stored_in[small_dof]);
423 for (
auto dfs : out[large_store].dofs) {
426 already_stored_in[dfs] = small_store;
427 for (
auto other_dfs : out[small_store].dofs) {
428 assert(dfs != other_dfs);
430 out[small_store].dofs.push_back(dfs);
433 out[large_store].dofs = {};
434 out[large_store].coef = 0;
437 }
else if (already_stored_in[large_dof] > -1) {
441 if (already_stored_in[small_dof] == -1) {
442 out[already_stored_in[large_dof]].dofs.push_back(small_dof);
443 already_stored_in[small_dof] = already_stored_in[large_dof];
445 if (!(already_stored_in[small_dof] ==
446 already_stored_in[large_dof])) {
451 const int small_store = std::min(already_stored_in[large_dof],
452 already_stored_in[small_dof]);
453 const int large_store = std::max(already_stored_in[large_dof],
454 already_stored_in[small_dof]);
455 for (
auto dfs : out[large_store].dofs) {
458 already_stored_in[dfs] = small_store;
459 for (
auto other_dfs : out[small_store].dofs) {
460 assert(dfs != other_dfs);
462 out[small_store].dofs.push_back(dfs);
465 out[large_store].dofs = {};
466 out[large_store].coef = 0;
472 d.dofs.push_back(small_dof);
473 d.dofs.push_back(large_dof);
474 already_stored_in[small_dof] = d_count;
475 already_stored_in[large_dof] = d_count;
485 for (
auto x = out.begin(); x != out.end(); ++x) {
486 if ((*x).dofs.size() == 0) {
501 template <
typename Derived>
503 static std::vector<dofIdentification> makeIdentification(
504 const std::vector<std::array<int, 4>> &edges_,
508 const int small_dim = superspace.get_polynomial_degree() +
510 ((1 << superspace.get_refinement_level()) - 1);
511 const int large_dim = superspace.get_polynomial_degree() + 1 +
513 ((1 << superspace.get_refinement_level()) - 1);
514 const int dofs_per_patch_per_component = small_dim * large_dim;
516 assert(dofs_per_patch_per_component ==
518 (superspace.get_number_of_patches() * 2)) &&
519 "The assembly of the glue matrix is highly specific to the space. "
520 "Something went wrong; is the discrete space correct?");
522 std::vector<GlueRoutines::dofIdentification> out;
523 out.reserve(edges_.size() * small_dim);
525 for (
auto edge : edges_) {
526 assert(edge[0] <= edge[1] || edge[1] == -1);
530 const int shift_e1 = edge[0] * dofs_per_patch_per_component +
531 (edgeToBeGluedInFirstComp(edge[2])
533 : (dofs_per_patch_per_component *
534 superspace.get_number_of_patches()));
535 const int shift_e2 = edge[1] * dofs_per_patch_per_component +
536 (edgeToBeGluedInFirstComp(edge[3])
538 : (dofs_per_patch_per_component *
539 superspace.get_number_of_patches()));
541 const int xdir_dim_1 =
542 edgeToBeGluedInFirstComp(edge[2]) ? large_dim : small_dim;
543 const int ydir_dim_1 =
544 edgeToBeGluedInFirstComp(edge[2]) ? small_dim : large_dim;
545 const int xdir_dim_2 =
546 edgeToBeGluedInFirstComp(edge[3]) ? large_dim : small_dim;
547 const int ydir_dim_2 =
548 edgeToBeGluedInFirstComp(edge[3]) ? small_dim : large_dim;
549 std::vector<int> dofs_e1 =
550 getEdgeDofIndices(edge[2], xdir_dim_1, ydir_dim_1, shift_e1);
551 std::vector<int> dofs_e2 =
552 getEdgeDofIndices(edge[3], xdir_dim_2, ydir_dim_2, shift_e2);
554 const int size_of_edge_dofs = dofs_e1.size();
555 assert(size_of_edge_dofs == dofs_e2.size() || edge[1] == -1);
557 const bool needReversion = reverseParametrized(edge);
558 const int coef = glueCoefficientDivergenceConforming(edge);
561 if (edge[1] > -1 && edge[0] > -1) {
562 for (
int i = 0; i < size_of_edge_dofs; ++i) {
564 const int j = needReversion ? size_of_edge_dofs - 1 - i : i;
565 assert(dofs_e1[i] != dofs_e2[j] &&
566 "If this happens something went horribly wrong.");
567 d.dofs.push_back(std::min(dofs_e1[i], dofs_e2[j]));
568 d.dofs.push_back(std::max(dofs_e1[i], dofs_e2[j]));
ElementTree & get_element_tree()
getter
std::vector< std::array< int, 4 > > patchTopologyInfo() const
Resolves neighborhood relations of the patches.
This class takes care of identifying DOFs on different edges, which must be identified with one anoth...
void init_Glue(const SuperSpace< Derived > &superspace, const Projector< Derived > &proj)
Initializes the AnsatzSpace object.
Eigen::SparseMatrix< double > get_glue_matrix() const
Returns Glue matrix to assemble continuous B-splines.
std::vector< GlueRoutines::dofIdentification > makeDofIdentificationList(const SuperSpace< Derived > &superspace, const Projector< Derived > &proj)
Generates a list of degrees of freedom (DOFs) identifications.
Glue(const SuperSpace< Derived > &superspace, const Projector< Derived > &proj)
Constructor for the Glue class.
Eigen::SparseMatrix< double > assembleGlueMatrix(const SuperSpace< Derived > &superspace, const Projector< Derived > &proj)
Assembles the Glue matrix according to the DOF identification list.
The Projector provides routines to assemble smooth B-Splines on each patch.
int get_knot_repetition() const
getter
int get_dofs_after_projector() const
Return number of smooth B-splines.
Routines for the evalutation of pointwise errors.
Helper struct for assembling the Glue for the different cases.
struct containing specifications on the linear operator has to be specialized or derived for any part...
The superspace manages local polynomial bases on each element of the mesh and provides an interface t...