GetFEM  5.5
getfem_interpolated_fem.cc
1 /*===========================================================================
2 
3  Copyright (C) 2004-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 
22 
23 namespace getfem {
24 
25  void interpolated_fem::build_rtree(void) const {
26  base_node min, max;
27  boxtree.clear();
28  box_to_convexes_map.clear();
29  for (dal::bv_visitor cv(mf.convex_index()); !cv.finished(); ++cv) {
30  bounding_box(min, max, mf.linked_mesh().points_of_convex(cv),
31  mf.linked_mesh().trans_of_convex(cv));
32  box_to_convexes_map[boxtree.add_box(min, max, cv)].push_back(cv);
33  }
34  boxtree.build_tree();
35  }
36 
37  bool interpolated_fem::find_a_point(base_node pt, base_node &ptr,
38  size_type &cv) const {
39  bool gt_invertible;
40  if (pif) { base_node ptreal = pt; pif->val(ptreal, pt); }
41  if (cv_stored != size_type(-1) && gic.invert(pt, ptr, gt_invertible))
42  { cv = cv_stored; if (gt_invertible) return true; }
43  boxtree.find_boxes_at_point(pt, boxlst);
44  bgeot::rtree::pbox_set::const_iterator it = boxlst.begin(),
45  ite = boxlst.end();
46  for (; it != ite; ++it) {
47  for (auto candidate : box_to_convexes_map.at((*it)->id)) {
49  (mf.linked_mesh().convex(candidate),
50  mf.linked_mesh().trans_of_convex(candidate));
51  cv_stored = candidate;
52  if (gic.invert(pt, ptr, gt_invertible)) {
53  cv = candidate; return true;
54  }
55  }
56  }
57  return false;
58  }
59 
61  fictx_cv = cv_stored = size_type(-1);
62  dim_ = dim_type(-1);
63  build_rtree();
64 
65  GMM_ASSERT1(!mf.is_reduced(),
66  "Interpolated fem works only on non reduced mesh_fems");
67 
68  std::vector<elt_interpolation_data> vv(mim.convex_index().last_true() + 1);
69  elements.swap(vv);
70  base_node gpt;
71  ind_dof.resize(mf.nb_basic_dof());
72  dal::bit_vector alldofs;
73  size_type max_dof = 0;
74  if (mim.convex_index().card() == 0) return;
75  for (dal::bv_visitor cv(mim.convex_index()); !cv.finished(); ++cv) {
76  if (dim_ == dim_type(-1))
77  dim_ = mim.linked_mesh().structure_of_convex(cv)->dim();
78 
79  GMM_ASSERT1(dim_ == mim.linked_mesh().structure_of_convex(cv)->dim(),
80  "Convexes of different dimension: to be done");
81  pintegration_method pim = mim.int_method_of_element(cv);
82  GMM_ASSERT1(pim->type() == IM_APPROX, "You have to use approximated "
83  "integration to interpolate a fem");
84  papprox_integration pai = pim->approx_method();
85  bgeot::pgeometric_trans pgt = mim.linked_mesh().trans_of_convex(cv);
86  elements[cv].gausspt.resize(pai->nb_points());
87  dal::bit_vector dofs;
88  size_type last_cv = size_type(-1);
89  for (size_type k = 0; k < pai->nb_points(); ++k) {
90  gausspt_interpolation_data &gpid = elements[cv].gausspt[k];
91  /* todo: use a geotrans_interpolation_context */
92  gpt = pgt->transform(pai->point(k),
93  mim.linked_mesh().points_of_convex(cv));
94  gpid.iflags = find_a_point(gpt, gpid.ptref, gpid.elt) ? 1 : 0;
95  if (gpid.iflags && last_cv != gpid.elt) {
96  size_type nbd = mf.nb_basic_dof_of_element(gpid.elt);
97  for (size_type i = 0; i < nbd; ++i) {
98  size_type idof = mf.ind_basic_dof_of_element(gpid.elt)[i];
99  if (!(blocked_dof[idof])) dofs.add(idof);
100  }
101  last_cv = gpid.elt;
102  }
103  }
104  elements[cv].nb_dof = dofs.card();
105  elements[cv].pim = pim;
106  max_dof = std::max(max_dof, dofs.card());
107  elements[cv].inddof.resize(dofs.card());
108  size_type cnt = 0;
109  for (dal::bv_visitor idof(dofs); !idof.finished(); ++idof)
110  { elements[cv].inddof[cnt] = idof; ind_dof[idof] = cnt++; }
111  for (size_type k = 0; k < pai->nb_points(); ++k) {
112  gausspt_interpolation_data &gpid = elements[cv].gausspt[k];
113  if (gpid.iflags) {
114  size_type nbd = mf.nb_basic_dof_of_element(gpid.elt);
115  gpid.local_dof.resize(nbd);
116  for (size_type i = 0; i < nbd; ++i) {
117  size_type ndof = mf.ind_basic_dof_of_element(gpid.elt)[i];
118  gpid.local_dof[i] = dofs.is_in(ndof) ? ind_dof[ndof]
119  : size_type(-1);
120  }
121  }
122  }
123  alldofs |= dofs;
124  }
125  /** setup global dofs, with dummy coordinates */
126  base_node P(dim()); gmm::fill(P,1./20);
127  std::vector<base_node> node_tab_(max_dof, P);
128  pspt_override = bgeot::store_point_tab(node_tab_);
129  pspt_valid = false;
130  dof_types_.resize(max_dof);
131  std::fill(dof_types_.begin(), dof_types_.end(),
132  global_dof(dim()));
133 
134  /* ind_dof should be kept full of -1 ( real_base_value and
135  grad_base_value expect that)
136  */
137  std::fill(ind_dof.begin(), ind_dof.end(), size_type(-1));
138  }
139 
141  { context_check();
142  if (mim.linked_mesh().convex_index().is_in(cv))
143  return elements[cv].nb_dof;
144  else GMM_ASSERT1(false, "Wrong convex number: " << cv);
145  }
146 
147  size_type interpolated_fem::index_of_global_dof
148  (size_type cv, size_type i) const
149  { return elements[cv].inddof[i]; }
150 
151  bgeot::pconvex_ref interpolated_fem::ref_convex(size_type cv) const
152  { return mim.int_method_of_element(cv)->approx_method()->ref_convex(); }
153 
155  (size_type cv) const
156  {
157  if (mim.linked_mesh().convex_index().is_in(cv))
159  (dim(), nb_dof(cv),
160  mim.linked_mesh().structure_of_convex(cv)->nb_faces()));
161  else GMM_ASSERT1(false, "Wrong convex number: " << cv);
162  }
163 
164  void interpolated_fem::base_value(const base_node &, base_tensor &) const
165  { GMM_ASSERT1(false, "No base values, real only element."); }
166  void interpolated_fem::grad_base_value(const base_node &,
167  base_tensor &) const
168  { GMM_ASSERT1(false, "No grad values, real only element."); }
169  void interpolated_fem::hess_base_value(const base_node &,
170  base_tensor &) const
171  { GMM_ASSERT1(false, "No hess values, real only element."); }
172 
173  inline void interpolated_fem::actualize_fictx(pfem pf, size_type cv,
174  const base_node &ptr) const {
175  if (fictx_cv != cv) {
176  mf.linked_mesh().points_of_convex(cv, G);
178  (mf.linked_mesh().trans_of_convex(cv), pf, base_node(), G, cv,
179  short_type(-1));
180  fictx_cv = cv;
181  }
182  fictx.set_xref(ptr);
183  }
184 
186  base_tensor &t, bool) const {
187  elt_interpolation_data &e = elements.at(c.convex_num());
188  size_type cv;
189 
190  mi2[1] = target_dim(); mi2[0] = short_type(e.nb_dof);
191  t.adjust_sizes(mi2);
192  std::fill(t.begin(), t.end(), scalar_type(0));
193  if (e.nb_dof == 0) return;
194 
195  if (c.have_pgp() &&
196  (c.pgp()->get_ppoint_tab()
197  == e.pim->approx_method()->pintegration_points())) {
198  gausspt_interpolation_data &gpid = e.gausspt.at(c.ii());
199  if (gpid.iflags & 1) {
200  cv = gpid.elt;
201  pfem pf = mf.fem_of_element(cv);
202  unsigned rdim = target_dim() / pf->target_dim(), mdim = (rdim==1) ? 0 : 1;
203  if (gpid.iflags & 2) { t = gpid.base_val; return; }
204  actualize_fictx(pf, cv, gpid.ptref);
205  pf->real_base_value(fictx, taux);
206  for (size_type i = 0; i < pf->nb_dof(cv); ++i)
207  if (gpid.local_dof[i*rdim] != size_type(-1))
208  for (size_type j = 0; j < target_dim(); ++j)
209  t(gpid.local_dof[i*rdim+j*mdim],j) = taux(i, j*(1-mdim));
210  if (store_values) { gpid.base_val = t; gpid.iflags |= 2; }
211  }
212  }
213  else {
214  if (find_a_point(c.xreal(), ptref, cv)) {
215  pfem pf = mf.fem_of_element(cv);
216  unsigned rdim = target_dim() / pf->target_dim(), mdim = (rdim==1) ? 0 : 1;
217  actualize_fictx(pf, cv, ptref);
218  pf->real_base_value(fictx, taux);
219  for (size_type i = 0; i < e.nb_dof; ++i) {
220  ind_dof.at(e.inddof[i]) = i;
221  }
222  for (size_type i = 0; i < pf->nb_dof(cv); ++i)
223  for (size_type j = 0; j < target_dim(); ++j)
224  if (ind_dof.at(mf.ind_basic_dof_of_element(cv)[i*rdim+j*mdim])
225  != size_type(-1)) {
226  t(ind_dof[mf.ind_basic_dof_of_element(cv)[i*rdim+j*mdim]], j)
227  = taux(i, j*(1-mdim));
228  }
229  for (size_type i = 0; i < elements[c.convex_num()].nb_dof; ++i)
230  ind_dof[e.inddof[i]] = size_type(-1);
231  }
232  }
233 
234  }
235 
237  (const fem_interpolation_context& c, base_tensor &t, bool) const {
238  size_type N0 = mf.linked_mesh().dim();
239  elt_interpolation_data &e = elements.at(c.convex_num());
240  size_type nbdof = e.nb_dof, cv;
241 
242  mi3[2] = short_type(N0); mi3[1] = target_dim(); mi3[0] = short_type(nbdof);
243  t.adjust_sizes(mi3);
244  std::fill(t.begin(), t.end(), scalar_type(0));
245  if (nbdof == 0) return;
246 
247  if (c.have_pgp() &&
248  (c.pgp()->get_ppoint_tab()
249  == e.pim->approx_method()->pintegration_points())) {
250  gausspt_interpolation_data &gpid = e.gausspt.at(c.ii());
251  if (gpid.iflags & 1) {
252  cv = gpid.elt;
253  pfem pf = mf.fem_of_element(cv);
254  unsigned rdim = target_dim() / pf->target_dim(), mdim = (rdim==1) ? 0 : 1;
255  if (gpid.iflags & 4) { t = gpid.grad_val; return; }
256  actualize_fictx(pf, cv, gpid.ptref);
257  pf->real_grad_base_value(fictx, taux);
258 
259  if (pif) {
260  pif->grad(c.xreal(), trans);
261  for (size_type i = 0; i < pf->nb_dof(cv); ++i)
262  if (gpid.local_dof[i*rdim] != size_type(-1))
263  for (size_type j = 0; j < target_dim(); ++j)
264  for (size_type k = 0; k < N0; ++k) {
265  scalar_type ee(0);
266  for (size_type l = 0; l < N0; ++l)
267  ee += trans(l, k) * taux(i, j*(1-mdim), l);
268  t(gpid.local_dof[i*rdim+j*mdim], j, k) = ee;
269  }
270  }
271  else {
272  for (size_type i = 0; i < pf->nb_dof(cv); ++i)
273  if (gpid.local_dof[i*rdim] != size_type(-1))
274  for (size_type j = 0; j < target_dim(); ++j)
275  for (size_type k = 0; k < N0; ++k)
276  t(gpid.local_dof[i*rdim+j*mdim], j, k)
277  = taux(i, j*(1-mdim), k);
278  if (store_values) { gpid.grad_val = t; gpid.iflags |= 4; }
279  }
280  }
281  }
282  else {
283  cerr << "NON PGP OU MAUVAIS PTS sz=" << elements.size() << ", cv="
284  << c.convex_num() << " ";
285  cerr << "ii=" << c.ii() << ", sz=" << e.gausspt.size() << "\n";
286 
287  if (find_a_point(c.xreal(), ptref, cv)) {
288  pfem pf = mf.fem_of_element(cv);
289  unsigned rdim = target_dim() / pf->target_dim(), mdim = (rdim==1) ? 0 : 1;
290  actualize_fictx(pf, cv, ptref);
291  pf->real_grad_base_value(fictx, taux);
292  for (size_type i = 0; i < nbdof; ++i)
293  ind_dof.at(e.inddof.at(i)) = i;
294  if (pif) {
295  pif->grad(c.xreal(), trans);
296  for (size_type i = 0; i < pf->nb_dof(cv); ++i)
297  for (size_type j = 0; j < target_dim(); ++j)
298  for (size_type k = 0; k < N0; ++k)
299  if (ind_dof[mf.ind_basic_dof_of_element(cv)[i*rdim+j*mdim]]
300  != size_type(-1)) {
301  scalar_type ee(0);
302  for (size_type l = 0; l < N0; ++l)
303  ee += trans(l, k) * taux(i, j*(1-mdim), l);
304  t(ind_dof[mf.ind_basic_dof_of_element(cv)[i*rdim+j*mdim]],j,k)=ee;
305  }
306  }
307  else {
308  for (size_type i = 0; i < pf->nb_dof(cv); ++i)
309  for (size_type j = 0; j < target_dim(); ++j)
310  for (size_type k = 0; k < N0; ++k)
311  if (ind_dof[mf.ind_basic_dof_of_element(cv)[i*rdim+j*mdim]]
312  != size_type(-1))
313  t(ind_dof[mf.ind_basic_dof_of_element(cv)[i*rdim+j*mdim]],j,k)
314  = taux(i,j*(1-mdim),k);
315  }
316  for (size_type i = 0; i < nbdof; ++i)
317  ind_dof[e.inddof[i]] = size_type(-1);
318  }
319  }
320  }
321 
323  (const fem_interpolation_context&, base_tensor &, bool) const
324  { GMM_ASSERT1(false, "Sorry, to be done."); }
325 
326 
327  dal::bit_vector interpolated_fem::interpolated_convexes() const {
328  dal::bit_vector bv;
329  for (dal::bv_visitor cv(mim.linked_mesh().convex_index()); !cv.finished();
330  ++cv) {
331  for (unsigned ii=0; ii < elements.at(cv).gausspt.size(); ++ii) {
332  if (elements[cv].gausspt[ii].iflags)
333  bv.add(elements[cv].gausspt[ii].elt);
334  }
335  }
336  return bv;
337  }
338 
339  void interpolated_fem::gauss_pts_stats(unsigned &ming, unsigned &maxg,
340  scalar_type &meang) const {
341  std::vector<unsigned> v(mf.linked_mesh().nb_allocated_convex());
342  for (dal::bv_visitor cv(mim.linked_mesh().convex_index());
343  !cv.finished(); ++cv) {
344  for (unsigned ii=0; ii < elements.at(cv).gausspt.size(); ++ii) {
345  if (elements[cv].gausspt[ii].iflags)
346  v[elements[cv].gausspt[ii].elt]++;
347  }
348  }
349  ming = 100000; maxg = 0; meang = 0;
350  for (dal::bv_visitor cv(mf.linked_mesh().convex_index());
351  !cv.finished(); ++cv) {
352  ming = std::min(ming, v[cv]);
353  maxg = std::max(maxg, v[cv]);
354  meang += v[cv];
355  }
356  meang /= scalar_type(mf.linked_mesh().convex_index().card());
357  }
358 
359  size_type interpolated_fem::memsize() const {
360  size_type sz = 0;
361  sz += blocked_dof.memsize();
362  sz += sizeof(*this);
363  sz += elements.capacity() * sizeof(elt_interpolation_data);
364  for (unsigned i=0; i < elements.size(); ++i) {
365  sz += elements[i].gausspt.capacity()*sizeof(gausspt_interpolation_data);
366  sz += elements[i].inddof.capacity() * sizeof(size_type);
367  for (unsigned j=0; j < elements[i].gausspt.size(); ++j) {
368  sz += elements[i].gausspt[j].local_dof.capacity() * sizeof(size_type);
369  }
370  }
371  return sz;
372  }
373 
374  interpolated_fem::interpolated_fem(const mesh_fem &mef,
375  const mesh_im &meim,
376  pinterpolated_func pif_,
377  dal::bit_vector blocked_dof_,
378  bool store_val)
379  : mf(mef), mim(meim), pif(pif_), store_values(store_val),
380  blocked_dof(blocked_dof_), boxtree{1E-13}, mi2(2), mi3(3) {
381  DAL_STORED_OBJECT_DEBUG_CREATED(this, "Interpolated fem");
382  this->add_dependency(mf);
383  this->add_dependency(mim);
384  is_pol = is_lag = is_standard_fem = false; es_degree = 5;
385  is_equiv = real_element_defined = true;
386  gmm::resize(trans, mf.linked_mesh().dim(), mf.linked_mesh().dim());
387  ntarget_dim = mef.get_qdim();
389  }
390 
391  DAL_SIMPLE_KEY(special_intfem_key, pfem);
392 
393  pfem new_interpolated_fem(const mesh_fem &mef, const mesh_im &mim,
394  pinterpolated_func pif,
395  dal::bit_vector blocked_dof, bool store_val) {
396  pfem pf = std::make_shared<interpolated_fem>
397  (mef, mim, pif, blocked_dof, store_val);
398  dal::pstatic_stored_object_key pk=std::make_shared<special_intfem_key>(pf);
399  dal::add_stored_object(pk, pf);
400  return pf;
401  }
402 
403 
404 } /* end of namespace getfem. */
405 
const base_node & xreal() const
coordinates of the current point, in the real convex.
void set_xref(const base_node &P)
change the current point (coordinates given in the reference convex)
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 ...
const dal::bit_vector & convex_index() const
Return the list of valid convex IDs.
pconvex_structure structure_of_convex(size_type ic) const
Return the pconvex_structure of the convex ic.
size_type nb_allocated_convex() const
The number of convex indexes from 0 to the index of the last convex.
bool context_check() const
return true if update_from_context was called
structure passed as the argument of fem interpolation functions.
Definition: getfem_fem.h:749
void real_grad_base_value(const fem_interpolation_context &c, base_tensor &t, bool=true) const
Give the gradient of all components of the base functions at the current point of the fem_interpolati...
void real_base_value(const fem_interpolation_context &c, base_tensor &t, bool=true) const
Give the value of all components of the base functions at the current point of the fem_interpolation_...
void grad_base_value(const base_node &, base_tensor &) const
Give the value of all gradients (on ref.
virtual size_type nb_dof(size_type cv) const
Number of degrees of freedom.
virtual void update_from_context(void) const
this function has to be defined and should update the object when the context is modified.
virtual const bgeot::convex< base_node > & node_convex(size_type cv) const
Gives the convex representing the nodes on the reference element.
void hess_base_value(const base_node &, base_tensor &) const
Give the value of all hessians (on ref.
void base_value(const base_node &, base_tensor &) const
Give the value of all components of the base functions at the point x of the reference element.
void gauss_pts_stats(unsigned &ming, unsigned &maxg, scalar_type &meang) const
return the min/max/mean number of gauss points in the convexes of the interpolated mesh_fem
void real_hess_base_value(const fem_interpolation_context &, base_tensor &, bool=true) const
Give the hessian of all components of the base functions at the current point of the fem_interpolatio...
dal::bit_vector interpolated_convexes() const
return the list of convexes of the interpolated mesh_fem which contain at least one gauss point (shou...
virtual bgeot::pconvex_ref ref_convex(size_type cv) const
Return the convex of the reference element.
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.
const mesh & linked_mesh() const
Return a reference to the underlying mesh.
virtual size_type nb_basic_dof() const
Return the total number of basic degrees of freedom (before the optional reduction).
virtual pfem fem_of_element(size_type cv) const
Return the basic fem associated with an element (if no fem is associated, the function will crash!...
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.
bool is_reduced() const
Return true if a reduction matrix is applied to the dofs.
Describe an integration method linked to a mesh.
virtual pintegration_method int_method_of_element(size_type cv) const
return the integration method associated with an element (in no integration is associated,...
const mesh & linked_mesh() const
Give a reference to the linked mesh of type mesh.
const dal::bit_vector & convex_index(void) const
Get the set of convexes where an integration method has been assigned.
ref_convex convex(size_type ic) const
return a bgeot::convex object for the convex number ic.
Definition: getfem_mesh.h:188
FEM which interpolates a mesh_fem on a different mesh.
void resize(V &v, size_type n)
*‍/
Definition: gmm_blas.h:209
pdof_description global_dof(dim_type d)
Description of a global dof, i.e.
Definition: getfem_fem.cc:542
dim_type target_dim() const
dimension of the target space.
Definition: getfem_fem.h:313
size_type convex_num() const
get the current convex number
Definition: getfem_fem.cc:52
std::shared_ptr< const getfem::virtual_fem > pfem
type of pointer on a fem description
Definition: getfem_fem.h:243
dim_type dim() const
dimension of the reference element.
Definition: getfem_fem.h:310
gmm::uint16_type short_type
used as the common short type integer in the library
Definition: bgeot_config.h:72
pconvex_ref generic_dummy_convex_ref(dim_type nc, size_type n, short_type nf)
generic convex with n global nodes
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
void add_stored_object(pstatic_stored_object_key k, pstatic_stored_object o, permanence perm)
Add an object with two optional dependencies.
GEneric Tool for Finite Element Methods.
pfem new_interpolated_fem(const mesh_fem &mef, const mesh_im &mim, pinterpolated_func pif=0, dal::bit_vector blocked_dof=dal::bit_vector(), bool store_val=true)
create a new interpolated FEM.