24 #include "getfem/getfem_generic_assembly_compile_and_exec.h" 25 #include "getfem/getfem_generic_assembly_functions_and_operators.h" 34 void ga_interpolation(ga_workspace &workspace,
35 ga_interpolation_context &gic) {
36 ga_instruction_set gis;
37 ga_compile_interpolation(workspace, gis);
38 ga_interpolation_exec(gis, workspace, gic);
42 struct ga_interpolation_context_fem_same_mesh
43 :
public ga_interpolation_context {
45 std::vector<int> dof_count;
51 virtual bgeot::pstored_point_tab
53 std::vector<size_type> &ind)
const {
54 pfem pf = mf.fem_of_element(cv);
55 GMM_ASSERT1(pf->is_lagrange(),
56 "Only Lagrange fems can be used in interpolation");
61 i < pf->node_convex(cv).structure()->nb_points_of_face(f); ++i)
63 (pf->node_convex(cv).structure()->ind_points_of_face(f)[i]);
65 for (
size_type i = 0; i < pf->node_convex(cv).nb_points(); ++i)
69 return pf->node_tab(cv);
72 virtual bool use_pgp(
size_type)
const {
return true; }
73 virtual bool use_mim()
const {
return false; }
85 size_type target_dim = mf.fem_of_element(cv)->target_dim();
86 GMM_ASSERT2(target_dim == 3,
"Invalid torus fem.");
89 if (!initialized) {init_(qdim, qdim, qdim);}
90 size_type idof = mf.ind_basic_dof_of_element(cv)[i];
91 result[idof] = t[idof%result_dim];
97 store_result_for_torus(cv, i, t);
103 GMM_ASSERT1( (si % q) == 0,
"Incompatibility between the mesh_fem and " 104 "the size of the expression to be interpolated");
105 if (!initialized) { init_(si, q, qmult); }
106 GMM_ASSERT1(s == si,
"Internal error");
107 size_type idof = mf.ind_basic_dof_of_element(cv)[i*q];
108 gmm::add(t.as_vector(),
109 gmm::sub_vector(result, gmm::sub_interval(qmult*idof, s)));
110 (dof_count[idof/q])++;
113 virtual void finalize() {
114 std::vector<size_type> data(3);
115 data[0] = initialized ? result.size() : 0;
116 data[1] = initialized ? dof_count.size() : 0;
117 data[2] = initialized ? s : 0;
118 MPI_MAX_VECTOR(data);
129 GMM_ASSERT1(gmm::vect_size(result) == data[0] && s == data[2] &&
130 gmm::vect_size(dof_count) == data[1],
"Incompatible sizes");
131 MPI_SUM_VECTOR(result);
132 MPI_SUM_VECTOR(dof_count);
133 for (
size_type i = 0; i < dof_count.size(); ++i)
135 gmm::scale(gmm::sub_vector(result, gmm::sub_interval(s*i, s)),
136 scalar_type(1) / scalar_type(dof_count[i]));
139 virtual const mesh &linked_mesh() {
return mf.linked_mesh(); }
141 ga_interpolation_context_fem_same_mesh(
const mesh_fem &mf_, base_vector &r)
142 : result(r), mf(mf_), initialized(false), is_torus(false) {
143 GMM_ASSERT1(!(mf.is_reduced()),
144 "Interpolation on reduced fem is not allowed");
145 if (dynamic_cast<const torus_mesh_fem*>(&mf)){
146 auto first_cv = mf.first_convex_of_basic_dof(0);
147 auto target_dim = mf.fem_of_element(first_cv)->target_dim();
148 if (target_dim == 3) is_torus =
true;
153 void ga_interpolation_Lagrange_fem
154 (ga_workspace &workspace,
const mesh_fem &mf, base_vector &result) {
155 ga_interpolation_context_fem_same_mesh gic(mf, result);
156 ga_interpolation(workspace, gic);
159 void ga_interpolation_Lagrange_fem
160 (
const getfem::model &md,
const std::string &expr,
const mesh_fem &mf,
161 base_vector &result,
const mesh_region &rg) {
162 ga_workspace workspace(md);
163 workspace.add_interpolation_expression(expr, mf.linked_mesh(), rg);
164 ga_interpolation_Lagrange_fem(workspace, mf, result);
168 struct ga_interpolation_context_mti
169 :
public ga_interpolation_context {
171 const mesh_trans_inv &mti;
176 virtual bgeot::pstored_point_tab
178 std::vector<size_type> &ind)
const {
179 std::vector<size_type> itab;
180 mti.points_on_convex(cv, itab);
181 std::vector<base_node> pt_tab(itab.size());
182 for (
size_type i = 0; i < itab.size(); ++i) {
183 pt_tab[i] = mti.reference_coords()[itab[i]];
186 return store_point_tab(pt_tab);
189 virtual bool use_pgp(
size_type)
const {
return false; }
190 virtual bool use_mim()
const {
return false; }
200 GMM_ASSERT1(s == si,
"Internal error");
201 size_type ipt = mti.point_on_convex(cv, i);
203 gmm::add(t.as_vector(),
204 gmm::sub_vector(result, gmm::sub_interval(s*dof_t, s)));
207 virtual void finalize() {
208 std::vector<size_type> data(2);
209 data[0] = initialized ? result.size() : 0;
210 data[1] = initialized ? s : 0;
211 MPI_MAX_VECTOR(data);
220 GMM_ASSERT1(gmm::vect_size(result) == data[0] && s == data[1],
221 "Incompatible sizes");
222 MPI_SUM_VECTOR(result);
225 virtual const mesh &linked_mesh() {
return mti.linked_mesh(); }
227 ga_interpolation_context_mti(
const mesh_trans_inv &mti_, base_vector &r,
229 : result(r), mti(mti_), initialized(false), nbdof(nbdof_) {
230 if (nbdof ==
size_type(-1)) nbdof = mti.nb_points();
235 void ga_interpolation_mti
236 (
const getfem::model &md,
const std::string &expr, mesh_trans_inv &mti,
237 base_vector &result,
const mesh_region &rg,
int extrapolation,
238 const mesh_region &rg_source,
size_type nbdof) {
240 ga_workspace workspace(md);
241 workspace.add_interpolation_expression(expr, mti.linked_mesh(), rg);
243 mti.distribute(extrapolation, rg_source);
244 ga_interpolation_context_mti gic(mti, result, nbdof);
245 ga_interpolation(workspace, gic);
250 struct ga_interpolation_context_im_data
251 :
public ga_interpolation_context {
257 virtual bgeot::pstored_point_tab
259 std::vector<size_type> &ind)
const {
260 pintegration_method pim =imd.linked_mesh_im().int_method_of_element(cv);
261 if (pim->type() == IM_NONE)
return bgeot::pstored_point_tab();
262 GMM_ASSERT1(pim->type() == IM_APPROX,
"Sorry, exact methods cannot " 263 "be used in high level generic assembly");
266 i_end = pim->approx_method()->nb_points_on_convex();
268 i_start = pim->approx_method()->ind_first_point_on_face(f);
269 i_end = i_start + pim->approx_method()->nb_points_on_face(f);
271 for (
size_type i = i_start; i < i_end; ++i) ind.push_back(i);
272 return pim->approx_method()->pintegration_points();
275 virtual bool use_pgp(
size_type cv)
const {
276 pintegration_method pim =imd.linked_mesh_im().int_method_of_element(cv);
277 if (pim->type() == IM_NONE)
return false;
278 GMM_ASSERT1(pim->type() == IM_APPROX,
"Sorry, exact methods cannot " 279 "be used in high level generic assembly");
280 return !(pim->approx_method()->is_built_on_the_fly());
282 virtual bool use_mim()
const {
return true; }
288 GMM_ASSERT1(imd.tensor_size() == t.sizes() ||
289 (imd.tensor_size().size() ==
size_type(1) &&
292 "Im_data tensor size " << imd.tensor_size() <<
293 " does not match the size of the interpolated " 294 "expression " << t.sizes() <<
".");
299 GMM_ASSERT1(s == si,
"Internal error");
300 size_type ipt = imd.filtered_index_of_point(cv, i);
302 "Im data with no data on the current integration point.");
303 gmm::add(t.as_vector(),
304 gmm::sub_vector(result, gmm::sub_interval(s*ipt, s)));
307 virtual void finalize() {
308 std::vector<size_type> data(2);
309 data[0] = initialized ? result.size() : 0;
310 data[1] = initialized ? s : 0;
311 MPI_MAX_VECTOR(data);
313 GMM_ASSERT1(gmm::vect_size(result) == data[0] && s == data[1],
314 "Incompatible sizes");
322 MPI_SUM_VECTOR(result);
325 virtual const mesh &linked_mesh() {
return imd.linked_mesh(); }
327 ga_interpolation_context_im_data(
const im_data &imd_, base_vector &r)
328 : result(r), imd(imd_), initialized(false) { }
331 void ga_interpolation_im_data
332 (ga_workspace &workspace,
const im_data &imd, base_vector &result) {
333 ga_interpolation_context_im_data gic(imd, result);
334 ga_interpolation(workspace, gic);
337 void ga_interpolation_im_data
338 (
const getfem::model &md,
const std::string &expr,
const im_data &imd,
339 base_vector &result,
const mesh_region &rg) {
340 ga_workspace workspace(md);
341 workspace.add_interpolation_expression
342 (expr, imd.linked_mesh_im(), rg);
344 ga_interpolation_im_data(workspace, imd, result);
349 struct ga_interpolation_context_mesh_slice
350 :
public ga_interpolation_context {
352 const stored_mesh_slice &sl;
355 std::vector<size_type> first_node;
357 virtual bgeot::pstored_point_tab
359 std::vector<size_type> &ind)
const {
360 GMM_ASSERT1(f ==
short_type(-1),
"No support for interpolation on faces" 361 " for a stored_mesh_slice yet.");
363 const mesh_slicer::cs_nodes_ct &nodes = sl.nodes(ic);
364 std::vector<base_node> pt_tab(nodes.size());
365 for (
size_type i=0; i < nodes.size(); ++i) {
366 pt_tab[i] = nodes[i].pt_ref;
369 return store_point_tab(pt_tab);
372 virtual bool use_pgp(
size_type )
const {
return false; }
373 virtual bool use_mim()
const {
return false; }
382 first_node.resize(sl.nb_convex());
383 for (
size_type ic=0; ic < sl.nb_convex()-1; ++ic)
384 first_node[ic+1] = first_node[ic] + sl.nodes(ic).size();
386 GMM_ASSERT1(s == si && result.size() == s * sl.nb_points(),
"Internal error");
389 gmm::add(t.as_vector(),
390 gmm::sub_vector(result, gmm::sub_interval(s*ipt, s)));
393 virtual void finalize() {
394 std::vector<size_type> data(2);
395 data[0] = initialized ? result.size() : 0;
396 data[1] = initialized ? s : 0;
397 MPI_MAX_VECTOR(data);
399 GMM_ASSERT1(gmm::vect_size(result) == data[0] && s == data[1],
400 "Incompatible sizes");
408 MPI_SUM_VECTOR(result);
411 virtual const mesh &linked_mesh() {
return sl.linked_mesh(); }
413 ga_interpolation_context_mesh_slice(
const stored_mesh_slice &sl_, base_vector &r)
414 : result(r), sl(sl_), initialized(false) { }
417 void ga_interpolation_mesh_slice
418 (ga_workspace &workspace,
const stored_mesh_slice &sl, base_vector &result) {
419 ga_interpolation_context_mesh_slice gic(sl, result);
420 ga_interpolation(workspace, gic);
423 void ga_interpolation_mesh_slice
424 (
const getfem::model &md,
const std::string &expr,
const stored_mesh_slice &sl,
425 base_vector &result,
const mesh_region &rg) {
426 ga_workspace workspace(md);
427 workspace.add_interpolation_expression(expr, sl.linked_mesh(), rg);
428 ga_interpolation_mesh_slice(workspace, sl, result);
437 const std::string &expr,
const mesh_fem &mf,
448 ga_workspace workspace(md);
450 gmm::sub_interval I(nbdof, mf.
nb_dof());
451 workspace.add_fem_variable(
"c__dummy_var_95_", mf, I, base_vector(nbdof));
452 if (mf.get_qdims().size() > 1)
453 workspace.add_expression(
"("+expr+
"):Test_c__dummy_var_95_",mim,region,2);
455 workspace.add_expression(
"("+expr+
").Test_c__dummy_var_95_",mim,region,2);
456 base_vector residual(nbdof+mf.
nb_dof());
457 workspace.set_assembled_vector(residual);
458 workspace.assembly(1);
459 getfem::base_vector F(mf.
nb_dof());
460 gmm::resize(result, mf.
nb_dof());
461 gmm::copy(gmm::sub_vector(residual, I), F);
463 getfem::base_matrix loc_M;
464 getfem::base_vector loc_U;
468 loc_M.base_resize(nd, nd); gmm::resize(loc_U, nd);
470 gmm::copy(gmm::sub_matrix(M, J, J), loc_M);
471 gmm::lu_solve(loc_M, loc_U, gmm::sub_vector(F, J));
472 gmm::copy(loc_U, gmm::sub_vector(result, J));
475 MPI_SUM_VECTOR(result);
482 class interpolate_transformation_expression
485 struct workspace_gis_pair :
public std::pair<ga_workspace, ga_instruction_set> {
486 inline ga_workspace &workspace() {
return this->first; }
487 inline ga_instruction_set &gis() {
return this->second; }
490 const mesh &source_mesh;
491 const mesh &target_mesh;
495 mutable bool recompute_elt_boxes;
496 mutable ga_workspace local_workspace;
497 mutable ga_instruction_set local_gis;
500 mutable std::set<var_trans_pair> used_vars;
501 mutable std::set<var_trans_pair> used_data;
502 mutable std::map<var_trans_pair,
503 workspace_gis_pair> compiled_derivatives;
504 mutable bool extract_variable_done;
505 mutable bool extract_data_done;
508 void update_from_context()
const {
509 recompute_elt_boxes =
true;
512 void extract_variables(
const ga_workspace &workspace,
513 std::set<var_trans_pair> &vars,
514 bool ignore_data,
const mesh &,
515 const std::string &)
const {
516 if ((ignore_data && !extract_variable_done) ||
517 (!ignore_data && !extract_data_done)) {
522 ga_workspace aux_workspace;
523 aux_workspace = ga_workspace(
true, workspace);
524 aux_workspace.clear_expressions();
525 aux_workspace.add_interpolation_expression(expr, source_mesh);
526 for (
size_type i = 0; i < aux_workspace.nb_trees(); ++i)
527 ga_extract_variables(aux_workspace.tree_info(i).ptree->root,
528 aux_workspace, source_mesh,
529 ignore_data ? used_vars : used_data,
532 extract_variable_done =
true;
534 extract_data_done =
true;
537 vars.insert(used_vars.begin(), used_vars.end());
539 vars.insert(used_data.begin(), used_data.end());
542 void init(
const ga_workspace &workspace)
const {
546 local_workspace = ga_workspace(
true, workspace);
547 local_workspace.clear_expressions();
549 local_workspace.add_interpolation_expression(expr, source_mesh);
550 local_gis = ga_instruction_set();
551 ga_compile_interpolation(local_workspace, local_gis);
554 for (
const std::string &transname : local_gis.transformations)
555 local_workspace.interpolate_transformation(transname)
556 ->init(local_workspace);
558 if (!extract_variable_done) {
559 std::set<var_trans_pair> vars;
560 extract_variables(workspace, vars,
true, source_mesh,
"");
563 for (
const var_trans_pair &var : used_vars) {
564 workspace_gis_pair &pwi = compiled_derivatives[var];
565 pwi.workspace() = local_workspace;
566 pwi.gis() = ga_instruction_set();
567 if (pwi.workspace().nb_trees()) {
568 ga_tree tree = *(pwi.workspace().tree_info(0).ptree);
569 ga_derivative(tree, pwi.workspace(), source_mesh,
570 var.varname, var.transname, 1);
572 ga_semantic_analysis(tree, local_workspace, *((
const mesh *)(0)),
574 ga_compile_interpolation(pwi.workspace(), pwi.gis());
579 if (recompute_elt_boxes) {
581 element_boxes.clear();
582 base_node bmin(N), bmax(N);
583 const dal::bit_vector&
587 for (dal::bv_visitor cv(convex_index); !cv.finished(); ++cv) {
593 gmm::copy(target_mesh.points_of_convex(cv)[0], bmin);
594 gmm::copy(bmin, bmax);
601 const base_node &pt = target_mesh.points_of_convex(cv)[ip];
604 bmin[k] = std::min(bmin[k], pt[k]);
605 bmax[k] = std::max(bmax[k], pt[k]);
609 scalar_type h = bmax[0] - bmin[0];
610 for (
size_type k = 1; k < N; ++k) h = std::max(h, bmax[k]-bmin[k]);
611 if (pgt->is_linear()) h *= 1E-4;
612 for (
auto &&val : bmin) val -= h*0.2;
613 for (
auto &&val : bmax) val += h*0.2;
615 element_boxes.add_box(bmin, bmax, cv);
617 recompute_elt_boxes =
false;
621 void finalize()
const {
622 for (
const std::string &transname : local_gis.transformations)
623 local_workspace.interpolate_transformation(transname)->finalize();
624 local_gis = ga_instruction_set();
627 std::string expression()
const {
return expr; }
629 int transform(
const ga_workspace &,
const mesh &m,
631 const base_small_vector &Normal,
636 std::map<var_trans_pair, base_tensor> &derivatives,
637 bool compute_derivatives)
const {
640 ga_interpolation_single_point_exec(local_gis, local_workspace, ctx_x,
643 GMM_ASSERT1(local_workspace.assembled_tensor().size() == m.dim(),
644 "Wrong dimension of the transformation expression");
646 gmm::copy(local_workspace.assembled_tensor().as_vector(), P);
648 bgeot::rtree::pbox_set bset;
649 element_boxes.find_boxes_at_point(P, bset);
652 while (bset.size()) {
653 bgeot::rtree::pbox_set::iterator it = bset.begin(), itmax = it;
655 if (bset.size() > 1) {
657 scalar_type rate_max = scalar_type(-1);
658 for (; it != bset.end(); ++it) {
660 scalar_type rate_box = scalar_type(1);
661 for (
size_type i = 0; i < m.dim(); ++i) {
662 scalar_type h = (*it)->max[i] - (*it)->min[i];
663 if (h > scalar_type(0)) {
665 = std::min((*it)->max[i] - P[i], P[i] - (*it)->min[i]) / h;
666 rate_box = std::min(rate, rate_box);
669 if (rate_box > rate_max) {
677 gic.init(target_mesh.points_of_convex(cv),
678 target_mesh.trans_of_convex(cv));
680 bool converged =
true;
681 bool is_in = gic.
invert(P, P_ref, converged, 1E-4);
688 if (is_in && converged) {
694 if (bset.size() == 1)
break;
703 if (compute_derivatives) {
704 for (
auto &&d : derivatives) {
705 workspace_gis_pair &pwi = compiled_derivatives[d.first];
707 gmm::clear(pwi.workspace().assembled_tensor().as_vector());
708 ga_function_exec(pwi.gis());
709 d.second = pwi.workspace().assembled_tensor();
715 interpolate_transformation_expression
717 : source_mesh(sm), target_mesh(tm), target_region(trg), expr(expr_),
718 recompute_elt_boxes(
true), extract_variable_done(
false),
719 extract_data_done(
false)
720 { this->add_dependency(tm); }
726 (ga_workspace &workspace,
const std::string &name,
const mesh &sm,
727 const mesh &tm,
const std::string &expr) {
729 (workspace, name, sm, tm,
size_type(-1), expr);
733 (ga_workspace &workspace,
const std::string &name,
const mesh &sm,
735 pinterpolate_transformation
736 p = std::make_shared<interpolate_transformation_expression>
738 workspace.add_interpolate_transformation(name, p);
742 (
model &md,
const std::string &name,
const mesh &sm,
const mesh &tm,
743 const std::string &expr) {
749 (
model &md,
const std::string &name,
const mesh &sm,
const mesh &tm,
750 size_type trg,
const std::string &expr) {
751 pinterpolate_transformation
752 p = std::make_shared<interpolate_transformation_expression>
761 class interpolate_transformation_neighbour
765 void update_from_context()
const {}
766 void extract_variables(
const ga_workspace &,
767 std::set<var_trans_pair> &,
769 const std::string &)
const {}
770 void init(
const ga_workspace &)
const {}
771 void finalize()
const {}
773 std::string expression(
void)
const {
return "X"; }
775 int transform(
const ga_workspace &,
const mesh &m_x,
777 const base_small_vector &,
const mesh **m_t,
781 std::map<var_trans_pair, base_tensor> &,
782 bool compute_derivatives)
const {
788 GMM_ASSERT1(face_x !=
short_type(-1),
"Neighbour transformation can " 789 "only be applied to internal faces");
795 gic.init(m_x.points_of_convex(adj_face.cv),
796 m_x.trans_of_convex(adj_face.cv));
797 bool converged =
true;
798 bool is_in = gic.
invert(ctx_x.
xreal(), P_ref, converged, 1E-4);
799 GMM_ASSERT1(is_in && converged,
"Geometric transformation inversion " 800 "has failed in neighbour transformation");
801 face_num = adj_face.f;
805 GMM_ASSERT1(!compute_derivatives,
806 "No derivative for this transformation");
810 interpolate_transformation_neighbour() { }
816 return (std::make_shared<interpolate_transformation_neighbour>());
823 class interpolate_transformation_element_extrapolation
827 std::map<size_type, size_type> elt_corr;
830 void update_from_context()
const {}
831 void extract_variables(
const ga_workspace &,
832 std::set<var_trans_pair> &,
834 const std::string &)
const {}
835 void init(
const ga_workspace &)
const {}
836 void finalize()
const {}
837 std::string expression(
void)
const {
return "X"; }
839 int transform(
const ga_workspace &,
const mesh &m_x,
841 const base_small_vector &,
const mesh **m_t,
845 std::map<var_trans_pair, base_tensor> &,
846 bool compute_derivatives)
const {
849 GMM_ASSERT1(&sm == &m_x,
"Bad mesh");
851 auto it = elt_corr.find(cv_x);
852 if (it != elt_corr.end()) cv_y = it->second;
856 gic.init(m_x.points_of_convex(cv_y),
857 m_x.trans_of_convex(cv_y));
858 bool converged =
true;
860 GMM_ASSERT1(converged,
"Geometric transformation inversion " 861 "has failed in element extrapolation transformation");
868 P_ref = ctx_x.
xref();
871 GMM_ASSERT1(!compute_derivatives,
872 "No derivative for this transformation");
876 void set_correspondance(
const std::map<size_type, size_type> &ec) {
880 interpolate_transformation_element_extrapolation
881 (
const mesh &sm_,
const std::map<size_type, size_type> &ec)
882 : sm(sm_), elt_corr(ec) { }
886 void add_element_extrapolation_transformation
887 (
model &md,
const std::string &name,
const mesh &sm,
888 std::map<size_type, size_type> &elt_corr) {
889 pinterpolate_transformation
890 p = std::make_shared<interpolate_transformation_element_extrapolation>
895 void add_element_extrapolation_transformation
896 (ga_workspace &workspace,
const std::string &name,
const mesh &sm,
897 std::map<size_type, size_type> &elt_corr) {
898 pinterpolate_transformation
899 p = std::make_shared<interpolate_transformation_element_extrapolation>
901 workspace.add_interpolate_transformation(name, p);
904 void set_element_extrapolation_correspondance
905 (ga_workspace &workspace,
const std::string &name,
906 std::map<size_type, size_type> &elt_corr) {
907 GMM_ASSERT1(workspace.interpolate_transformation_exists(name),
908 "Unknown transformation");
909 auto pit=workspace.interpolate_transformation(name).get();
911 =
dynamic_cast<const interpolate_transformation_element_extrapolation *
> 914 "The transformation is not of element extrapolation type");
915 const_cast<interpolate_transformation_element_extrapolation *
>(cpext)
916 ->set_correspondance(elt_corr);
919 void set_element_extrapolation_correspondance
920 (
model &md,
const std::string &name,
921 std::map<size_type, size_type> &elt_corr) {
923 "Unknown transformation");
926 =
dynamic_cast<const interpolate_transformation_element_extrapolation *
> 929 "The transformation is not of element extrapolation type");
930 const_cast<interpolate_transformation_element_extrapolation *
>(cpext)
931 ->set_correspondance(elt_corr);
virtual size_type nb_basic_dof_of_element(size_type cv) const
Return the number of degrees of freedom attached to a given convex.
virtual size_type nb_dof() const
Return the total number of degrees of freedom.
structure used to hold a set of convexes and/or convex faces.
size_type nb_dof() const
Total number of degrees of freedom in the model.
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...
short_type face_num() const
get the current face number
static mesh_region all_convexes()
provide a default value for the mesh_region parameters of assembly procedures etc.
const dal::bit_vector & index() const
Index of the region convexes, or the convexes from the partition on the current thread.
size_type convex_num() const
get the current convex number
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) ...
Describe a mesh (collection of convexes (elements) and points).
pinterpolate_transformation interpolate_transformation_neighbour_instance()
Create a new instance of a transformation corresponding to the interpolation on the neighbour element...
convex_face adjacent_face(size_type ic, short_type f) const
Return a face of the neighbouring element that is adjacent to the given face.
does the inversion of the geometric transformation for a given convex
Describe an integration method linked to a mesh.
Semantic analysis of assembly trees and semantic manipulations.
const base_node & xreal() const
coordinates of the current point, in the real convex.
pinterpolate_transformation interpolate_transformation(const std::string &name) const
Get a pointer to the interpolate transformation name.
``Model'' variables store the variables, the data and the description of a model. ...
size_t size_type
used as the common size type in the library
bool interpolate_transformation_exists(const std::string &name) const
Tests if name correpsonds to an interpolate transformation.
bool invert(const base_node &n, base_node &n_ref, scalar_type IN_EPS=1e-12)
given the node on the real element, returns the node on the reference element (even if it is outside ...
structure passed as the argument of fem interpolation functions.
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.
"iterator" class for regions.
void add_interpolate_transformation(const std::string &name, pinterpolate_transformation ptrans)
Add a interpolate transformation to the model to be used with the generic assembly.
Deal with interdependencies of objects.
GEneric Tool for Finite Element Methods.
std::shared_ptr< const getfem::virtual_fem > pfem
type of pointer on a fem description
Describe a finite element method linked to a mesh.
void clear(L &l)
clear (fill with zeros) a vector or matrix.
gmm::uint16_type short_type
used as the common short type integer in the library
const base_node & xref() const
coordinates of the current point, in the reference convex.
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.
const mesh_region region(size_type id) const
Return the region of index 'id'.
Compilation and execution operations.
virtual ind_dof_ct ind_basic_dof_of_element(size_type cv) const
Give an array of the dof numbers a of convex.
Balanced tree of n-dimensional rectangles.
void resize(M &v, size_type m, size_type n)
*/
const dal::bit_vector & convex_index() const
Return the list of valid convex IDs.
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation