23 #include "getfem/getfem_generic_assembly_compile_and_exec.h"
24 #include "getfem/getfem_generic_assembly_functions_and_operators.h"
33 void ga_interpolation(ga_workspace &workspace,
34 ga_interpolation_context &gic) {
35 ga_instruction_set gis;
36 ga_compile_interpolation(workspace, gis);
37 ga_interpolation_exec(gis, workspace, gic);
41 struct ga_interpolation_context_fem_same_mesh
42 :
public ga_interpolation_context {
44 std::vector<int> dof_count;
50 virtual bgeot::pstored_point_tab
52 std::vector<size_type> &ind)
const {
53 pfem pf = mf.fem_of_element(cv);
54 GMM_ASSERT1(pf->is_lagrange(),
55 "Only Lagrange fems can be used in interpolation");
58 for (
size_type i = 0; i < pf->node_convex(cv).structure()
59 ->nb_points_of_face(f); ++i)
60 ind.push_back(pf->node_convex(cv).structure()
61 ->ind_points_of_face(f)[i]);
63 for (
size_type i = 0; i < pf->node_convex(cv).nb_points(); ++i)
66 return pf->node_tab(cv);
69 virtual bool use_pgp(
size_type)
const {
return true; }
70 virtual bool use_mim()
const {
return false; }
83 size_type target_dim = mf.fem_of_element(cv)->target_dim();
84 GMM_ASSERT2(target_dim == 3,
"Invalid torus fem.");
87 if (!initialized) {init_(qdim, qdim, qdim);}
88 size_type idof = mf.ind_basic_dof_of_element(cv)[i];
89 result[idof] = t[idof%result_dim];
95 GMM_ASSERT1( (si % q) == 0,
"Incompatibility between the mesh_fem and "
96 "the size of the expression to be interpolated");
97 if (!initialized) { init_(si, q, qmult); }
98 GMM_ASSERT1(s == si,
"Internal error");
99 size_type idof = mf.ind_basic_dof_of_element(cv)[i*q];
101 gmm::sub_vector(result, gmm::sub_interval(qmult*idof, s)));
102 (dof_count[idof/q])++;
106 virtual void finalize() {
107 std::vector<size_type> data(3);
108 data[0] = initialized ? result.size() : 0;
109 data[1] = initialized ? dof_count.size() : 0;
110 data[2] = initialized ? s : 0;
111 MPI_MAX_VECTOR(data);
122 GMM_ASSERT1(gmm::vect_size(result) == data[0] && s == data[2] &&
123 gmm::vect_size(dof_count) == data[1],
"Incompatible sizes");
124 MPI_SUM_VECTOR(result);
125 MPI_SUM_VECTOR(dof_count);
126 for (
size_type i = 0; i < dof_count.size(); ++i)
128 gmm::scale(gmm::sub_vector(result, gmm::sub_interval(s*i, s)),
129 scalar_type(1) / scalar_type(dof_count[i]));
132 virtual const mesh &linked_mesh() {
return mf.linked_mesh(); }
134 ga_interpolation_context_fem_same_mesh(
const mesh_fem &mf_, base_vector &r)
135 : result(r), mf(mf_), initialized(false), is_torus(false) {
136 GMM_ASSERT1(!(mf.is_reduced()),
137 "Interpolation on reduced fem is not allowed");
138 if (
dynamic_cast<const torus_mesh_fem*
>(&mf)){
139 auto first_cv = mf.first_convex_of_basic_dof(0);
140 auto target_dim = mf.fem_of_element(first_cv)->target_dim();
141 if (target_dim == 3) is_torus =
true;
146 void ga_interpolation_Lagrange_fem
147 (ga_workspace &workspace,
const mesh_fem &mf, base_vector &result) {
148 ga_interpolation_context_fem_same_mesh gic(mf, result);
149 ga_interpolation(workspace, gic);
152 void ga_interpolation_Lagrange_fem
153 (
const getfem::model &md,
const std::string &expr,
const mesh_fem &mf,
154 base_vector &result,
const mesh_region &rg) {
155 ga_workspace workspace(md);
156 workspace.add_interpolation_expression(expr, mf.linked_mesh(), rg);
157 ga_interpolation_Lagrange_fem(workspace, mf, result);
161 struct ga_interpolation_context_mti
162 :
public ga_interpolation_context {
164 const mesh_trans_inv &mti;
168 virtual bgeot::pstored_point_tab
170 std::vector<size_type> &ind)
const {
171 std::vector<size_type> itab;
172 mti.points_on_convex(cv, itab);
173 std::vector<base_node> pt_tab(itab.size());
174 for (
size_type i = 0; i < itab.size(); ++i) {
175 pt_tab[i] = mti.reference_coords()[itab[i]];
178 return store_point_tab(pt_tab);
181 virtual bool use_pgp(
size_type)
const {
return false; }
182 virtual bool use_mim()
const {
return false; }
192 GMM_ASSERT1(s == si,
"Internal error");
193 size_type ipt = mti.point_on_convex(cv, i);
196 gmm::sub_vector(result, gmm::sub_interval(s*dof_t, s)));
199 virtual void finalize() {
200 std::vector<size_type> data(2);
201 data[0] = initialized ? result.size() : 0;
202 data[1] = initialized ? s : 0;
203 MPI_MAX_VECTOR(data);
212 GMM_ASSERT1(gmm::vect_size(result) == data[0] && s == data[1],
213 "Incompatible sizes");
214 MPI_SUM_VECTOR(result);
217 virtual const mesh &linked_mesh() {
return mti.linked_mesh(); }
219 ga_interpolation_context_mti(
const mesh_trans_inv &mti_, base_vector &r,
221 : result(r), mti(mti_), initialized(false), nbdof(nbdof_) {
222 if (nbdof ==
size_type(-1)) nbdof = mti.nb_points();
227 void ga_interpolation_mti
228 (
const getfem::model &md,
const std::string &expr, mesh_trans_inv &mti,
229 base_vector &result,
const mesh_region &rg,
int extrapolation,
230 const mesh_region &rg_source,
size_type nbdof) {
232 ga_workspace workspace(md);
233 workspace.add_interpolation_expression(expr, mti.linked_mesh(), rg);
235 mti.distribute(extrapolation, rg_source);
236 ga_interpolation_context_mti gic(mti, result, nbdof);
237 ga_interpolation(workspace, gic);
242 struct ga_interpolation_context_im_data
243 :
public ga_interpolation_context {
249 virtual bgeot::pstored_point_tab
251 std::vector<size_type> &ind)
const {
252 pintegration_method pim =imd.linked_mesh_im().int_method_of_element(cv);
253 if (pim->type() == IM_NONE)
return bgeot::pstored_point_tab();
254 GMM_ASSERT1(pim->type() == IM_APPROX,
"Sorry, exact methods cannot "
255 "be used in high level generic assembly");
258 i_end = pim->approx_method()->nb_points_on_convex();
260 i_start = pim->approx_method()->ind_first_point_on_face(f);
261 i_end = i_start + pim->approx_method()->nb_points_on_face(f);
263 for (
size_type i = i_start; i < i_end; ++i) ind.push_back(i);
264 return pim->approx_method()->pintegration_points();
267 virtual bool use_pgp(
size_type cv)
const {
268 pintegration_method pim =imd.linked_mesh_im().int_method_of_element(cv);
269 if (pim->type() == IM_NONE)
return false;
270 GMM_ASSERT1(pim->type() == IM_APPROX,
"Sorry, exact methods cannot "
271 "be used in high level generic assembly");
272 return !(pim->approx_method()->is_built_on_the_fly());
274 virtual bool use_mim()
const {
return true; }
280 GMM_ASSERT1(imd.tensor_size() == t.sizes() ||
281 (imd.tensor_size().size() ==
size_type(1) &&
284 "Im_data tensor size " << imd.tensor_size() <<
285 " does not match the size of the interpolated "
286 "expression " << t.sizes() <<
".");
291 GMM_ASSERT1(s == si,
"Internal error");
292 size_type ipt = imd.filtered_index_of_point(cv, i);
294 "Im data with no data on the current integration point.");
296 gmm::sub_vector(result, gmm::sub_interval(s*ipt, s)));
299 virtual void finalize() {
300 std::vector<size_type> data(2);
301 data[0] = initialized ? result.size() : 0;
302 data[1] = initialized ? s : 0;
303 MPI_MAX_VECTOR(data);
305 GMM_ASSERT1(gmm::vect_size(result) == data[0] && s == data[1],
306 "Incompatible sizes");
314 MPI_SUM_VECTOR(result);
317 virtual const mesh &linked_mesh() {
return imd.linked_mesh(); }
319 ga_interpolation_context_im_data(
const im_data &imd_, base_vector &r)
320 : result(r), imd(imd_), initialized(false) { }
323 void ga_interpolation_im_data
324 (ga_workspace &workspace,
const im_data &imd, base_vector &result) {
325 ga_interpolation_context_im_data gic(imd, result);
326 ga_interpolation(workspace, gic);
329 void ga_interpolation_im_data
330 (
const getfem::model &md,
const std::string &expr,
const im_data &imd,
331 base_vector &result,
const mesh_region &rg) {
332 ga_workspace workspace(md);
333 workspace.add_interpolation_expression
334 (expr, imd.linked_mesh_im(), rg);
336 ga_interpolation_im_data(workspace, imd, result);
341 struct ga_interpolation_context_mesh_slice
342 :
public ga_interpolation_context {
344 const stored_mesh_slice &sl;
347 std::vector<size_type> first_node;
349 virtual bgeot::pstored_point_tab
351 std::vector<size_type> &ind)
const {
352 GMM_ASSERT1(f ==
short_type(-1),
"No support for interpolation on faces"
353 " for a stored_mesh_slice yet.");
355 const mesh_slicer::cs_nodes_ct &nodes = sl.nodes(ic);
356 std::vector<base_node> pt_tab(nodes.size());
357 for (
size_type i=0; i < nodes.size(); ++i) {
358 pt_tab[i] = nodes[i].pt_ref;
361 return store_point_tab(pt_tab);
364 virtual bool use_pgp(
size_type )
const {
return false; }
365 virtual bool use_mim()
const {
return false; }
374 first_node.resize(sl.nb_convex());
375 for (
size_type ic=0; ic < sl.nb_convex()-1; ++ic)
376 first_node[ic+1] = first_node[ic] + sl.nodes(ic).size();
378 GMM_ASSERT1(s == si && result.size() == s * sl.nb_points(),
"Internal error");
382 gmm::sub_vector(result, gmm::sub_interval(s*ipt, s)));
385 virtual void finalize() {
386 std::vector<size_type> data(2);
387 data[0] = initialized ? result.size() : 0;
388 data[1] = initialized ? s : 0;
389 MPI_MAX_VECTOR(data);
391 GMM_ASSERT1(gmm::vect_size(result) == data[0] && s == data[1],
392 "Incompatible sizes");
400 MPI_SUM_VECTOR(result);
403 virtual const mesh &linked_mesh() {
return sl.linked_mesh(); }
405 ga_interpolation_context_mesh_slice(
const stored_mesh_slice &sl_,
407 : result(r), sl(sl_), initialized(false) { }
410 void ga_interpolation_mesh_slice
411 (ga_workspace &workspace,
const stored_mesh_slice &sl, base_vector &result) {
412 ga_interpolation_context_mesh_slice gic(sl, result);
413 ga_interpolation(workspace, gic);
416 void ga_interpolation_mesh_slice
417 (
const getfem::model &md,
const std::string &expr,
const stored_mesh_slice &sl,
418 base_vector &result,
const mesh_region &rg) {
419 ga_workspace workspace(md);
420 workspace.add_interpolation_expression(expr, sl.linked_mesh(), rg);
421 ga_interpolation_mesh_slice(workspace, sl, result);
430 const std::string &expr,
const mesh_fem &mf,
439 model_real_sparse_matrix M(nbdof, nbdof);
443 ga_workspace workspace(md, ga_workspace::inherit::NONE);
444 gmm::sub_interval I(0,nbdof);
445 workspace.add_fem_variable(
"c__dummy_var_95_", mf, I, base_vector(nbdof));
446 if (mf.get_qdims().size() > 1)
447 workspace.add_expression(
"("+expr+
"):Test_c__dummy_var_95_",mim,region,2);
449 workspace.add_expression(
"("+expr+
").Test_c__dummy_var_95_",mim,region,2);
450 getfem::base_vector F(nbdof);
451 workspace.set_assembled_vector(F);
452 workspace.assembly(1);
453 gmm::resize(result, nbdof);
455 getfem::base_matrix loc_M;
456 getfem::base_vector loc_U;
460 loc_M.base_resize(nd, nd); gmm::resize(loc_U, nd);
462 gmm::copy(gmm::sub_matrix(M, J, J), loc_M);
463 gmm::lu_solve(loc_M, loc_U, gmm::sub_vector(F, J));
464 gmm::copy(loc_U, gmm::sub_vector(result, J));
467 MPI_SUM_VECTOR(result);
474 struct rated_box_index_compare {
476 const std::pair<scalar_type, const bgeot::box_index*> &x,
477 const std::pair<scalar_type, const bgeot::box_index*> &y)
const {
478 GMM_ASSERT2(x.second !=
nullptr,
"x contains nullptr");
479 GMM_ASSERT2(y.second !=
nullptr,
"y contains nullptr");
480 return (x.first < y.first) || (!(y.first < x.first) && (x.second->id < y.second->id));
484 class interpolate_transformation_expression
485 :
public virtual_interpolate_transformation,
public context_dependencies {
487 struct workspace_gis_pair :
public std::pair<ga_workspace, ga_instruction_set> {
488 inline ga_workspace &workspace() {
return this->first; }
489 inline ga_instruction_set &gis() {
return this->second; }
492 const mesh &source_mesh;
493 const mesh &target_mesh;
497 mutable bool recompute_elt_boxes;
498 mutable ga_workspace local_workspace;
499 mutable ga_instruction_set local_gis;
502 mutable std::set<var_trans_pair> used_vars;
503 mutable std::set<var_trans_pair> used_data;
504 mutable std::map<var_trans_pair,
505 workspace_gis_pair> compiled_derivatives;
506 mutable bool extract_variable_done;
507 mutable bool extract_data_done;
510 mutable std::map<size_type, std::vector<size_type>> box_to_convexes;
513 void update_from_context()
const {
514 recompute_elt_boxes =
true;
517 void extract_variables(
const ga_workspace &workspace,
518 std::set<var_trans_pair> &vars,
519 bool ignore_data,
const mesh &,
520 const std::string &)
const {
521 if ((ignore_data && !extract_variable_done) ||
522 (!ignore_data && !extract_data_done)) {
527 ga_workspace aux_workspace(workspace, ga_workspace::inherit::ALL);
528 aux_workspace.clear_expressions();
529 aux_workspace.add_interpolation_expression(expr, source_mesh);
530 for (
size_type i = 0; i < aux_workspace.nb_trees(); ++i)
531 ga_extract_variables(aux_workspace.tree_info(i).ptree->root,
532 aux_workspace, source_mesh,
533 ignore_data ? used_vars : used_data,
536 extract_variable_done =
true;
538 extract_data_done =
true;
541 vars.insert(used_vars.begin(), used_vars.end());
543 vars.insert(used_data.begin(), used_data.end());
546 void init(
const ga_workspace &workspace)
const {
550 local_workspace = ga_workspace(workspace, ga_workspace::inherit::ALL);
551 local_workspace.clear_expressions();
553 local_workspace.add_interpolation_expression(expr, source_mesh);
554 local_gis = ga_instruction_set();
555 ga_compile_interpolation(local_workspace, local_gis);
558 for (
const std::string &transname : local_gis.transformations)
559 local_workspace.interpolate_transformation(transname)
560 ->init(local_workspace);
562 if (!extract_variable_done) {
563 std::set<var_trans_pair> vars;
564 extract_variables(workspace, vars,
true, source_mesh,
"");
567 for (
const var_trans_pair &var : used_vars) {
568 workspace_gis_pair &pwi = compiled_derivatives[var];
569 pwi.workspace() = local_workspace;
570 pwi.gis() = ga_instruction_set();
571 if (pwi.workspace().nb_trees()) {
572 ga_tree der_tree = *(pwi.workspace().tree_info(0).ptree);
573 ga_derivative(der_tree, pwi.workspace(), source_mesh,
574 var.varname, var.transname, 1);
575 pwi.workspace().clear_expressions();
576 const std::string der_expr=ga_tree_to_string(der_tree);
577 pwi.workspace().add_interpolation_expression(der_expr, target_mesh);
578 ga_compile_interpolation(pwi.workspace(), pwi.gis());
583 if (recompute_elt_boxes) {
585 box_to_convexes.clear();
586 element_boxes.clear();
587 base_node bmin(N), bmax(N);
588 const mesh_region &mr = target_mesh.region(target_region);
589 const dal::bit_vector&
591 ? target_mesh.convex_index() : mr.index();
592 for (dal::bv_visitor cv(convex_index); !cv.finished(); ++cv) {
598 gmm::copy(target_mesh.points_of_convex(cv)[0], bmin);
604 mesh::ref_mesh_pt_ct cvpts = target_mesh.points_of_convex(cv);
607 const base_node &pt = cvpts[ip];
609 bmin[k] = std::min(bmin[k], pt[k]);
610 bmax[k] = std::max(bmax[k], pt[k]);
614 scalar_type h = bmax[0] - bmin[0];
615 for (
size_type k = 1; k < N; ++k) h = std::max(h, bmax[k]-bmin[k]);
616 if (pgt->is_linear()) h *= 1E-4;
617 for (
auto &&val : bmin) val -= h*0.2;
618 for (
auto &&val : bmax) val += h*0.2;
620 box_to_convexes[element_boxes.add_box(bmin, bmax)].push_back(cv);
622 element_boxes.build_tree();
623 recompute_elt_boxes =
false;
627 void finalize()
const {
628 for (
const std::string &transname : local_gis.transformations)
629 local_workspace.interpolate_transformation(transname)->finalize();
630 local_gis = ga_instruction_set();
633 std::string expression()
const {
return expr; }
635 int transform(
const ga_workspace &,
const mesh &m,
636 fem_interpolation_context &ctx_x,
637 const base_small_vector &Normal,
642 std::map<var_trans_pair, base_tensor> &derivatives,
643 bool compute_derivatives)
const {
646 local_gis.ctx = ctx_x;
647 local_gis.Normal = Normal;
651 gmm::clear(local_workspace.assembled_tensor().as_vector());
653 for (
auto &&instr : local_gis.all_instructions) {
654 GMM_ASSERT1(instr.second.m == &m,
655 "Incompatibility of meshes in interpolation");
656 auto &gilb = instr.second.begin_instructions;
657 for (
size_type j = 0; j < gilb.size(); ++j) j += gilb[j]->exec();
658 auto &gile = instr.second.elt_instructions;
659 for (
size_type j = 0; j < gile.size(); ++j) j+=gile[j]->exec();
660 auto &gil = instr.second.instructions;
661 for (
size_type j = 0; j < gil.size(); ++j) j += gil[j]->exec();
664 GMM_ASSERT1(local_workspace.assembled_tensor().size()==target_mesh.dim(),
665 "Wrong dimension of the transformation expression");
666 P.resize(target_mesh.dim());
667 gmm::copy(local_workspace.assembled_tensor().as_vector(), P);
671 bgeot::rtree::pbox_cont boxes;
673 bgeot::rtree::pbox_set bset;
674 element_boxes.find_boxes_at_point(P, bset);
677 std::set<std::pair<scalar_type, const bgeot::box_index*>,
678 rated_box_index_compare> rated_boxes;
679 for (
const auto &box : bset) {
680 scalar_type rating = scalar_type(1);
681 for (
size_type i = 0; i < target_mesh.dim(); ++i) {
682 scalar_type h = box->max->at(i) - box->min->at(i);
683 if (h > scalar_type(0)) {
684 scalar_type r = std::min(box->max->at(i) - P[i],
685 P[i] - box->min->at(i)) / h;
686 rating = std::min(r, rating);
689 rated_boxes.insert(std::make_pair(rating, box));
693 for (
const auto &p : rated_boxes)
694 boxes.push_back(p.second);
697 scalar_type best_dist(1e10);
699 base_node best_P_ref;
700 for (
size_type i = boxes.size(); i > 0 && !ret_type; --i) {
701 for (
auto convex : box_to_convexes.at(boxes[i-1]->id)) {
702 gic.init(target_mesh.points_of_convex(convex),
703 target_mesh.trans_of_convex(convex));
706 bool is_in = gic.
invert(P, P_ref, converged, 1E-4);
721 = target_mesh.trans_of_convex(cv)->convex_ref()->is_in(P_ref);
722 if (dist < best_dist) {
732 if (ret_type == 0 && best_dist < 5e-3) {
744 if (compute_derivatives) {
745 for (
auto &&d : derivatives) {
746 workspace_gis_pair &pwi = compiled_derivatives[d.first];
747 pwi.gis().ctx = ctx_x;
748 pwi.gis().Normal = Normal;
752 gmm::clear(pwi.workspace().assembled_tensor().as_vector());
753 ga_function_exec(pwi.gis());
754 d.second = pwi.workspace().assembled_tensor();
761 interpolate_transformation_expression
762 (
const mesh &sm,
const mesh &tm,
size_type trg,
const std::string &expr_)
763 : source_mesh(sm), target_mesh(tm), target_region(trg), expr(expr_),
764 recompute_elt_boxes(true), extract_variable_done(false),
765 extract_data_done(false)
766 { this->add_dependency(tm); }
772 (ga_workspace &workspace,
const std::string &name,
const mesh &sm,
773 const mesh &tm,
const std::string &expr) {
775 (workspace, name, sm, tm,
size_type(-1), expr);
779 (ga_workspace &workspace,
const std::string &name,
const mesh &sm,
780 const mesh &tm,
size_type trg,
const std::string &expr) {
781 pinterpolate_transformation
782 p = std::make_shared<interpolate_transformation_expression>
784 workspace.add_interpolate_transformation(name, p);
788 (model &md,
const std::string &name,
const mesh &sm,
const mesh &tm,
789 const std::string &expr) {
795 (model &md,
const std::string &name,
const mesh &sm,
const mesh &tm,
796 size_type trg,
const std::string &expr) {
797 pinterpolate_transformation
798 p = std::make_shared<interpolate_transformation_expression>
800 md.add_interpolate_transformation(name, p);
807 class interpolate_transformation_neighbor
808 :
public virtual_interpolate_transformation,
public context_dependencies {
811 void update_from_context()
const {}
812 void extract_variables(
const ga_workspace &,
813 std::set<var_trans_pair> &,
815 const std::string &)
const {}
816 void init(
const ga_workspace &)
const {}
817 void finalize()
const {}
819 std::string expression()
const {
return "X"; }
821 int transform(
const ga_workspace &,
const mesh &m_x,
822 fem_interpolation_context &ctx_x,
823 const base_small_vector &,
const mesh **m_t,
827 std::map<var_trans_pair, base_tensor> &,
828 bool compute_derivatives)
const {
834 GMM_ASSERT1(face_x !=
short_type(-1),
"Neighbor transformation can "
835 "only be applied to internal faces");
837 auto adj_face = m_x.adjacent_face(cv_x, face_x);
841 gic.init(m_x.points_of_convex(adj_face.cv),
842 m_x.trans_of_convex(adj_face.cv));
843 bool converged =
true;
844 gic.
invert(ctx_x.xreal(), P_ref, converged);
845 bool is_in = (ctx_x.pgt()->convex_ref()->is_in(P_ref) < 1E-4);
846 GMM_ASSERT1(is_in && converged,
"Geometric transformation inversion "
847 "has failed in neighbor transformation");
848 face_num = adj_face.f;
852 GMM_ASSERT1(!compute_derivatives,
853 "No derivative for this transformation");
857 interpolate_transformation_neighbor() { }
863 return (std::make_shared<interpolate_transformation_neighbor>());
870 class interpolate_transformation_element_extrapolation
871 :
public virtual_interpolate_transformation,
public context_dependencies {
874 std::map<size_type, size_type> elt_corr;
877 void update_from_context()
const {}
878 void extract_variables(
const ga_workspace &,
879 std::set<var_trans_pair> &,
881 const std::string &)
const {}
882 void init(
const ga_workspace &)
const {}
883 void finalize()
const {}
884 std::string expression()
const {
return "X"; }
886 int transform(
const ga_workspace &,
const mesh &m_x,
887 fem_interpolation_context &ctx_x,
888 const base_small_vector &,
const mesh **m_t,
892 std::map<var_trans_pair, base_tensor> &,
893 bool compute_derivatives)
const {
896 GMM_ASSERT1(&sm == &m_x,
"Bad mesh");
897 size_type cv_x = ctx_x.convex_num(), cv_y = cv_x;
898 auto it = elt_corr.find(cv_x);
899 if (it != elt_corr.end()) cv_y = it->second;
903 gic.init(m_x.points_of_convex(cv_y),
904 m_x.trans_of_convex(cv_y));
905 bool converged =
true;
906 gic.
invert(ctx_x.xreal(), P_ref, converged, 1E-4);
907 GMM_ASSERT1(converged,
"Geometric transformation inversion "
908 "has failed in element extrapolation transformation");
915 P_ref = ctx_x.xref();
918 GMM_ASSERT1(!compute_derivatives,
919 "No derivative for this transformation");
923 void set_correspondence(
const std::map<size_type, size_type> &ec) {
927 interpolate_transformation_element_extrapolation
928 (
const mesh &sm_,
const std::map<size_type, size_type> &ec)
929 : sm(sm_), elt_corr(ec) { }
933 void add_element_extrapolation_transformation
934 (model &md,
const std::string &name,
const mesh &sm,
935 std::map<size_type, size_type> &elt_corr) {
936 pinterpolate_transformation
937 p = std::make_shared<interpolate_transformation_element_extrapolation>
939 md.add_interpolate_transformation(name, p);
942 void add_element_extrapolation_transformation
943 (ga_workspace &workspace,
const std::string &name,
const mesh &sm,
944 std::map<size_type, size_type> &elt_corr) {
945 pinterpolate_transformation
946 p = std::make_shared<interpolate_transformation_element_extrapolation>
948 workspace.add_interpolate_transformation(name, p);
951 void set_element_extrapolation_correspondence
952 (ga_workspace &workspace,
const std::string &name,
953 std::map<size_type, size_type> &elt_corr) {
954 GMM_ASSERT1(workspace.interpolate_transformation_exists(name),
955 "Unknown transformation");
956 auto pit=workspace.interpolate_transformation(name).get();
958 =
dynamic_cast<const interpolate_transformation_element_extrapolation *
>
961 "The transformation is not of element extrapolation type");
962 const_cast<interpolate_transformation_element_extrapolation *
>(cpext)
963 ->set_correspondence(elt_corr);
966 void set_element_extrapolation_correspondence
967 (model &md,
const std::string &name,
968 std::map<size_type, size_type> &elt_corr) {
969 GMM_ASSERT1(md.interpolate_transformation_exists(name),
970 "Unknown transformation");
971 auto pit=md.interpolate_transformation(name).get();
973 =
dynamic_cast<const interpolate_transformation_element_extrapolation *
>
976 "The transformation is not of element extrapolation type");
977 const_cast<interpolate_transformation_element_extrapolation *
>(cpext)
978 ->set_correspondence(elt_corr);
987 class standard_secondary_domain :
public virtual_secondary_domain {
991 virtual const mesh_region &give_region(
const mesh &,
997 standard_secondary_domain(
const mesh_im &mim__,
const mesh_region ®ion_)
998 : virtual_secondary_domain(mim__, region_) {}
1001 void add_standard_secondary_domain
1002 (model &md,
const std::string &name,
const mesh_im &mim,
1003 const mesh_region &rg) {
1004 psecondary_domain p = std::make_shared<standard_secondary_domain>(mim, rg);
1005 md.add_secondary_domain(name, p);
1008 void add_standard_secondary_domain
1009 (ga_workspace &workspace,
const std::string &name,
const mesh_im &mim,
1010 const mesh_region &rg) {
1011 psecondary_domain p = std::make_shared<standard_secondary_domain>(mim, rg);
1012 workspace.add_secondary_domain(name, p);
does the inversion of the geometric transformation for a given convex
bool invert(const base_node &n, base_node &n_ref, scalar_type IN_EPS=1e-12, bool project_into_element=false)
given the node on the real element, returns the node on the reference element (even if it is outside ...
Balanced tree of n-dimensional rectangles.
Describe a finite element method linked to a mesh.
virtual ind_dof_ct ind_basic_dof_of_element(size_type cv) const
Give an array of the dof numbers a of convex.
virtual size_type nb_dof() const
Return the total number of degrees of freedom.
const mesh & linked_mesh() const
Return a reference to the underlying mesh.
const dal::bit_vector & convex_index() const
Get the set of convexes where a finite element has been assigned.
virtual size_type nb_basic_dof_of_element(size_type cv) const
Return the number of degrees of freedom attached to a given convex.
Describe an integration method linked to a mesh.
"iterator" class for regions.
structure used to hold a set of convexes and/or convex faces.
static mesh_region all_convexes()
provide a default value for the mesh_region parameters of assembly procedures etc.
Describe a mesh (collection of convexes (elements) and points).
`‘Model’' variables store the variables, the data and the description of a model.
Semantic analysis of assembly trees and semantic manipulations.
Compilation and execution operations.
void copy(const L1 &l1, L2 &l2)
*/
void clear(L &l)
clear (fill with zeros) a vector or matrix.
void resize(V &v, size_type n)
*/
void add(const L1 &l1, L2 &l2)
*/
void asm_mass_matrix(const MAT &M, const mesh_im &mim, const mesh_fem &mf1, const mesh_region &rg=mesh_region::all_convexes())
generic mass matrix assembly (on the whole mesh or on the specified convex set or boundary)
std::shared_ptr< const getfem::virtual_fem > pfem
type of pointer on a fem description
gmm::uint16_type short_type
used as the common short type integer in the library
size_t size_type
used as the common size type in the library
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation
GEneric Tool for Finite Element Methods.
pinterpolate_transformation interpolate_transformation_neighbor_instance()
Create a new instance of a transformation corresponding to the interpolation on the neighbor element.
void add_interpolate_transformation_from_expression(ga_workspace &workspace, const std::string &transname, const mesh &source_mesh, const mesh &target_mesh, const std::string &expr)
Add a transformation to a workspace workspace or a model md mapping point in mesh source_mesh to mesh...
void ga_local_projection(const getfem::model &md, const mesh_im &mim, const std::string &expr, const mesh_fem &mf, base_vector &result, const mesh_region &rg=mesh_region::all_convexes())
Make an elementwise L2 projection of an expression with respect to the mesh_fem mf.