GetFEM  5.5
getfem_generic_assembly_interpolation.cc
1 /*===========================================================================
2 
3  Copyright (C) 2013-2026 Yves Renard
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 
23 #include "getfem/getfem_generic_assembly_compile_and_exec.h"
24 #include "getfem/getfem_generic_assembly_functions_and_operators.h"
25 
26 namespace getfem {
27 
28  //=========================================================================
29  // Interpolation functions
30  //=========================================================================
31 
32  // general Interpolation
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);
38  }
39 
40  // Interpolation on a Lagrange fem on the same mesh
41  struct ga_interpolation_context_fem_same_mesh
42  : public ga_interpolation_context {
43  base_vector &result;
44  std::vector<int> dof_count;
45  const mesh_fem &mf;
46  bool initialized;
47  bool is_torus;
48  size_type s;
49 
50  virtual bgeot::pstored_point_tab
51  ppoints_for_element(size_type cv, short_type f,
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");
56 
57  if (f != short_type(-1))
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]);
62  else
63  for (size_type i = 0; i < pf->node_convex(cv).nb_points(); ++i)
64  ind.push_back(i);
65 
66  return pf->node_tab(cv);
67  }
68 
69  virtual bool use_pgp(size_type) const { return true; }
70  virtual bool use_mim() const { return false; }
71 
72  void init_(size_type si, size_type q, size_type qmult) {
73  s = si;
74  gmm::resize(result, qmult * mf.nb_basic_dof());
75  gmm::clear(result);
76  gmm::resize(dof_count, mf.nb_basic_dof()/q);
77  gmm::clear(dof_count);
78  initialized = true;
79  }
80 
81  virtual void store_result(size_type cv, size_type i, base_tensor &t) {
82  if (is_torus) {
83  size_type target_dim = mf.fem_of_element(cv)->target_dim();
84  GMM_ASSERT2(target_dim == 3, "Invalid torus fem.");
85  size_type qdim = 1;
86  size_type result_dim = 2;
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];
90  ++dof_count[idof];
91  } else {
92  size_type si = t.size();
93  size_type q = mf.get_qdim();
94  size_type qmult = si / q;
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];
100  gmm::add(t.as_vector(),
101  gmm::sub_vector(result, gmm::sub_interval(qmult*idof, s)));
102  (dof_count[idof/q])++;
103  }
104  }
105 
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);
112  if (!initialized) {
113  if (data[0]) {
114  gmm::resize(result, data[0]);
115  gmm::resize(dof_count, data[1]);
116  gmm::clear(dof_count);
117  s = data[2];
118  }
119  gmm::clear(result);
120  }
121  if (initialized)
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)
127  if (dof_count[i])
128  gmm::scale(gmm::sub_vector(result, gmm::sub_interval(s*i, s)),
129  scalar_type(1) / scalar_type(dof_count[i]));
130  }
131 
132  virtual const mesh &linked_mesh() { return mf.linked_mesh(); }
133 
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;
142  }
143  }
144  };
145 
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);
150  }
151 
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);
158  }
159 
160  // Interpolation on a cloud of points
161  struct ga_interpolation_context_mti
162  : public ga_interpolation_context {
163  base_vector &result;
164  const mesh_trans_inv &mti;
165  bool initialized;
166  size_type s, nbdof;
167 
168  virtual bgeot::pstored_point_tab
169  ppoints_for_element(size_type cv, short_type,
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]];
176  ind.push_back(i);
177  }
178  return store_point_tab(pt_tab);
179  }
180 
181  virtual bool use_pgp(size_type) const { return false; }
182  virtual bool use_mim() const { return false; }
183 
184  virtual void store_result(size_type cv, size_type i, base_tensor &t) {
185  size_type si = t.size();
186  if (!initialized) {
187  s = si;
188  gmm::resize(result, s * nbdof);
189  gmm::clear(result);
190  initialized = true;
191  }
192  GMM_ASSERT1(s == si, "Internal error");
193  size_type ipt = mti.point_on_convex(cv, i);
194  size_type dof_t = mti.id_of_point(ipt);
195  gmm::add(t.as_vector(),
196  gmm::sub_vector(result, gmm::sub_interval(s*dof_t, s)));
197  }
198 
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);
204  if (!initialized) {
205  if (data[0]) {
206  gmm::resize(result, data[0]);
207  s = data[1];
208  }
209  gmm::clear(result);
210  }
211  if (initialized)
212  GMM_ASSERT1(gmm::vect_size(result) == data[0] && s == data[1],
213  "Incompatible sizes");
214  MPI_SUM_VECTOR(result);
215  }
216 
217  virtual const mesh &linked_mesh() { return mti.linked_mesh(); }
218 
219  ga_interpolation_context_mti(const mesh_trans_inv &mti_, base_vector &r,
220  size_type nbdof_ = size_type(-1))
221  : result(r), mti(mti_), initialized(false), nbdof(nbdof_) {
222  if (nbdof == size_type(-1)) nbdof = mti.nb_points();
223  }
224  };
225 
226  // Distribute to be parallelized
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) {
231 
232  ga_workspace workspace(md);
233  workspace.add_interpolation_expression(expr, mti.linked_mesh(), rg);
234 
235  mti.distribute(extrapolation, rg_source);
236  ga_interpolation_context_mti gic(mti, result, nbdof);
237  ga_interpolation(workspace, gic);
238  }
239 
240 
241  // Interpolation on a im_data
242  struct ga_interpolation_context_im_data
243  : public ga_interpolation_context {
244  base_vector &result;
245  const im_data &imd;
246  bool initialized;
247  size_type s;
248 
249  virtual bgeot::pstored_point_tab
250  ppoints_for_element(size_type cv, short_type f,
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");
256  size_type i_start(0), i_end(0);
257  if (f == short_type(-1))
258  i_end = pim->approx_method()->nb_points_on_convex();
259  else {
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);
262  }
263  for (size_type i = i_start; i < i_end; ++i) ind.push_back(i);
264  return pim->approx_method()->pintegration_points();
265  }
266 
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());
273  }
274  virtual bool use_mim() const { return true; }
275 
276  virtual void store_result(size_type cv, size_type i, base_tensor &t) {
277  size_type si = t.size();
278  if (!initialized) {
279  s = si;
280  GMM_ASSERT1(imd.tensor_size() == t.sizes() ||
281  (imd.tensor_size().size() == size_type(1) &&
282  imd.tensor_size()[0] == size_type(1) &&
283  si == size_type(1)),
284  "Im_data tensor size " << imd.tensor_size() <<
285  " does not match the size of the interpolated "
286  "expression " << t.sizes() << ".");
287  gmm::resize(result, s * imd.nb_filtered_index());
288  gmm::clear(result);
289  initialized = true;
290  }
291  GMM_ASSERT1(s == si, "Internal error");
292  size_type ipt = imd.filtered_index_of_point(cv, i);
293  GMM_ASSERT1(ipt != size_type(-1),
294  "Im data with no data on the current integration point.");
295  gmm::add(t.as_vector(),
296  gmm::sub_vector(result, gmm::sub_interval(s*ipt, s)));
297  }
298 
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);
304  if (initialized) {
305  GMM_ASSERT1(gmm::vect_size(result) == data[0] && s == data[1],
306  "Incompatible sizes");
307  } else {
308  if (data[0]) {
309  gmm::resize(result, data[0]);
310  s = data[1];
311  }
312  gmm::clear(result);
313  }
314  MPI_SUM_VECTOR(result);
315  }
316 
317  virtual const mesh &linked_mesh() { return imd.linked_mesh(); }
318 
319  ga_interpolation_context_im_data(const im_data &imd_, base_vector &r)
320  : result(r), imd(imd_), initialized(false) { }
321  };
322 
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);
327  }
328 
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);
335 
336  ga_interpolation_im_data(workspace, imd, result);
337  }
338 
339 
340  // Interpolation on a stored_mesh_slice
341  struct ga_interpolation_context_mesh_slice
342  : public ga_interpolation_context {
343  base_vector &result;
344  const stored_mesh_slice &sl;
345  bool initialized;
346  size_type s;
347  std::vector<size_type> first_node;
348 
349  virtual bgeot::pstored_point_tab
350  ppoints_for_element(size_type cv, short_type f,
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.");
354  size_type ic = sl.convex_pos(cv);
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;
359  ind.push_back(i);
360  }
361  return store_point_tab(pt_tab);
362  }
363 
364  virtual bool use_pgp(size_type /* cv */) const { return false; } // why not?
365  virtual bool use_mim() const { return false; }
366 
367  virtual void store_result(size_type cv, size_type i, base_tensor &t) {
368  size_type si = t.size();
369  if (!initialized) {
370  s = si;
371  gmm::resize(result, s * sl.nb_points());
372  gmm::clear(result);
373  initialized = true;
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();
377  }
378  GMM_ASSERT1(s == si && result.size() == s * sl.nb_points(), "Internal error");
379  size_type ic = sl.convex_pos(cv);
380  size_type ipt = first_node[ic] + i;
381  gmm::add(t.as_vector(),
382  gmm::sub_vector(result, gmm::sub_interval(s*ipt, s)));
383  }
384 
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);
390  if (initialized) {
391  GMM_ASSERT1(gmm::vect_size(result) == data[0] && s == data[1],
392  "Incompatible sizes");
393  } else {
394  if (data[0]) {
395  gmm::resize(result, data[0]);
396  s = data[1];
397  }
398  gmm::clear(result);
399  }
400  MPI_SUM_VECTOR(result);
401  }
402 
403  virtual const mesh &linked_mesh() { return sl.linked_mesh(); }
404 
405  ga_interpolation_context_mesh_slice(const stored_mesh_slice &sl_,
406  base_vector &r)
407  : result(r), sl(sl_), initialized(false) { }
408  };
409 
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);
414  }
415 
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);
422  }
423 
424 
425  //=========================================================================
426  // Local projection functions
427  //=========================================================================
428 
429  void ga_local_projection(const getfem::model &md, const mesh_im &mim,
430  const std::string &expr, const mesh_fem &mf,
431  base_vector &result, const mesh_region &region) {
432 
433  // The case where the expression is a vector one and mf a scalar fem is
434  // not taken into account for the moment.
435 
436  // Could be improved by not performing the assembly of the global mass matrix
437  // working locally. This means a specific assembly.
438  size_type nbdof = mf.nb_dof();
439  model_real_sparse_matrix M(nbdof, nbdof);
440  asm_mass_matrix(M, mim, mf, region);
441  // FIXME: M should be cached for performance
442 
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);
448  else
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);
454 
455  getfem::base_matrix loc_M;
456  getfem::base_vector loc_U;
457  for (mr_visitor v(region, mf.linked_mesh(), true); !v.finished(); ++v) {
458  if (mf.convex_index().is_in(v.cv())) {
459  size_type nd = mf.nb_basic_dof_of_element(v.cv());
460  loc_M.base_resize(nd, nd); gmm::resize(loc_U, nd);
461  gmm::sub_index J(mf.ind_basic_dof_of_element(v.cv()));
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));
465  }
466  }
467  MPI_SUM_VECTOR(result);
468  }
469 
470  //=========================================================================
471  // Interpolate transformation with an expression
472  //=========================================================================
473 
474  struct rated_box_index_compare {
475  bool operator()(
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));
481  }
482  };
483 
484  class interpolate_transformation_expression
485  : public virtual_interpolate_transformation, public context_dependencies {
486 
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; }
490  };
491 
492  const mesh &source_mesh;
493  const mesh &target_mesh;
494  const size_type target_region;
495  std::string expr;
496  mutable bgeot::rtree element_boxes;
497  mutable bool recompute_elt_boxes;
498  mutable ga_workspace local_workspace;
499  mutable ga_instruction_set local_gis;
500  mutable bgeot::geotrans_inv_convex gic;
501  mutable base_node P;
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;
508 
509  private:
510  mutable std::map<size_type, std::vector<size_type>> box_to_convexes;
511 
512  public:
513  void update_from_context() const {
514  recompute_elt_boxes = true;
515  }
516 
517  void extract_variables(const ga_workspace &workspace,
518  std::set<var_trans_pair> &vars,
519  bool ignore_data, const mesh &/* m */,
520  const std::string &/* interpolate_name */) const {
521  if ((ignore_data && !extract_variable_done) ||
522  (!ignore_data && !extract_data_done)) {
523  if (ignore_data)
524  used_vars.clear();
525  else
526  used_data.clear();
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,
534  ignore_data);
535  if (ignore_data)
536  extract_variable_done = true;
537  else
538  extract_data_done = true;
539  }
540  if (ignore_data)
541  vars.insert(used_vars.begin(), used_vars.end());
542  else
543  vars.insert(used_data.begin(), used_data.end());
544  }
545 
546  void init(const ga_workspace &workspace) const {
547  size_type N = target_mesh.dim();
548 
549  // Expression compilation
550  local_workspace = ga_workspace(workspace, ga_workspace::inherit::ALL);
551  local_workspace.clear_expressions();
552 
553  local_workspace.add_interpolation_expression(expr, source_mesh);
554  local_gis = ga_instruction_set();
555  ga_compile_interpolation(local_workspace, local_gis);
556 
557  // In fact, transformations are not allowed ... for future compatibility
558  for (const std::string &transname : local_gis.transformations)
559  local_workspace.interpolate_transformation(transname)
560  ->init(local_workspace);
561 
562  if (!extract_variable_done) {
563  std::set<var_trans_pair> vars;
564  extract_variables(workspace, vars, true, source_mesh, "");
565  }
566 
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());
579  }
580  }
581 
582  // Element_boxes update (if necessary)
583  if (recompute_elt_boxes) {
584 
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&
590  convex_index = (target_region == mesh_region::all_convexes().id())
591  ? target_mesh.convex_index() : mr.index();
592  for (dal::bv_visitor cv(convex_index); !cv.finished(); ++cv) {
593 
594  bgeot::pgeometric_trans pgt = target_mesh.trans_of_convex(cv);
595 
596  size_type nbd_t = pgt->nb_points();
597  if (nbd_t) {
598  gmm::copy(target_mesh.points_of_convex(cv)[0], bmin);
599  gmm::copy(bmin, bmax);
600  } else {
601  gmm::clear(bmin);
602  gmm::clear(bmax);
603  }
604  mesh::ref_mesh_pt_ct cvpts = target_mesh.points_of_convex(cv);
605  for (short_type ip = 1; ip < nbd_t; ++ip) {
606  // size_type ind = target_mesh.ind_points_of_convex(cv)[ip];
607  const base_node &pt = cvpts[ip];
608  for (size_type k = 0; k < N; ++k) {
609  bmin[k] = std::min(bmin[k], pt[k]);
610  bmax[k] = std::max(bmax[k], pt[k]);
611  }
612  }
613 
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;
619 
620  box_to_convexes[element_boxes.add_box(bmin, bmax)].push_back(cv);
621  }
622  element_boxes.build_tree();
623  recompute_elt_boxes = false;
624  }
625  }
626 
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();
631  }
632 
633  std::string expression() const { return expr; }
634 
635  int transform(const ga_workspace &/*workspace*/, const mesh &m,
636  fem_interpolation_context &ctx_x,
637  const base_small_vector &Normal,
638  const mesh **m_t,
639  size_type &cv, short_type &face_num,
640  base_node &P_ref,
641  base_small_vector &/*N_y*/,
642  std::map<var_trans_pair, base_tensor> &derivatives,
643  bool compute_derivatives) const {
644  int ret_type = 0;
645 
646  local_gis.ctx = ctx_x;
647  local_gis.Normal = Normal;
648  local_gis.nbpt = 1;
649  local_gis.ipt = 0;
650  local_gis.pai = 0;
651  gmm::clear(local_workspace.assembled_tensor().as_vector());
652 
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();
662  }
663 
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);
668 
669  *m_t = &target_mesh;
670 
671  bgeot::rtree::pbox_cont boxes;
672  {
673  bgeot::rtree::pbox_set bset;
674  element_boxes.find_boxes_at_point(P, bset);
675 
676  // using a std::set as a sorter
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);
687  }
688  }
689  rated_boxes.insert(std::make_pair(rating, box));
690  }
691 
692  // boxes should now be ordered in increasing rating order
693  for (const auto &p : rated_boxes)
694  boxes.push_back(p.second);
695  }
696 
697  scalar_type best_dist(1e10);
698  size_type best_cv(-1);
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));
704 
705  bool converged;
706  bool is_in = gic.invert(P, P_ref, converged, 1E-4);
707  // cout << "cv = " << convex << " P = " << P << " P_ref = " << P_ref;
708  // cout << " is_in = " << int(is_in) << endl;
709  // for (size_type iii = 0;
710  // iii < target_mesh.points_of_convex(cv).size(); ++iii)
711  // cout << target_mesh.points_of_convex(cv)[iii] << endl;
712 
713  if (converged) {
714  cv = convex;
715  if (is_in) {
716  face_num = short_type(-1); // Should detect potential faces ?
717  ret_type = 1;
718  break;
719  } else {
720  scalar_type dist
721  = target_mesh.trans_of_convex(cv)->convex_ref()->is_in(P_ref);
722  if (dist < best_dist) {
723  best_dist = dist;
724  best_cv = cv;
725  best_P_ref = P_ref;
726  }
727  }
728  }
729  }
730  }
731 
732  if (ret_type == 0 && best_dist < 5e-3) {
733  cv = best_cv;
734  P_ref = best_P_ref;
735  face_num = short_type(-1); // Should detect potential faces ?
736  ret_type = 1;
737  }
738 
739  // Note on derivatives of the transformation : for efficiency and
740  // simplicity reasons, the derivative should be computed with
741  // the value of corresponding test functions. This means that
742  // for a transformation F(u) the computed derivative is F'(u).Test_u
743  // including the Test_u.
744  if (compute_derivatives) { // Could be optimized ?
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;
749  pwi.gis().nbpt = 1;
750  pwi.gis().ipt = 0;
751  pwi.gis().pai = 0;
752  gmm::clear(pwi.workspace().assembled_tensor().as_vector());
753  ga_function_exec(pwi.gis());
754  d.second = pwi.workspace().assembled_tensor();
755  }
756  }
757 
758  return ret_type;
759  }
760 
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); }
767 
768  };
769 
770 
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);
776  }
777 
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>
783  (sm, tm, trg, expr);
784  workspace.add_interpolate_transformation(name, p);
785  }
786 
788  (model &md, const std::string &name, const mesh &sm, const mesh &tm,
789  const std::string &expr) {
791  (md, name, sm, tm, size_type(-1), expr);
792  }
793 
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>
799  (sm, tm, trg, expr);
800  md.add_interpolate_transformation(name, p);
801  }
802 
803  //=========================================================================
804  // Interpolate transformation on neighbor element (for internal faces)
805  //=========================================================================
806 
807  class interpolate_transformation_neighbor
808  : public virtual_interpolate_transformation, public context_dependencies {
809 
810  public:
811  void update_from_context() const {}
812  void extract_variables(const ga_workspace &/* workspace */,
813  std::set<var_trans_pair> &/* vars */,
814  bool /* ignore_data */, const mesh &/* m */,
815  const std::string &/* interpolate_name */) const {}
816  void init(const ga_workspace &/* workspace */) const {}
817  void finalize() const {}
818 
819  std::string expression() const { return "X"; }
820 
821  int transform(const ga_workspace &/*workspace*/, const mesh &m_x,
822  fem_interpolation_context &ctx_x,
823  const base_small_vector &/*Normal*/, const mesh **m_t,
824  size_type &cv, short_type &face_num,
825  base_node &P_ref,
826  base_small_vector &/*N_y*/,
827  std::map<var_trans_pair, base_tensor> &/*derivatives*/,
828  bool compute_derivatives) const {
829 
830  int ret_type = 0;
831  *m_t = &m_x;
832  size_type cv_x = ctx_x.convex_num();
833  short_type face_x = ctx_x.face_num();
834  GMM_ASSERT1(face_x != short_type(-1), "Neighbor transformation can "
835  "only be applied to internal faces");
836 
837  auto adj_face = m_x.adjacent_face(cv_x, face_x);
838 
839  if (adj_face.cv != size_type(-1)) {
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;
849  cv = adj_face.cv;
850  ret_type = 1;
851  }
852  GMM_ASSERT1(!compute_derivatives,
853  "No derivative for this transformation");
854  return ret_type;
855  }
856 
857  interpolate_transformation_neighbor() { }
858 
859  };
860 
861 
862  pinterpolate_transformation interpolate_transformation_neighbor_instance() {
863  return (std::make_shared<interpolate_transformation_neighbor>());
864  }
865 
866  //=========================================================================
867  // Interpolate transformation on neighbor element (for extrapolation)
868  //=========================================================================
869 
870  class interpolate_transformation_element_extrapolation
871  : public virtual_interpolate_transformation, public context_dependencies {
872 
873  const mesh &sm;
874  std::map<size_type, size_type> elt_corr;
875 
876  public:
877  void update_from_context() const {}
878  void extract_variables(const ga_workspace &/* workspace */,
879  std::set<var_trans_pair> &/* vars */,
880  bool /* ignore_data */, const mesh &/* m */,
881  const std::string &/* interpolate_name */) const {}
882  void init(const ga_workspace &/* workspace */) const {}
883  void finalize() const {}
884  std::string expression() const { return "X"; }
885 
886  int transform(const ga_workspace &/*workspace*/, const mesh &m_x,
887  fem_interpolation_context &ctx_x,
888  const base_small_vector &/*Normal*/, const mesh **m_t,
889  size_type &cv, short_type &face_num,
890  base_node &P_ref,
891  base_small_vector &/*N_y*/,
892  std::map<var_trans_pair, base_tensor> &/*derivatives*/,
893  bool compute_derivatives) const {
894  int ret_type = 0;
895  *m_t = &m_x;
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;
900 
901  if (cv_x != cv_y) {
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");
909  face_num = short_type(-1);
910  cv = cv_y;
911  ret_type = 1;
912  } else {
913  cv = cv_x;
914  face_num = short_type(-1);
915  P_ref = ctx_x.xref();
916  ret_type = 1;
917  }
918  GMM_ASSERT1(!compute_derivatives,
919  "No derivative for this transformation");
920  return ret_type;
921  }
922 
923  void set_correspondence(const std::map<size_type, size_type> &ec) {
924  elt_corr = ec;
925  }
926 
927  interpolate_transformation_element_extrapolation
928  (const mesh &sm_, const std::map<size_type, size_type> &ec)
929  : sm(sm_), elt_corr(ec) { }
930  };
931 
932 
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>
938  (sm, elt_corr);
939  md.add_interpolate_transformation(name, p);
940  }
941 
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>
947  (sm, elt_corr);
948  workspace.add_interpolate_transformation(name, p);
949  }
950 
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();
957  auto cpext
958  = dynamic_cast<const interpolate_transformation_element_extrapolation *>
959  (pit);
960  GMM_ASSERT1(cpext,
961  "The transformation is not of element extrapolation type");
962  const_cast<interpolate_transformation_element_extrapolation *>(cpext)
963  ->set_correspondence(elt_corr);
964  }
965 
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();
972  auto cpext
973  = dynamic_cast<const interpolate_transformation_element_extrapolation *>
974  (pit);
975  GMM_ASSERT1(cpext,
976  "The transformation is not of element extrapolation type");
977  const_cast<interpolate_transformation_element_extrapolation *>(cpext)
978  ->set_correspondence(elt_corr);
979  }
980 
981 
982  //=========================================================================
983  // Secondary domains
984  //=========================================================================
985 
986 
987  class standard_secondary_domain : public virtual_secondary_domain {
988 
989  public:
990 
991  virtual const mesh_region &give_region(const mesh &,
992  size_type, short_type) const
993  { return region; }
994  // virtual void init(const ga_workspace &workspace) const = 0;
995  // virtual void finalize() const = 0;
996 
997  standard_secondary_domain(const mesh_im &mim__, const mesh_region &region_)
998  : virtual_secondary_domain(mim__, region_) {}
999  };
1000 
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);
1006  }
1007 
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);
1013  }
1014 
1015 
1016 } /* end of namespace */
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.
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).
Definition: getfem_mesh.h:98
`‘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)
*‍/
Definition: gmm_blas.h:976
void clear(L &l)
clear (fill with zeros) a vector or matrix.
Definition: gmm_blas.h:58
void resize(V &v, size_type n)
*‍/
Definition: gmm_blas.h:209
void add(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:1275
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
Definition: getfem_fem.h:243
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
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.