28 struct contact_boundary {
 
   33   mutable const model_real_plain_vector *U;
 
   34   mutable model_real_plain_vector U_unred; 
 
   36   contact_boundary(
size_type r, 
const mesh_fem *mf, 
const std::string &dn)
 
   37     : region(r), mfu(mf), dispname(dn)
 
   42 base_small_vector element_U(
const contact_boundary &cb, 
size_type cv)
 
   44   auto U_elm = base_small_vector{};
 
   50 auto most_central_box(
const bgeot::rtree::pbox_set &bset,
 
   55   auto itmax = begin(bset);
 
   58   if (bset.size() > 1) {
 
   59     auto rate_max = scalar_type{-1};
 
   60     for (; it != end(bset); ++it) {
 
   61       auto rate_box = scalar_type{1};
 
   62       for (
size_type i = 0; i < pt.size(); ++i) {
 
   63         auto h = (*it)->max->at(i) - (*it)->min->at(i);
 
   65           auto rate = min((*it)->max->at(i) - pt[i], pt[i] - (*it)->min->at(i)) / h;
 
   66           rate_box = min(rate, rate_box);
 
   69       if (rate_box > rate_max) {
 
   81 class  interpolate_transformation_on_deformed_domains
 
   82   : 
public virtual_interpolate_transformation {
 
   84   contact_boundary master;
 
   85   contact_boundary slave; 
 
   88   mutable std::map<size_type, std::vector<size_type>> box_to_convex; 
 
   91   mutable fem_precomp_pool fppool;
 
   94   void compute_element_boxes()
 const { 
 
   96     model_real_plain_vector Uelm; 
 
   97     element_boxes.clear();
 
   99     auto bnum = master.region;
 
  100     auto &mfu = *(master.mfu);
 
  101     auto &U   = *(master.U);
 
  102     auto &m   = mfu.linked_mesh();
 
  105     base_node Xdeformed(N), bmin(N), bmax(N);
 
  106     auto region = m.region(bnum);
 
  111     region.prohibit_partitioning();
 
  113     GMM_ASSERT1(mfu.get_qdim() == N, 
"Wrong mesh_fem qdim");
 
  115     dal::bit_vector points_already_interpolated;
 
  116     std::vector<base_node> transformed_points(m.nb_max_points());
 
  117     box_to_convex.clear();
 
  121       auto pgt  = m.trans_of_convex(cv);
 
  122       auto pf_s = mfu.fem_of_element(cv);
 
  123       auto pfp  = fppool(pf_s, pgt->pgeometric_nodes());
 
  126       mfu.linked_mesh().points_of_convex(cv, G);
 
  128       auto ctx   = fem_interpolation_context{pgt, pfp, 
size_type(-1), G, cv};
 
  129       auto nb_pt = pgt->structure()->nb_points();
 
  132         auto ind = m.ind_points_of_convex(cv)[k];
 
  136         if (points_already_interpolated.is_in(ind)) {
 
  137           Xdeformed = transformed_points[ind];
 
  139           pf_s->interpolation(ctx, Uelm, Xdeformed, dim_type{N});
 
  140           Xdeformed += ctx.xreal(); 
 
  141           transformed_points[ind] = Xdeformed;
 
  142           points_already_interpolated.add(ind);
 
  146           bmin = bmax = Xdeformed;
 
  149             bmin[l] = std::min(bmin[l], Xdeformed[l]);
 
  150             bmax[l] = std::max(bmax[l], Xdeformed[l]);
 
  156       box_to_convex[element_boxes.add_box(bmin, bmax)].push_back(cv);
 
  158     element_boxes.build_tree();
 
  161   fem_interpolation_context deformed_master_context(
size_type cv)
 const 
  163     auto &mfu  = *(master.mfu);
 
  164     auto G     = base_matrix{};
 
  165     auto pfu   = mfu.fem_of_element(cv);
 
  166     auto pgt   = master.mfu->linked_mesh().trans_of_convex(cv);
 
  167     auto pfp   = fppool(pfu, pgt->pgeometric_nodes());
 
  168     master.mfu->linked_mesh().points_of_convex(cv, G);
 
  172   std::vector<bgeot::base_node> deformed_master_nodes(
size_type cv)
 const {
 
  173     using namespace bgeot;
 
  176     auto nodes = vector<base_node>{};
 
  178     auto U_elm = element_U(master, cv);
 
  179     auto &mfu  = *(master.mfu);
 
  180     auto G     = base_matrix{};
 
  181     auto pfu   = mfu.fem_of_element(cv);
 
  182     auto pgt   = master.mfu->linked_mesh().trans_of_convex(cv);
 
  183     auto pfp   = fppool(pfu, pgt->pgeometric_nodes());
 
  184     auto N     = mfu.linked_mesh().dim();
 
  185     auto pt    = base_node(N);
 
  187     master.mfu->linked_mesh().points_of_convex(cv, G);
 
  188     auto ctx = fem_interpolation_context{pgt, pfp, 
size_type(-1), G, cv};
 
  189     auto nb_pt = pgt->structure()->nb_points();
 
  190     nodes.reserve(nb_pt);
 
  193       pfu->interpolation(ctx, U_elm, U, dim_type{N});
 
  203   interpolate_transformation_on_deformed_domains(
 
  206     const std::string      &source_displacements,
 
  209     const std::string      &target_displacements)
 
  211       slave{source_region, &mf_source, source_displacements},
 
  212       master{target_region, &mf_target, target_displacements}
 
  216   void extract_variables(
const ga_workspace           &workspace,
 
  217                          std::set<var_trans_pair>     &vars,
 
  220                          const std::string            &interpolate_name)
 const override {
 
  221     if (!ignore_data || !(workspace.is_constant(master.dispname))){
 
  222       vars.emplace(master.dispname, interpolate_name);
 
  223       vars.emplace(slave.dispname, 
"");
 
  227   void init(
const ga_workspace &workspace)
 const override {
 
  229     for (
auto pcb : std::list<const contact_boundary*>{&master, &slave}) {
 
  230       auto &mfu = *(pcb->mfu);
 
  231       if (mfu.is_reduced()) {
 
  233         mfu.extend_vector(workspace.value(pcb->dispname), pcb->U_unred);
 
  234         pcb->U = &(pcb->U_unred);
 
  236         pcb->U = &(workspace.value(pcb->dispname));
 
  239     compute_element_boxes();
 
  242   void finalize()
 const override {
 
  243     element_boxes.clear();
 
  244     box_to_convex.clear();
 
  245     master.U_unred.clear();
 
  246     slave.U_unred.clear();
 
  250   int transform(
const ga_workspace                    &workspace,
 
  252                 fem_interpolation_context             &ctx_x,
 
  259                 std::map<var_trans_pair, base_tensor> &derivatives,
 
  260                 bool                                  compute_derivatives)
 const override {
 
  262     auto &target_mesh = master.mfu->linked_mesh();
 
  264     auto transformation_success = 
false;
 
  267     using namespace bgeot;
 
  271     auto cv_x    = ctx_x.convex_num();
 
  272     auto U_elm_x = element_U(slave, cv_x);
 
  273     auto &mfu_x  = *(slave.mfu);
 
  274     auto pfu_x   = mfu_x.fem_of_element(cv_x);
 
  275     auto N       = mfu_x.linked_mesh().dim();
 
  277     auto G_x     = base_matrix{}; 
 
  278     m_x.points_of_convex(cv_x, G_x);
 
  280     pfu_x->interpolation(ctx_x, U_elm_x, U_x, dim_type{N});
 
  282     add(ctx_x.xreal(), U_x, pt_x);
 
  290     auto bset = rtree::pbox_set{};
 
  291     element_boxes.find_boxes_at_point(pt_x, bset);
 
  292     while (!bset.empty())
 
  294       auto itmax = most_central_box(bset, pt_x);
 
  296       for (
auto i : box_to_convex.at((*itmax)->id))
 
  298         auto deformed_nodes_y = deformed_master_nodes(i);
 
  299         gic.init(deformed_nodes_y, target_mesh.trans_of_convex(i));
 
  300         auto converged = 
true;
 
  301         auto is_in = gic.
invert(pt_x, P_ref, converged);
 
  302         if (is_in && converged) {
 
  305           transformation_success = 
true;
 
  309       if (transformation_success || (bset.size() == 1)) 
break;
 
  318     if (compute_derivatives && transformation_success) {
 
  319       GMM_ASSERT2(derivatives.size() == 2,
 
  320                   "Expecting to return derivatives only for Umaster and Uslave");
 
  322       for (
auto &pair : derivatives)
 
  324         if (pair.first.varname == slave.dispname)
 
  326           auto base_ux = base_tensor{};
 
  327           auto vbase_ux = base_matrix{} ;
 
  328           ctx_x.base_value(base_ux);
 
  329           auto qdim_ux = pfu_x->target_dim();
 
  330           auto ndof_ux = pfu_x->nb_dof(cv_x) * N / qdim_ux;
 
  331           vectorize_base_tensor(base_ux, vbase_ux, ndof_ux, qdim_ux, N);
 
  332           pair.second.adjust_sizes(ndof_ux, N);
 
  333           copy(vbase_ux.as_vector(), pair.second.as_vector());
 
  336         if (pair.first.varname == master.dispname)
 
  338           auto ctx_y = deformed_master_context(cv);
 
  339           ctx_y.set_xref(P_ref);
 
  340           auto base_uy = base_tensor{};
 
  341           auto vbase_uy = base_matrix{} ;
 
  342           ctx_y.base_value(base_uy);
 
  343           auto pfu_y   = master.mfu->fem_of_element(cv);
 
  344           auto dim_y = master.mfu->linked_mesh().dim();
 
  345           auto qdim_uy = pfu_y->target_dim();
 
  346           auto ndof_uy = pfu_y->nb_dof(cv) * dim_y / qdim_uy;
 
  347           vectorize_base_tensor(base_uy, vbase_uy, ndof_uy, qdim_uy, dim_y);
 
  348           pair.second.adjust_sizes(ndof_uy, dim_y);
 
  349           copy(vbase_uy.as_vector(), pair.second.as_vector());
 
  350           scale(pair.second.as_vector(), -1.);
 
  352         else GMM_ASSERT2(
false, 
"unexpected derivative variable");
 
  356     return transformation_success ? 1 : 0;
 
  362   (ga_workspace &workspace, 
const std::string &transname,
 
  363    const mesh &source_mesh, 
const std::string &source_displacements,
 
  365    const std::string &target_displacements, 
const mesh_region &target_region)
 
  367     auto pmf_source = workspace.associated_mf(source_displacements);
 
  368     auto pmf_target = workspace.associated_mf(target_displacements);
 
  369     auto p_transformation
 
  370       = std::make_shared<interpolate_transformation_on_deformed_domains>(source_region.id(),
 
  372                                                                          source_displacements,
 
  375                                                                          target_displacements);
 
  376     workspace.add_interpolate_transformation(transname, p_transformation);
 
  380   (
model &md, 
const std::string &transname,
 
  381    const mesh &source_mesh, 
const std::string &source_displacements,
 
  383    const std::string &target_displacements, 
const mesh_region &target_region)
 
  387     auto p_transformation
 
  388       = std::make_shared<interpolate_transformation_on_deformed_domains>(source_region.id(),
 
  390                                                                          source_displacements,
 
  393                                                                          target_displacements);
 
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.
 
"iterator" class for regions.
 
structure used to hold a set of convexes and/or convex faces.
 
Describe a mesh (collection of convexes (elements) and points).
 
`‘Model’' variables store the variables, the data and the description of a model.
 
const mesh_fem & mesh_fem_of_variable(const std::string &name) const
Gives the access to the mesh_fem of a variable if any.
 
void add_interpolate_transformation(const std::string &name, pinterpolate_transformation ptrans)
Add an interpolate transformation to the model to be used with the generic assembly.
 
A language for generic assembly of pde boundary value problems.
 
Model representation in Getfem.
 
void copy(const L1 &l1, L2 &l2)
*/
 
void resize(V &v, size_type n)
*/
 
void add(const L1 &l1, L2 &l2)
*/
 
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
 
GEneric Tool for Finite Element Methods.
 
void add_interpolate_transformation_on_deformed_domains(ga_workspace &workspace, const std::string &transname, const mesh &source_mesh, const std::string &source_displacements, const mesh_region &source_region, const mesh &target_mesh, const std::string &target_displacements, const mesh_region &target_region)
Add a transformation to the workspace that creates an identity mapping between two meshes in deformed...
 
void slice_vector_on_basic_dof_of_element(const mesh_fem &mf, const VEC1 &vec, size_type cv, VEC2 &coeff, size_type qmult1=size_type(-1), size_type qmult2=size_type(-1))
Given a mesh_fem.