GetFEM  5.5
getfem_interpolation_on_deformed_domains.cc
1 /*===========================================================================
2 
3  Copyright (C) 2013-2026 Yves Renard, Konstantinos Poulios and Andriy Andreykiv.
4 
5  This file is a part of GetFEM
6 
7  GetFEM is free software; you can redistribute it and/or modify it
8  under the terms of the GNU Lesser General Public License as published
9  by the Free Software Foundation; either version 3 of the License, or
10  (at your option) any later version along with the GCC Runtime Library
11  Exception either version 3.1 or (at your option) any later version.
12  This program is distributed in the hope that it will be useful, but
13  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
14  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
15  License and GCC Runtime Library Exception for more details.
16  You should have received a copy of the GNU Lesser General Public License
17  along with this program. If not, see https://www.gnu.org/licenses/.
18 
19 ===========================================================================*/
20 
22 #include "getfem/getfem_models.h"
23 
24 namespace getfem {
25 
26 // Structure describing a contact boundary (or contact body)
27 struct contact_boundary {
28  size_type region; // boundary region for the slave (source)
29  // and volume region for the master (target)
30  const getfem::mesh_fem *mfu; // F.e.m. for the displacement.
31  std::string dispname; // Variable name for the displacement
32  mutable const model_real_plain_vector *U;// Displacement
33  mutable model_real_plain_vector U_unred; // Unreduced displacement
34 
35  contact_boundary(size_type r, const mesh_fem *mf, const std::string &dn)
36  : region(r), mfu(mf), dispname(dn)
37  {}
38 };
39 
40 //extract element displacements from a contact boundary object
41 base_small_vector element_U(const contact_boundary &cb, size_type cv)
42 {
43  auto U_elm = base_small_vector{};
44  slice_vector_on_basic_dof_of_element(*(cb.mfu), *cb.U, cv, U_elm);
45  return U_elm;
46 }
47 
48 //Returns an iterator of a box which centre is closest to the given point
49 auto most_central_box(const bgeot::rtree::pbox_set &bset,
50  const bgeot::base_node &pt) -> decltype(begin(bset))
51 {
52  using namespace std;
53 
54  auto itmax = begin(bset);
55 
56  auto it = itmax;
57  if (bset.size() > 1) {
58  auto rate_max = scalar_type{-1};
59  for (; it != end(bset); ++it) {
60  auto rate_box = scalar_type{1};
61  for (size_type i = 0; i < pt.size(); ++i) {
62  auto h = (*it)->max->at(i) - (*it)->min->at(i);
63  if (h > 0.) {
64  auto rate = min((*it)->max->at(i) - pt[i], pt[i] - (*it)->min->at(i)) / h;
65  rate_box = min(rate, rate_box);
66  }
67  }
68  if (rate_box > rate_max) {
69  itmax = it;
70  rate_max = rate_box;
71  }
72  }
73  }
74 
75  return itmax;
76 }
77 
78 //Transformation that creates identity mapping between two contact boundaries,
79 //deformed with provided displacement fields
80 class interpolate_transformation_on_deformed_domains
81  : public virtual_interpolate_transformation {
82 
83  contact_boundary master;//also marked with a target or Y prefix/suffix
84  contact_boundary slave; //also marked with a source or X prefix/suffix
85 
86  mutable bgeot::rtree element_boxes;
87  mutable std::map<size_type, std::vector<size_type>> box_to_convex; //index to obtain a convex
88  //number from a box number
89  mutable bgeot::geotrans_inv_convex gic;
90  mutable fem_precomp_pool fppool;
91 
92  //Create a box tree based on the deformed elements of the master (target)
93  void compute_element_boxes() const { // called by init
94  base_matrix G;
95  model_real_plain_vector Uelm; //element displacement
96  element_boxes.clear();
97 
98  auto bnum = master.region;
99  auto &mfu = *(master.mfu);
100  auto &U = *(master.U);
101  auto &m = mfu.linked_mesh();
102  auto N = m.dim();
103 
104  base_node Xdeformed(N), bmin(N), bmax(N);
105  auto region = m.region(bnum);
106 
107  //the box tree creation and subsequent transformation inversion
108  //should be done for all elements of the master, while integration
109  //will be performed only on a thread partition of the slave
110  region.prohibit_partitioning();
111 
112  GMM_ASSERT1(mfu.get_qdim() == N, "Wrong mesh_fem qdim");
113 
114  dal::bit_vector points_already_interpolated;
115  std::vector<base_node> transformed_points(m.nb_max_points());
116  box_to_convex.clear();
117 
118  for (getfem::mr_visitor v(region, m); !v.finished(); ++v) {
119  auto cv = v.cv();
120  auto pgt = m.trans_of_convex(cv);
121  auto pf_s = mfu.fem_of_element(cv);
122  auto pfp = fppool(pf_s, pgt->pgeometric_nodes());
123 
124  slice_vector_on_basic_dof_of_element(mfu, U, cv, Uelm);
125  mfu.linked_mesh().points_of_convex(cv, G);
126 
127  auto ctx = fem_interpolation_context{pgt, pfp, size_type(-1), G, cv};
128  auto nb_pt = pgt->structure()->nb_points();
129 
130  for (size_type k = 0; k < nb_pt; ++k) {
131  auto ind = m.ind_points_of_convex(cv)[k];
132 
133  // computation of a transformed vertex
134  ctx.set_ii(k);
135  if (points_already_interpolated.is_in(ind)) {
136  Xdeformed = transformed_points[ind];
137  } else {
138  pf_s->interpolation(ctx, Uelm, Xdeformed, dim_type{N});
139  Xdeformed += ctx.xreal(); //Xdeformed = U + Xo
140  transformed_points[ind] = Xdeformed;
141  points_already_interpolated.add(ind);
142  }
143 
144  if (k == 0) // computation of bounding box
145  bmin = bmax = Xdeformed;
146  else {
147  for (size_type l = 0; l < N; ++l) {
148  bmin[l] = std::min(bmin[l], Xdeformed[l]);
149  bmax[l] = std::max(bmax[l], Xdeformed[l]);
150  }
151  }
152  }
153 
154  // Store the bounding box and additional information.
155  box_to_convex[element_boxes.add_box(bmin, bmax)].push_back(cv);
156  }
157  element_boxes.build_tree();
158  }
159 
160  fem_interpolation_context deformed_master_context(size_type cv) const
161  {
162  auto &mfu = *(master.mfu);
163  auto G = base_matrix{};
164  auto pfu = mfu.fem_of_element(cv);
165  auto pgt = master.mfu->linked_mesh().trans_of_convex(cv);
166  auto pfp = fppool(pfu, pgt->pgeometric_nodes());
167  master.mfu->linked_mesh().points_of_convex(cv, G);
168  return {pgt, pfp, size_type(-1), G, cv};
169  }
170 
171  std::vector<bgeot::base_node> deformed_master_nodes(size_type cv) const {
172  using namespace bgeot;
173  using namespace std;
174 
175  auto nodes = vector<base_node>{};
176 
177  auto U_elm = element_U(master, cv);
178  auto &mfu = *(master.mfu);
179  auto G = base_matrix{};
180  auto pfu = mfu.fem_of_element(cv);
181  auto pgt = master.mfu->linked_mesh().trans_of_convex(cv);
182  auto pfp = fppool(pfu, pgt->pgeometric_nodes());
183  auto N = mfu.linked_mesh().dim();
184  auto pt = base_node(N);
185  auto U = base_small_vector(N);
186  master.mfu->linked_mesh().points_of_convex(cv, G);
187  auto ctx = fem_interpolation_context{pgt, pfp, size_type(-1), G, cv};
188  auto nb_pt = pgt->structure()->nb_points();
189  nodes.reserve(nb_pt);
190  for (size_type k = 0; k < nb_pt; ++k) {
191  ctx.set_ii(k);
192  pfu->interpolation(ctx, U_elm, U, dim_type{N});
193  gmm::add(ctx.xreal(), U, pt);
194  nodes.push_back(pt);
195  }
196 
197  return nodes;
198  }
199 
200 public:
201 
202  interpolate_transformation_on_deformed_domains(
203  size_type source_region,
204  const getfem::mesh_fem &mf_source,
205  const std::string &source_displacements,
206  size_type target_region,
207  const getfem::mesh_fem &mf_target,
208  const std::string &target_displacements)
209  :
210  slave{source_region, &mf_source, source_displacements},
211  master{target_region, &mf_target, target_displacements}
212 {}
213 
214 
215  void extract_variables(const ga_workspace &workspace,
216  std::set<var_trans_pair> &vars,
217  bool ignore_data,
218  const mesh &m_x,
219  const std::string &interpolate_name) const override {
220  if (!ignore_data || !(workspace.is_constant(master.dispname))){
221  vars.emplace(master.dispname, interpolate_name);
222  vars.emplace(slave.dispname, "");
223  }
224  }
225 
226  void init(const ga_workspace &workspace) const override {
227 
228  for (auto pcb : std::list<const contact_boundary*>{&master, &slave}) {
229  auto &mfu = *(pcb->mfu);
230  if (mfu.is_reduced()) {
231  gmm::resize(pcb->U_unred, mfu.nb_basic_dof());
232  mfu.extend_vector(workspace.value(pcb->dispname), pcb->U_unred);
233  pcb->U = &(pcb->U_unred);
234  } else {
235  pcb->U = &(workspace.value(pcb->dispname));
236  }
237  }
238  compute_element_boxes();
239  };
240 
241  void finalize() const override {
242  element_boxes.clear();
243  box_to_convex.clear();
244  master.U_unred.clear();
245  slave.U_unred.clear();
246  fppool.clear();
247  }
248 
249  int transform(const ga_workspace &workspace,
250  const mesh &m_x,
251  fem_interpolation_context &ctx_x,
252  const base_small_vector &/*Normal*/,
253  const mesh **m_t,
254  size_type &cv,
255  short_type &face_num,
256  base_node &P_ref,
257  base_small_vector &N_y,
258  std::map<var_trans_pair, base_tensor> &derivatives,
259  bool compute_derivatives) const override {
260 
261  auto &target_mesh = master.mfu->linked_mesh();
262  *m_t = &target_mesh;
263  auto transformation_success = false;
264 
265  using namespace gmm;
266  using namespace bgeot;
267  using namespace std;
268 
269  //compute a deformed point of the slave
270  auto cv_x = ctx_x.convex_num();
271  auto U_elm_x = element_U(slave, cv_x);
272  auto &mfu_x = *(slave.mfu);
273  auto pfu_x = mfu_x.fem_of_element(cv_x);
274  auto N = mfu_x.linked_mesh().dim();
275  auto U_x = base_small_vector(N);
276  auto G_x = base_matrix{}; //coordinates of the source element nodes
277  m_x.points_of_convex(cv_x, G_x);
278  ctx_x.set_pf(pfu_x);
279  pfu_x->interpolation(ctx_x, U_elm_x, U_x, dim_type{N});
280  auto pt_x = base_small_vector(N); //deformed point of the slave
281  add(ctx_x.xreal(), U_x, pt_x);
282 
283  //Find the best box from the master (target) that
284  //corresponds to this point (The box which centre is the closest to the point).
285  //Obtain the corresponding element number using the box id and box_to_convex
286  //indices. Compute deformed nodes of the target element. Invert the geometric
287  //transformation of the target element with deformed nodes, obtaining this way
288  //reference coordinates of the target element
289  auto bset = rtree::pbox_set{};
290  element_boxes.find_boxes_at_point(pt_x, bset);
291  while (!bset.empty())
292  {
293  auto itmax = most_central_box(bset, pt_x);
294 
295  for (auto i : box_to_convex.at((*itmax)->id))
296  {
297  auto deformed_nodes_y = deformed_master_nodes(i);
298  gic.init(deformed_nodes_y, target_mesh.trans_of_convex(i));
299  auto converged = true;
300  auto is_in = gic.invert(pt_x, P_ref, converged);
301  if (is_in && converged) {
302  cv = i;
303  face_num = static_cast<short_type>(-1);
304  transformation_success = true;
305  break;
306  }
307  }
308  if (transformation_success || (bset.size() == 1)) break;
309  bset.erase(itmax);
310  }
311 
312  //Since this transformation can be seen as Xsource + Usource - Utarget,
313  //the corresponding stiffnesses are identity matrix for Usource and
314  //minus identity for Utarget. The required answer in this function is
315  //stiffness X shape function. Hence, returning shape function for Usource
316  //and min shape function for Utarget
317  if (compute_derivatives && transformation_success) {
318  GMM_ASSERT2(derivatives.size() == 2,
319  "Expecting to return derivatives only for Umaster and Uslave");
320 
321  for (auto &pair : derivatives)
322  {
323  if (pair.first.varname == slave.dispname)
324  {
325  auto base_ux = base_tensor{};
326  auto vbase_ux = base_matrix{} ;
327  ctx_x.base_value(base_ux);
328  auto qdim_ux = pfu_x->target_dim();
329  auto ndof_ux = pfu_x->nb_dof(cv_x) * N / qdim_ux;
330  vectorize_base_tensor(base_ux, vbase_ux, ndof_ux, qdim_ux, N);
331  pair.second.adjust_sizes(ndof_ux, N);
332  copy(vbase_ux.as_vector(), pair.second.as_vector());
333  }
334  else
335  if (pair.first.varname == master.dispname)
336  {
337  auto ctx_y = deformed_master_context(cv);
338  ctx_y.set_xref(P_ref);
339  auto base_uy = base_tensor{};
340  auto vbase_uy = base_matrix{} ;
341  ctx_y.base_value(base_uy);
342  auto pfu_y = master.mfu->fem_of_element(cv);
343  auto dim_y = master.mfu->linked_mesh().dim();
344  auto qdim_uy = pfu_y->target_dim();
345  auto ndof_uy = pfu_y->nb_dof(cv) * dim_y / qdim_uy;
346  vectorize_base_tensor(base_uy, vbase_uy, ndof_uy, qdim_uy, dim_y);
347  pair.second.adjust_sizes(ndof_uy, dim_y);
348  copy(vbase_uy.as_vector(), pair.second.as_vector());
349  scale(pair.second.as_vector(), -1.);
350  }
351  else GMM_ASSERT2(false, "unexpected derivative variable");
352  }
353  }
354 
355  return transformation_success ? 1 : 0;
356  }
357 
358 };
359 
361  (ga_workspace &workspace, const std::string &transname,
362  const mesh &source_mesh, const std::string &source_displacements,
363  const mesh_region &source_region, const mesh &target_mesh,
364  const std::string &target_displacements, const mesh_region &target_region)
365  {
366  auto pmf_source = workspace.associated_mf(source_displacements);
367  auto pmf_target = workspace.associated_mf(target_displacements);
368  auto p_transformation
369  = std::make_shared<interpolate_transformation_on_deformed_domains>(source_region.id(),
370  *pmf_source,
371  source_displacements,
372  target_region.id(),
373  *pmf_target,
374  target_displacements);
375  workspace.add_interpolate_transformation(transname, p_transformation);
376  }
377 
379  (model &md, const std::string &transname,
380  const mesh &source_mesh, const std::string &source_displacements,
381  const mesh_region &source_region, const mesh &target_mesh,
382  const std::string &target_displacements, const mesh_region &target_region)
383  {
384  auto &mf_source = md.mesh_fem_of_variable(source_displacements);
385  auto mf_target = md.mesh_fem_of_variable(target_displacements);
386  auto p_transformation
387  = std::make_shared<interpolate_transformation_on_deformed_domains>(source_region.id(),
388  mf_source,
389  source_displacements,
390  target_region.id(),
391  mf_target,
392  target_displacements);
393  md.add_interpolate_transformation(transname, p_transformation);
394  }
395 
396 } /* end of namespace getfem. */
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.
Definition: bgeot_rtree.h:97
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).
Definition: getfem_mesh.h:98
`‘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)
*‍/
Definition: gmm_blas.h:976
void resize(V &v, size_type n)
*‍/
Definition: gmm_blas.h:209
void add(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:1275
Basic Geometric Tools.
gmm::uint16_type short_type
used as the common short type integer in the library
Definition: bgeot_config.h:72
size_t size_type
used as the common size type in the library
Definition: bgeot_poly.h:48
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.