GetFEM  5.4.2
getfem_generic_assembly_interpolation.cc
1 /*===========================================================================
2 
3  Copyright (C) 2013-2020 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, write to the Free Software Foundation,
18  Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301, USA.
19 
20 ===========================================================================*/
21 
24 #include "getfem/getfem_generic_assembly_compile_and_exec.h"
25 #include "getfem/getfem_generic_assembly_functions_and_operators.h"
26 
27 namespace getfem {
28 
29  //=========================================================================
30  // Interpolation functions
31  //=========================================================================
32 
33  // general Interpolation
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);
39  }
40 
41  // Interpolation on a Lagrange fem on the same mesh
42  struct ga_interpolation_context_fem_same_mesh
43  : public ga_interpolation_context {
44  base_vector &result;
45  std::vector<int> dof_count;
46  const mesh_fem &mf;
47  bool initialized;
48  bool is_torus;
49  size_type s;
50 
51  virtual bgeot::pstored_point_tab
52  ppoints_for_element(size_type cv, short_type f,
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");
57 
58  if (f != short_type(-1)) {
59 
60  for (size_type i = 0;
61  i < pf->node_convex(cv).structure()->nb_points_of_face(f); ++i)
62  ind.push_back
63  (pf->node_convex(cv).structure()->ind_points_of_face(f)[i]);
64  } else {
65  for (size_type i = 0; i < pf->node_convex(cv).nb_points(); ++i)
66  ind.push_back(i);
67  }
68 
69  return pf->node_tab(cv);
70  }
71 
72  virtual bool use_pgp(size_type) const { return true; }
73  virtual bool use_mim() const { return false; }
74 
75  void init_(size_type si, size_type q, size_type qmult) {
76  s = si;
77  gmm::resize(result, qmult * mf.nb_basic_dof());
78  gmm::clear(result);
79  gmm::resize(dof_count, mf.nb_basic_dof()/q);
80  gmm::clear(dof_count);
81  initialized = true;
82  }
83 
84  void store_result_for_torus(size_type cv, size_type i, base_tensor &t) {
85  size_type target_dim = mf.fem_of_element(cv)->target_dim();
86  GMM_ASSERT2(target_dim == 3, "Invalid torus fem.");
87  size_type qdim = 1;
88  size_type result_dim = 2;
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];
92  ++dof_count[idof];
93  }
94 
95  virtual void store_result(size_type cv, size_type i, base_tensor &t) {
96  if (is_torus){
97  store_result_for_torus(cv, i, t);
98  return;
99  }
100  size_type si = t.size();
101  size_type q = mf.get_qdim();
102  size_type qmult = si / q;
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])++;
111  }
112 
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);
119  if (!initialized) {
120  if (data[0]) {
121  gmm::resize(result, data[0]);
122  gmm::resize(dof_count, data[1]);
123  gmm::clear(dof_count);
124  s = data[2];
125  }
126  gmm::clear(result);
127  }
128  if (initialized)
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)
134  if (dof_count[i])
135  gmm::scale(gmm::sub_vector(result, gmm::sub_interval(s*i, s)),
136  scalar_type(1) / scalar_type(dof_count[i]));
137  }
138 
139  virtual const mesh &linked_mesh() { return mf.linked_mesh(); }
140 
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;
149  }
150  }
151  };
152 
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);
157  }
158 
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);
165  }
166 
167  // Interpolation on a cloud of points
168  struct ga_interpolation_context_mti
169  : public ga_interpolation_context {
170  base_vector &result;
171  const mesh_trans_inv &mti;
172  bool initialized;
173  size_type s, nbdof;
174 
175 
176  virtual bgeot::pstored_point_tab
177  ppoints_for_element(size_type cv, short_type,
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]];
184  ind.push_back(i);
185  }
186  return store_point_tab(pt_tab);
187  }
188 
189  virtual bool use_pgp(size_type) const { return false; }
190  virtual bool use_mim() const { return false; }
191 
192  virtual void store_result(size_type cv, size_type i, base_tensor &t) {
193  size_type si = t.size();
194  if (!initialized) {
195  s = si;
196  gmm::resize(result, s * nbdof);
197  gmm::clear(result);
198  initialized = true;
199  }
200  GMM_ASSERT1(s == si, "Internal error");
201  size_type ipt = mti.point_on_convex(cv, i);
202  size_type dof_t = mti.id_of_point(ipt);
203  gmm::add(t.as_vector(),
204  gmm::sub_vector(result, gmm::sub_interval(s*dof_t, s)));
205  }
206 
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);
212  if (!initialized) {
213  if (data[0]) {
214  gmm::resize(result, data[0]);
215  s = data[1];
216  }
217  gmm::clear(result);
218  }
219  if (initialized)
220  GMM_ASSERT1(gmm::vect_size(result) == data[0] && s == data[1],
221  "Incompatible sizes");
222  MPI_SUM_VECTOR(result);
223  }
224 
225  virtual const mesh &linked_mesh() { return mti.linked_mesh(); }
226 
227  ga_interpolation_context_mti(const mesh_trans_inv &mti_, base_vector &r,
228  size_type nbdof_ = size_type(-1))
229  : result(r), mti(mti_), initialized(false), nbdof(nbdof_) {
230  if (nbdof == size_type(-1)) nbdof = mti.nb_points();
231  }
232  };
233 
234  // Distribute to be parallelized
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) {
239 
240  ga_workspace workspace(md);
241  workspace.add_interpolation_expression(expr, mti.linked_mesh(), rg);
242 
243  mti.distribute(extrapolation, rg_source);
244  ga_interpolation_context_mti gic(mti, result, nbdof);
245  ga_interpolation(workspace, gic);
246  }
247 
248 
249  // Interpolation on a im_data
250  struct ga_interpolation_context_im_data
251  : public ga_interpolation_context {
252  base_vector &result;
253  const im_data &imd;
254  bool initialized;
255  size_type s;
256 
257  virtual bgeot::pstored_point_tab
258  ppoints_for_element(size_type cv, short_type f,
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");
264  size_type i_start(0), i_end(0);
265  if (f == short_type(-1))
266  i_end = pim->approx_method()->nb_points_on_convex();
267  else {
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);
270  }
271  for (size_type i = i_start; i < i_end; ++i) ind.push_back(i);
272  return pim->approx_method()->pintegration_points();
273  }
274 
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());
281  }
282  virtual bool use_mim() const { return true; }
283 
284  virtual void store_result(size_type cv, size_type i, base_tensor &t) {
285  size_type si = t.size();
286  if (!initialized) {
287  s = si;
288  GMM_ASSERT1(imd.tensor_size() == t.sizes() ||
289  (imd.tensor_size().size() == size_type(1) &&
290  imd.tensor_size()[0] == size_type(1) &&
291  si == size_type(1)),
292  "Im_data tensor size " << imd.tensor_size() <<
293  " does not match the size of the interpolated "
294  "expression " << t.sizes() << ".");
295  gmm::resize(result, s * imd.nb_filtered_index());
296  gmm::clear(result);
297  initialized = true;
298  }
299  GMM_ASSERT1(s == si, "Internal error");
300  size_type ipt = imd.filtered_index_of_point(cv, i);
301  GMM_ASSERT1(ipt != size_type(-1),
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)));
305  }
306 
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);
312  if (initialized) {
313  GMM_ASSERT1(gmm::vect_size(result) == data[0] && s == data[1],
314  "Incompatible sizes");
315  } else {
316  if (data[0]) {
317  gmm::resize(result, data[0]);
318  s = data[1];
319  }
320  gmm::clear(result);
321  }
322  MPI_SUM_VECTOR(result);
323  }
324 
325  virtual const mesh &linked_mesh() { return imd.linked_mesh(); }
326 
327  ga_interpolation_context_im_data(const im_data &imd_, base_vector &r)
328  : result(r), imd(imd_), initialized(false) { }
329  };
330 
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);
335  }
336 
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);
343 
344  ga_interpolation_im_data(workspace, imd, result);
345  }
346 
347 
348  // Interpolation on a stored_mesh_slice
349  struct ga_interpolation_context_mesh_slice
350  : public ga_interpolation_context {
351  base_vector &result;
352  const stored_mesh_slice &sl;
353  bool initialized;
354  size_type s;
355  std::vector<size_type> first_node;
356 
357  virtual bgeot::pstored_point_tab
358  ppoints_for_element(size_type cv, short_type f,
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.");
362  size_type ic = sl.convex_pos(cv);
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;
367  ind.push_back(i);
368  }
369  return store_point_tab(pt_tab);
370  }
371 
372  virtual bool use_pgp(size_type /* cv */) const { return false; } // why not?
373  virtual bool use_mim() const { return false; }
374 
375  virtual void store_result(size_type cv, size_type i, base_tensor &t) {
376  size_type si = t.size();
377  if (!initialized) {
378  s = si;
379  gmm::resize(result, s * sl.nb_points());
380  gmm::clear(result);
381  initialized = true;
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();
385  }
386  GMM_ASSERT1(s == si && result.size() == s * sl.nb_points(), "Internal error");
387  size_type ic = sl.convex_pos(cv);
388  size_type ipt = first_node[ic] + i;
389  gmm::add(t.as_vector(),
390  gmm::sub_vector(result, gmm::sub_interval(s*ipt, s)));
391  }
392 
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);
398  if (initialized) {
399  GMM_ASSERT1(gmm::vect_size(result) == data[0] && s == data[1],
400  "Incompatible sizes");
401  } else {
402  if (data[0]) {
403  gmm::resize(result, data[0]);
404  s = data[1];
405  }
406  gmm::clear(result);
407  }
408  MPI_SUM_VECTOR(result);
409  }
410 
411  virtual const mesh &linked_mesh() { return sl.linked_mesh(); }
412 
413  ga_interpolation_context_mesh_slice(const stored_mesh_slice &sl_,
414  base_vector &r)
415  : result(r), sl(sl_), initialized(false) { }
416  };
417 
418  void ga_interpolation_mesh_slice
419  (ga_workspace &workspace, const stored_mesh_slice &sl, base_vector &result) {
420  ga_interpolation_context_mesh_slice gic(sl, result);
421  ga_interpolation(workspace, gic);
422  }
423 
424  void ga_interpolation_mesh_slice
425  (const getfem::model &md, const std::string &expr, const stored_mesh_slice &sl,
426  base_vector &result, const mesh_region &rg) {
427  ga_workspace workspace(md);
428  workspace.add_interpolation_expression(expr, sl.linked_mesh(), rg);
429  ga_interpolation_mesh_slice(workspace, sl, result);
430  }
431 
432 
433  //=========================================================================
434  // Local projection functions
435  //=========================================================================
436 
437  void ga_local_projection(const getfem::model &md, const mesh_im &mim,
438  const std::string &expr, const mesh_fem &mf,
439  base_vector &result, const mesh_region &region) {
440 
441  // The case where the expression is a vector one and mf a scalar fem is
442  // not taken into account for the moment.
443 
444  // Could be improved by not performing the assembly of the global mass matrix
445  // working locally. This means a specific assembly.
446  size_type nbdof = mf.nb_dof();
447  model_real_sparse_matrix M(nbdof, nbdof);
448  asm_mass_matrix(M, mim, mf, region);
449  // FIXME: M should be cached for performance
450 
451  ga_workspace workspace(md, ga_workspace::inherit::NONE);
452  gmm::sub_interval I(0,nbdof);
453  workspace.add_fem_variable("c__dummy_var_95_", mf, I, base_vector(nbdof));
454  if (mf.get_qdims().size() > 1)
455  workspace.add_expression("("+expr+"):Test_c__dummy_var_95_",mim,region,2);
456  else
457  workspace.add_expression("("+expr+").Test_c__dummy_var_95_",mim,region,2);
458  getfem::base_vector F(nbdof);
459  workspace.set_assembled_vector(F);
460  workspace.assembly(1);
461  gmm::resize(result, nbdof);
462 
463  getfem::base_matrix loc_M;
464  getfem::base_vector loc_U;
465  for (mr_visitor v(region, mf.linked_mesh(), true); !v.finished(); ++v) {
466  if (mf.convex_index().is_in(v.cv())) {
467  size_type nd = mf.nb_basic_dof_of_element(v.cv());
468  loc_M.base_resize(nd, nd); gmm::resize(loc_U, nd);
469  gmm::sub_index J(mf.ind_basic_dof_of_element(v.cv()));
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));
473  }
474  }
475  MPI_SUM_VECTOR(result);
476  }
477 
478  //=========================================================================
479  // Interpolate transformation with an expression
480  //=========================================================================
481 
482  struct rated_box_index_compare {
483  bool operator()(
484  const std::pair<scalar_type, const bgeot::box_index*> &x,
485  const std::pair<scalar_type, const bgeot::box_index*> &y) const {
486  GMM_ASSERT2(x.second != nullptr, "x contains nullptr");
487  GMM_ASSERT2(y.second != nullptr, "y contains nullptr");
488  return (x.first < y.first) || (!(y.first < x.first) && (x.second->id < y.second->id));
489  }
490  };
491 
492  class interpolate_transformation_expression
493  : public virtual_interpolate_transformation, public context_dependencies {
494 
495  struct workspace_gis_pair : public std::pair<ga_workspace, ga_instruction_set> {
496  inline ga_workspace &workspace() { return this->first; }
497  inline ga_instruction_set &gis() { return this->second; }
498  };
499 
500  const mesh &source_mesh;
501  const mesh &target_mesh;
502  const size_type target_region;
503  std::string expr;
504  mutable bgeot::rtree element_boxes;
505  mutable bool recompute_elt_boxes;
506  mutable ga_workspace local_workspace;
507  mutable ga_instruction_set local_gis;
508  mutable bgeot::geotrans_inv_convex gic;
509  mutable base_node P;
510  mutable std::set<var_trans_pair> used_vars;
511  mutable std::set<var_trans_pair> used_data;
512  mutable std::map<var_trans_pair,
513  workspace_gis_pair> compiled_derivatives;
514  mutable bool extract_variable_done;
515  mutable bool extract_data_done;
516 
517  private:
518  mutable std::map<size_type, std::vector<size_type>> box_to_convexes;
519 
520  public:
521  void update_from_context() const {
522  recompute_elt_boxes = true;
523  }
524 
525  void extract_variables(const ga_workspace &workspace,
526  std::set<var_trans_pair> &vars,
527  bool ignore_data, const mesh &/* m */,
528  const std::string &/* interpolate_name */) const {
529  if ((ignore_data && !extract_variable_done) ||
530  (!ignore_data && !extract_data_done)) {
531  if (ignore_data)
532  used_vars.clear();
533  else
534  used_data.clear();
535  ga_workspace aux_workspace(workspace, ga_workspace::inherit::ALL);
536  aux_workspace.clear_expressions();
537  aux_workspace.add_interpolation_expression(expr, source_mesh);
538  for (size_type i = 0; i < aux_workspace.nb_trees(); ++i)
539  ga_extract_variables(aux_workspace.tree_info(i).ptree->root,
540  aux_workspace, source_mesh,
541  ignore_data ? used_vars : used_data,
542  ignore_data);
543  if (ignore_data)
544  extract_variable_done = true;
545  else
546  extract_data_done = true;
547  }
548  if (ignore_data)
549  vars.insert(used_vars.begin(), used_vars.end());
550  else
551  vars.insert(used_data.begin(), used_data.end());
552  }
553 
554  void init(const ga_workspace &workspace) const {
555  size_type N = target_mesh.dim();
556 
557  // Expression compilation
558  local_workspace = ga_workspace(workspace, ga_workspace::inherit::ALL);
559  local_workspace.clear_expressions();
560 
561  local_workspace.add_interpolation_expression(expr, source_mesh);
562  local_gis = ga_instruction_set();
563  ga_compile_interpolation(local_workspace, local_gis);
564 
565  // In fact, transformations are not allowed ... for future compatibility
566  for (const std::string &transname : local_gis.transformations)
567  local_workspace.interpolate_transformation(transname)
568  ->init(local_workspace);
569 
570  if (!extract_variable_done) {
571  std::set<var_trans_pair> vars;
572  extract_variables(workspace, vars, true, source_mesh, "");
573  }
574 
575  for (const var_trans_pair &var : used_vars) {
576  workspace_gis_pair &pwi = compiled_derivatives[var];
577  pwi.workspace() = local_workspace;
578  pwi.gis() = ga_instruction_set();
579  if (pwi.workspace().nb_trees()) {
580  ga_tree tree = *(pwi.workspace().tree_info(0).ptree);
581  ga_derivative(tree, pwi.workspace(), source_mesh,
582  var.varname, var.transname, 1);
583  if (tree.root)
584  ga_semantic_analysis(tree, local_workspace, dummy_mesh(),
585  1, false, true);
586  ga_compile_interpolation(pwi.workspace(), pwi.gis());
587  }
588  }
589 
590  // Element_boxes update (if necessary)
591  if (recompute_elt_boxes) {
592 
593  box_to_convexes.clear();
594  element_boxes.clear();
595  base_node bmin(N), bmax(N);
596  const dal::bit_vector&
597  convex_index = (target_region == mesh_region::all_convexes().id())
598  ? target_mesh.convex_index()
599  : target_mesh.region(target_region).index();
600  for (dal::bv_visitor cv(convex_index); !cv.finished(); ++cv) {
601 
602  bgeot::pgeometric_trans pgt = target_mesh.trans_of_convex(cv);
603 
604  size_type nbd_t = pgt->nb_points();
605  if (nbd_t) {
606  gmm::copy(target_mesh.points_of_convex(cv)[0], bmin);
607  gmm::copy(bmin, bmax);
608  } else {
609  gmm::clear(bmin);
610  gmm::clear(bmax);
611  }
612  for (short_type ip = 1; ip < nbd_t; ++ip) {
613  // size_type ind = target_mesh.ind_points_of_convex(cv)[ip];
614  const base_node &pt = target_mesh.points_of_convex(cv)[ip];
615 
616  for (size_type k = 0; k < N; ++k) {
617  bmin[k] = std::min(bmin[k], pt[k]);
618  bmax[k] = std::max(bmax[k], pt[k]);
619  }
620  }
621 
622  scalar_type h = bmax[0] - bmin[0];
623  for (size_type k = 1; k < N; ++k) h = std::max(h, bmax[k]-bmin[k]);
624  if (pgt->is_linear()) h *= 1E-4;
625  for (auto &&val : bmin) val -= h*0.2;
626  for (auto &&val : bmax) val += h*0.2;
627 
628  box_to_convexes[element_boxes.add_box(bmin, bmax)].push_back(cv);
629  }
630  element_boxes.build_tree();
631  recompute_elt_boxes = false;
632  }
633  }
634 
635  void finalize() const {
636  for (const std::string &transname : local_gis.transformations)
637  local_workspace.interpolate_transformation(transname)->finalize();
638  local_gis = ga_instruction_set();
639  }
640 
641  std::string expression() const { return expr; }
642 
643  int transform(const ga_workspace &/*workspace*/, const mesh &m,
644  fem_interpolation_context &ctx_x,
645  const base_small_vector &Normal,
646  const mesh **m_t,
647  size_type &cv, short_type &face_num,
648  base_node &P_ref,
649  base_small_vector &/*N_y*/,
650  std::map<var_trans_pair, base_tensor> &derivatives,
651  bool compute_derivatives) const {
652  int ret_type = 0;
653 
654  local_gis.ctx = ctx_x;
655  local_gis.Normal = Normal;
656  local_gis.nbpt = 1;
657  local_gis.ipt = 0;
658  local_gis.pai = 0;
659  gmm::clear(local_workspace.assembled_tensor().as_vector());
660 
661  for (auto &&instr : local_gis.all_instructions) {
662  GMM_ASSERT1(instr.second.m == &m,
663  "Incompatibility of meshes in interpolation");
664  auto &gilb = instr.second.begin_instructions;
665  for (size_type j = 0; j < gilb.size(); ++j) j += gilb[j]->exec();
666  auto &gile = instr.second.elt_instructions;
667  for (size_type j = 0; j < gile.size(); ++j) j+=gile[j]->exec();
668  auto &gil = instr.second.instructions;
669  for (size_type j = 0; j < gil.size(); ++j) j += gil[j]->exec();
670  }
671 
672  GMM_ASSERT1(local_workspace.assembled_tensor().size()==target_mesh.dim(),
673  "Wrong dimension of the transformation expression");
674  P.resize(target_mesh.dim());
675  gmm::copy(local_workspace.assembled_tensor().as_vector(), P);
676 
677  *m_t = &target_mesh;
678 
679  bgeot::rtree::pbox_cont boxes;
680  {
681  bgeot::rtree::pbox_set bset;
682  element_boxes.find_boxes_at_point(P, bset);
683 
684  // using a std::set as a sorter
685  std::set<std::pair<scalar_type, const bgeot::box_index*>,
686  rated_box_index_compare> rated_boxes;
687  for (const auto &box : bset) {
688  scalar_type rating = scalar_type(1);
689  for (size_type i = 0; i < m.dim(); ++i) {
690  scalar_type h = box->max->at(i) - box->min->at(i);
691  if (h > scalar_type(0)) {
692  scalar_type r = std::min(box->max->at(i) - P[i],
693  P[i] - box->min->at(i)) / h;
694  rating = std::min(r, rating);
695  }
696  }
697  rated_boxes.insert(std::make_pair(rating, box));
698  }
699 
700  // boxes should now be ordered in increasing rating order
701  for (const auto &p : rated_boxes)
702  boxes.push_back(p.second);
703  }
704 
705 
706  scalar_type best_dist(1e10);
707  size_type best_cv(-1);
708  base_node best_P_ref;
709  for (size_type i = boxes.size(); i > 0 && !ret_type; --i) {
710  for (auto convex : box_to_convexes.at(boxes[i-1]->id)) {
711  gic.init(target_mesh.points_of_convex(convex),
712  target_mesh.trans_of_convex(convex));
713 
714  bool converged;
715  bool is_in = gic.invert(P, P_ref, converged, 1E-4);
716  // cout << "cv = " << convex << " P = " << P << " P_ref = " << P_ref;
717  // cout << " is_in = " << int(is_in) << endl;
718  // for (size_type iii = 0;
719  // iii < target_mesh.points_of_convex(cv).size(); ++iii)
720  // cout << target_mesh.points_of_convex(cv)[iii] << endl;
721 
722  if (converged) {
723  cv = convex;
724  if (is_in) {
725  face_num = short_type(-1); // Should detect potential faces ?
726  ret_type = 1;
727  break;
728  } else {
729  scalar_type dist
730  = target_mesh.trans_of_convex(cv)->convex_ref()->is_in(P_ref);
731  if (dist < best_dist) {
732  best_dist = dist;
733  best_cv = cv;
734  best_P_ref = P_ref;
735  }
736  }
737  }
738  }
739  }
740 
741  if (ret_type == 0 && best_dist < 5e-3) {
742  cv = best_cv;
743  P_ref = best_P_ref;
744  face_num = short_type(-1); // Should detect potential faces ?
745  ret_type = 1;
746  }
747 
748  // Note on derivatives of the transformation : for efficiency and
749  // simplicity reasons, the derivative should be computed with
750  // the value of corresponding test functions. This means that
751  // for a transformation F(u) the computed derivative is F'(u).Test_u
752  // including the Test_u.
753  if (compute_derivatives) { // Could be optimized ?
754  for (auto &&d : derivatives) {
755  workspace_gis_pair &pwi = compiled_derivatives[d.first];
756 
757  gmm::clear(pwi.workspace().assembled_tensor().as_vector());
758  ga_function_exec(pwi.gis());
759  d.second = pwi.workspace().assembled_tensor();
760  }
761  }
762  return ret_type;
763  }
764 
765  interpolate_transformation_expression
766  (const mesh &sm, const mesh &tm, size_type trg, const std::string &expr_)
767  : source_mesh(sm), target_mesh(tm), target_region(trg), expr(expr_),
768  recompute_elt_boxes(true), extract_variable_done(false),
769  extract_data_done(false)
770  { this->add_dependency(tm); }
771 
772  };
773 
774 
776  (ga_workspace &workspace, const std::string &name, const mesh &sm,
777  const mesh &tm, const std::string &expr) {
779  (workspace, name, sm, tm, size_type(-1), expr);
780  }
781 
783  (ga_workspace &workspace, const std::string &name, const mesh &sm,
784  const mesh &tm, size_type trg, const std::string &expr) {
785  pinterpolate_transformation
786  p = std::make_shared<interpolate_transformation_expression>
787  (sm, tm, trg, expr);
788  workspace.add_interpolate_transformation(name, p);
789  }
790 
792  (model &md, const std::string &name, const mesh &sm, const mesh &tm,
793  const std::string &expr) {
795  (md, name, sm, tm, size_type(-1), expr);
796  }
797 
799  (model &md, const std::string &name, const mesh &sm, const mesh &tm,
800  size_type trg, const std::string &expr) {
801  pinterpolate_transformation
802  p = std::make_shared<interpolate_transformation_expression>
803  (sm, tm, trg, expr);
804  md.add_interpolate_transformation(name, p);
805  }
806 
807  //=========================================================================
808  // Interpolate transformation on neighbor element (for internal faces)
809  //=========================================================================
810 
811  class interpolate_transformation_neighbor
812  : public virtual_interpolate_transformation, public context_dependencies {
813 
814  public:
815  void update_from_context() const {}
816  void extract_variables(const ga_workspace &/* workspace */,
817  std::set<var_trans_pair> &/* vars */,
818  bool /* ignore_data */, const mesh &/* m */,
819  const std::string &/* interpolate_name */) const {}
820  void init(const ga_workspace &/* workspace */) const {}
821  void finalize() const {}
822 
823  std::string expression(void) const { return "X"; }
824 
825  int transform(const ga_workspace &/*workspace*/, const mesh &m_x,
826  fem_interpolation_context &ctx_x,
827  const base_small_vector &/*Normal*/, const mesh **m_t,
828  size_type &cv, short_type &face_num,
829  base_node &P_ref,
830  base_small_vector &/*N_y*/,
831  std::map<var_trans_pair, base_tensor> &/*derivatives*/,
832  bool compute_derivatives) const {
833 
834  int ret_type = 0;
835  *m_t = &m_x;
836  size_type cv_x = ctx_x.convex_num();
837  short_type face_x = ctx_x.face_num();
838  GMM_ASSERT1(face_x != short_type(-1), "Neighbor transformation can "
839  "only be applied to internal faces");
840 
841  auto adj_face = m_x.adjacent_face(cv_x, face_x);
842 
843  if (adj_face.cv != size_type(-1)) {
845  gic.init(m_x.points_of_convex(adj_face.cv),
846  m_x.trans_of_convex(adj_face.cv));
847  bool converged = true;
848  gic.invert(ctx_x.xreal(), P_ref, converged);
849  bool is_in = (ctx_x.pgt()->convex_ref()->is_in(P_ref) < 1E-4);
850  GMM_ASSERT1(is_in && converged, "Geometric transformation inversion "
851  "has failed in neighbor transformation");
852  face_num = adj_face.f;
853  cv = adj_face.cv;
854  ret_type = 1;
855  }
856  GMM_ASSERT1(!compute_derivatives,
857  "No derivative for this transformation");
858  return ret_type;
859  }
860 
861  interpolate_transformation_neighbor() { }
862 
863  };
864 
865 
866  pinterpolate_transformation interpolate_transformation_neighbor_instance() {
867  return (std::make_shared<interpolate_transformation_neighbor>());
868  }
869 
870  //=========================================================================
871  // Interpolate transformation on neighbor element (for extrapolation)
872  //=========================================================================
873 
874  class interpolate_transformation_element_extrapolation
875  : public virtual_interpolate_transformation, public context_dependencies {
876 
877  const mesh &sm;
878  std::map<size_type, size_type> elt_corr;
879 
880  public:
881  void update_from_context() const {}
882  void extract_variables(const ga_workspace &/* workspace */,
883  std::set<var_trans_pair> &/* vars */,
884  bool /* ignore_data */, const mesh &/* m */,
885  const std::string &/* interpolate_name */) const {}
886  void init(const ga_workspace &/* workspace */) const {}
887  void finalize() const {}
888  std::string expression(void) const { return "X"; }
889 
890  int transform(const ga_workspace &/*workspace*/, const mesh &m_x,
891  fem_interpolation_context &ctx_x,
892  const base_small_vector &/*Normal*/, const mesh **m_t,
893  size_type &cv, short_type &face_num,
894  base_node &P_ref,
895  base_small_vector &/*N_y*/,
896  std::map<var_trans_pair, base_tensor> &/*derivatives*/,
897  bool compute_derivatives) const {
898  int ret_type = 0;
899  *m_t = &m_x;
900  GMM_ASSERT1(&sm == &m_x, "Bad mesh");
901  size_type cv_x = ctx_x.convex_num(), cv_y = cv_x;
902  auto it = elt_corr.find(cv_x);
903  if (it != elt_corr.end()) cv_y = it->second;
904 
905  if (cv_x != cv_y) {
907  gic.init(m_x.points_of_convex(cv_y),
908  m_x.trans_of_convex(cv_y));
909  bool converged = true;
910  gic.invert(ctx_x.xreal(), P_ref, converged, 1E-4);
911  GMM_ASSERT1(converged, "Geometric transformation inversion "
912  "has failed in element extrapolation transformation");
913  face_num = short_type(-1);
914  cv = cv_y;
915  ret_type = 1;
916  } else {
917  cv = cv_x;
918  face_num = short_type(-1);
919  P_ref = ctx_x.xref();
920  ret_type = 1;
921  }
922  GMM_ASSERT1(!compute_derivatives,
923  "No derivative for this transformation");
924  return ret_type;
925  }
926 
927  void set_correspondence(const std::map<size_type, size_type> &ec) {
928  elt_corr = ec;
929  }
930 
931  interpolate_transformation_element_extrapolation
932  (const mesh &sm_, const std::map<size_type, size_type> &ec)
933  : sm(sm_), elt_corr(ec) { }
934  };
935 
936 
937  void add_element_extrapolation_transformation
938  (model &md, const std::string &name, const mesh &sm,
939  std::map<size_type, size_type> &elt_corr) {
940  pinterpolate_transformation
941  p = std::make_shared<interpolate_transformation_element_extrapolation>
942  (sm, elt_corr);
943  md.add_interpolate_transformation(name, p);
944  }
945 
946  void add_element_extrapolation_transformation
947  (ga_workspace &workspace, const std::string &name, const mesh &sm,
948  std::map<size_type, size_type> &elt_corr) {
949  pinterpolate_transformation
950  p = std::make_shared<interpolate_transformation_element_extrapolation>
951  (sm, elt_corr);
952  workspace.add_interpolate_transformation(name, p);
953  }
954 
955  void set_element_extrapolation_correspondence
956  (ga_workspace &workspace, const std::string &name,
957  std::map<size_type, size_type> &elt_corr) {
958  GMM_ASSERT1(workspace.interpolate_transformation_exists(name),
959  "Unknown transformation");
960  auto pit=workspace.interpolate_transformation(name).get();
961  auto cpext
962  = dynamic_cast<const interpolate_transformation_element_extrapolation *>
963  (pit);
964  GMM_ASSERT1(cpext,
965  "The transformation is not of element extrapolation type");
966  const_cast<interpolate_transformation_element_extrapolation *>(cpext)
967  ->set_correspondence(elt_corr);
968  }
969 
970  void set_element_extrapolation_correspondence
971  (model &md, const std::string &name,
972  std::map<size_type, size_type> &elt_corr) {
973  GMM_ASSERT1(md.interpolate_transformation_exists(name),
974  "Unknown transformation");
975  auto pit=md.interpolate_transformation(name).get();
976  auto cpext
977  = dynamic_cast<const interpolate_transformation_element_extrapolation *>
978  (pit);
979  GMM_ASSERT1(cpext,
980  "The transformation is not of element extrapolation type");
981  const_cast<interpolate_transformation_element_extrapolation *>(cpext)
982  ->set_correspondence(elt_corr);
983  }
984 
985 
986  //=========================================================================
987  // Secondary domains
988  //=========================================================================
989 
990 
991  class standard_secondary_domain : public virtual_secondary_domain {
992 
993  public:
994 
995  virtual const mesh_region &give_region(const mesh &,
996  size_type, short_type) const
997  { return region; }
998  // virtual void init(const ga_workspace &workspace) const = 0;
999  // virtual void finalize() const = 0;
1000 
1001  standard_secondary_domain(const mesh_im &mim__, const mesh_region &region_)
1002  : virtual_secondary_domain(mim__, region_) {}
1003  };
1004 
1005  void add_standard_secondary_domain
1006  (model &md, const std::string &name, const mesh_im &mim,
1007  const mesh_region &rg) {
1008  psecondary_domain p = std::make_shared<standard_secondary_domain>(mim, rg);
1009  md.add_secondary_domain(name, p);
1010  }
1011 
1012  void add_standard_secondary_domain
1013  (ga_workspace &workspace, const std::string &name, const mesh_im &mim,
1014  const mesh_region &rg) {
1015  psecondary_domain p = std::make_shared<standard_secondary_domain>(mim, rg);
1016  workspace.add_secondary_domain(name, p);
1017  }
1018 
1019 
1020 } /* end of namespace */
bgeot::geotrans_inv_convex
does the inversion of the geometric transformation for a given convex
Definition: bgeot_geotrans_inv.h:64
gmm::resize
void resize(M &v, size_type m, size_type n)
*‍/
Definition: gmm_blas.h:231
bgeot::size_type
size_t size_type
used as the common size type in the library
Definition: bgeot_poly.h:49
gmm::clear
void clear(L &l)
clear (fill with zeros) a vector or matrix.
Definition: gmm_blas.h:59
getfem::mesh_im
Describe an integration method linked to a mesh.
Definition: getfem_mesh_im.h:47
getfem_generic_assembly_semantic.h
Semantic analysis of assembly trees and semantic manipulations.
getfem_generic_assembly_tree.h
Compilation and execution operations.
getfem::asm_mass_matrix
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)
Definition: getfem_assembling.h:697
getfem::interpolate_transformation_neighbor_instance
pinterpolate_transformation interpolate_transformation_neighbor_instance()
Create a new instance of a transformation corresponding to the interpolation on the neighbor element.
Definition: getfem_generic_assembly_interpolation.cc:866
getfem::mesh_fem
Describe a finite element method linked to a mesh.
Definition: getfem_mesh_fem.h:148
getfem::mesh_region::all_convexes
static mesh_region all_convexes()
provide a default value for the mesh_region parameters of assembly procedures etc.
Definition: getfem_mesh_region.h:142
getfem::model
`‘Model’' variables store the variables, the data and the description of a model.
Definition: getfem_models.h:114
bgeot::short_type
gmm::uint16_type short_type
used as the common short type integer in the library
Definition: bgeot_config.h:72
bgeot::rtree
Balanced tree of n-dimensional rectangles.
Definition: bgeot_rtree.h:98
getfem
GEneric Tool for Finite Element Methods.
Definition: getfem_accumulated_distro.h:46
getfem::mesh_fem::convex_index
const dal::bit_vector & convex_index() const
Get the set of convexes where a finite element has been assigned.
Definition: getfem_mesh_fem.h:191
getfem::mesh_fem::nb_basic_dof_of_element
virtual size_type nb_basic_dof_of_element(size_type cv) const
Return the number of degrees of freedom attached to a given convex.
Definition: getfem_mesh_fem.h:494
getfem::mesh_region::visitor
"iterator" class for regions.
Definition: getfem_mesh_region.h:237
getfem::ga_local_projection
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.
Definition: getfem_generic_assembly_interpolation.cc:437
getfem::pfem
std::shared_ptr< const getfem::virtual_fem > pfem
type of pointer on a fem description
Definition: getfem_fem.h:244
getfem::mesh_fem::linked_mesh
const mesh & linked_mesh() const
Return a reference to the underlying mesh.
Definition: getfem_mesh_fem.h:279
getfem::mesh
Describe a mesh (collection of convexes (elements) and points).
Definition: getfem_mesh.h:95
bgeot::pgeometric_trans
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation
Definition: bgeot_geometric_trans.h:186
getfem::mesh_fem::ind_basic_dof_of_element
virtual ind_dof_ct ind_basic_dof_of_element(size_type cv) const
Give an array of the dof numbers a of convex.
Definition: getfem_mesh_fem.h:450
getfem::mesh_region
structure used to hold a set of convexes and/or convex faces.
Definition: getfem_mesh_region.h:55
getfem::add_interpolate_transformation_from_expression
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...
Definition: getfem_generic_assembly_interpolation.cc:776
bgeot::geotrans_inv_convex::invert
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 ...
Definition: bgeot_geotrans_inv.cc:61
getfem::mesh_fem::nb_dof
virtual size_type nb_dof() const
Return the total number of degrees of freedom.
Definition: getfem_mesh_fem.h:562

Rabisu Mirror Service We provide mirrors to support Open source communities. Our mirror server is located in Istanbul/Turkey region.

Please do not hesitate to contact mirror@rabisu.com for new open source mirror submissions.