27 struct contact_boundary {
32 mutable const model_real_plain_vector *U;
33 mutable model_real_plain_vector U_unred;
35 contact_boundary(
size_type r,
const mesh_fem *mf,
const std::string &dn)
36 : region(r), mfu(mf), dispname(dn)
41 base_small_vector element_U(
const contact_boundary &cb,
size_type cv)
43 auto U_elm = base_small_vector{};
49 auto most_central_box(
const bgeot::rtree::pbox_set &bset,
54 auto itmax = begin(bset);
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);
64 auto rate = min((*it)->max->at(i) - pt[i], pt[i] - (*it)->min->at(i)) / h;
65 rate_box = min(rate, rate_box);
68 if (rate_box > rate_max) {
80 class interpolate_transformation_on_deformed_domains
81 :
public virtual_interpolate_transformation {
83 contact_boundary master;
84 contact_boundary slave;
87 mutable std::map<size_type, std::vector<size_type>> box_to_convex;
90 mutable fem_precomp_pool fppool;
93 void compute_element_boxes()
const {
95 model_real_plain_vector Uelm;
96 element_boxes.clear();
98 auto bnum = master.region;
99 auto &mfu = *(master.mfu);
100 auto &U = *(master.U);
101 auto &m = mfu.linked_mesh();
104 base_node Xdeformed(N), bmin(N), bmax(N);
105 auto region = m.region(bnum);
110 region.prohibit_partitioning();
112 GMM_ASSERT1(mfu.get_qdim() == N,
"Wrong mesh_fem qdim");
114 dal::bit_vector points_already_interpolated;
115 std::vector<base_node> transformed_points(m.nb_max_points());
116 box_to_convex.clear();
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());
125 mfu.linked_mesh().points_of_convex(cv, G);
127 auto ctx = fem_interpolation_context{pgt, pfp,
size_type(-1), G, cv};
128 auto nb_pt = pgt->structure()->nb_points();
131 auto ind = m.ind_points_of_convex(cv)[k];
135 if (points_already_interpolated.is_in(ind)) {
136 Xdeformed = transformed_points[ind];
138 pf_s->interpolation(ctx, Uelm, Xdeformed, dim_type{N});
139 Xdeformed += ctx.xreal();
140 transformed_points[ind] = Xdeformed;
141 points_already_interpolated.add(ind);
145 bmin = bmax = Xdeformed;
148 bmin[l] = std::min(bmin[l], Xdeformed[l]);
149 bmax[l] = std::max(bmax[l], Xdeformed[l]);
155 box_to_convex[element_boxes.add_box(bmin, bmax)].push_back(cv);
157 element_boxes.build_tree();
160 fem_interpolation_context deformed_master_context(
size_type cv)
const
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);
171 std::vector<bgeot::base_node> deformed_master_nodes(
size_type cv)
const {
172 using namespace bgeot;
175 auto nodes = vector<base_node>{};
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);
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);
192 pfu->interpolation(ctx, U_elm, U, dim_type{N});
202 interpolate_transformation_on_deformed_domains(
205 const std::string &source_displacements,
208 const std::string &target_displacements)
210 slave{source_region, &mf_source, source_displacements},
211 master{target_region, &mf_target, target_displacements}
215 void extract_variables(
const ga_workspace &workspace,
216 std::set<var_trans_pair> &vars,
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,
"");
226 void init(
const ga_workspace &workspace)
const override {
228 for (
auto pcb : std::list<const contact_boundary*>{&master, &slave}) {
229 auto &mfu = *(pcb->mfu);
230 if (mfu.is_reduced()) {
232 mfu.extend_vector(workspace.value(pcb->dispname), pcb->U_unred);
233 pcb->U = &(pcb->U_unred);
235 pcb->U = &(workspace.value(pcb->dispname));
238 compute_element_boxes();
241 void finalize()
const override {
242 element_boxes.clear();
243 box_to_convex.clear();
244 master.U_unred.clear();
245 slave.U_unred.clear();
249 int transform(
const ga_workspace &workspace,
251 fem_interpolation_context &ctx_x,
258 std::map<var_trans_pair, base_tensor> &derivatives,
259 bool compute_derivatives)
const override {
261 auto &target_mesh = master.mfu->linked_mesh();
263 auto transformation_success =
false;
266 using namespace bgeot;
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();
276 auto G_x = base_matrix{};
277 m_x.points_of_convex(cv_x, G_x);
279 pfu_x->interpolation(ctx_x, U_elm_x, U_x, dim_type{N});
281 add(ctx_x.xreal(), U_x, pt_x);
289 auto bset = rtree::pbox_set{};
290 element_boxes.find_boxes_at_point(pt_x, bset);
291 while (!bset.empty())
293 auto itmax = most_central_box(bset, pt_x);
295 for (
auto i : box_to_convex.at((*itmax)->id))
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) {
304 transformation_success =
true;
308 if (transformation_success || (bset.size() == 1))
break;
317 if (compute_derivatives && transformation_success) {
318 GMM_ASSERT2(derivatives.size() == 2,
319 "Expecting to return derivatives only for Umaster and Uslave");
321 for (
auto &pair : derivatives)
323 if (pair.first.varname == slave.dispname)
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());
335 if (pair.first.varname == master.dispname)
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.);
351 else GMM_ASSERT2(
false,
"unexpected derivative variable");
355 return transformation_success ? 1 : 0;
361 (ga_workspace &workspace,
const std::string &transname,
362 const mesh &source_mesh,
const std::string &source_displacements,
364 const std::string &target_displacements,
const mesh_region &target_region)
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(),
371 source_displacements,
374 target_displacements);
375 workspace.add_interpolate_transformation(transname, p_transformation);
379 (
model &md,
const std::string &transname,
380 const mesh &source_mesh,
const std::string &source_displacements,
382 const std::string &target_displacements,
const mesh_region &target_region)
386 auto p_transformation
387 = std::make_shared<interpolate_transformation_on_deformed_domains>(source_region.id(),
389 source_displacements,
392 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.