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