26 void interpolated_fem::build_rtree(
void)
const {
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),
33 box_to_convexes_map[boxtree.add_box(min, max, cv)].push_back(cv);
38 bool interpolated_fem::find_a_point(base_node pt, base_node &ptr,
41 if (pif) { base_node ptreal = pt; pif->val(ptreal, pt); }
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(),
47 for (; it != ite; ++it) {
48 for (
auto candidate : box_to_convexes_map.at((*it)->id)) {
52 cv_stored = candidate;
53 if (gic.
invert(pt, ptr, gt_invertible)) {
54 cv = candidate;
return true;
67 "Interpolated fem works only on non reduced mesh_fems");
69 std::vector<elt_interpolation_data> vv(mim.
convex_index().last_true() + 1);
73 dal::bit_vector alldofs;
76 for (dal::bv_visitor cv(mim.
convex_index()); !cv.finished(); ++cv) {
77 if (dim_ == dim_type(-1))
81 "Convexes of different dimension: to be done");
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();
87 elements[cv].gausspt.resize(pai->nb_points());
90 for (
size_type k = 0; k < pai->nb_points(); ++k) {
91 gausspt_interpolation_data &gpid = elements[cv].gausspt[k];
93 gpt = pgt->transform(pai->point(k),
95 gpid.iflags = find_a_point(gpt, gpid.ptref, gpid.elt) ? 1 : 0;
96 if (gpid.iflags && last_cv != gpid.elt) {
100 if (!(blocked_dof[idof])) dofs.add(idof);
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());
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];
116 gpid.local_dof.resize(nbd);
119 gpid.local_dof[i] = dofs.is_in(ndof) ? ind_dof[ndof]
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_);
131 dof_types_.resize(max_dof);
132 std::fill(dof_types_.begin(), dof_types_.end(),
138 std::fill(ind_dof.begin(), ind_dof.end(),
size_type(-1));
144 return elements[cv].nb_dof;
145 else GMM_ASSERT1(
false,
"Wrong convex number: " << cv);
148 size_type interpolated_fem::index_of_global_dof
150 {
return elements[cv].inddof[i]; }
158 if (mim.linked_mesh().convex_index().is_in(cv))
161 mim.linked_mesh().structure_of_convex(cv)->nb_faces()));
162 else GMM_ASSERT1(
false,
"Wrong convex number: " << cv);
166 { GMM_ASSERT1(
false,
"No base values, real only element."); }
169 { GMM_ASSERT1(
false,
"No grad values, real only element."); }
172 { GMM_ASSERT1(
false,
"No hess values, real only element."); }
174 inline void interpolated_fem::actualize_fictx(
pfem pf,
size_type cv,
175 const base_node &ptr)
const {
176 if (fictx_cv != cv) {
179 (mf.
linked_mesh().trans_of_convex(cv), pf, base_node(), G, cv,
187 base_tensor &t,
bool)
const {
188 elt_interpolation_data &e = elements.at(c.
convex_num());
193 std::fill(t.begin(), t.end(), scalar_type(0));
194 if (e.nb_dof == 0)
return;
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) {
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))
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; }
215 if (find_a_point(c.
xreal(), ptref, 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;
223 for (
size_type i = 0; i < pf->nb_dof(cv); ++i)
228 = taux(i, j*(1-mdim));
240 elt_interpolation_data &e = elements.at(c.
convex_num());
245 std::fill(t.begin(), t.end(), scalar_type(0));
246 if (nbdof == 0)
return;
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) {
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);
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)
268 ee += trans(l, k) * taux(i, j*(1-mdim), l);
269 t(gpid.local_dof[i*rdim+j*mdim], j, k) = ee;
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)
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; }
284 cerr <<
"NON PGP OU MAUVAIS PTS sz=" << elements.size() <<
", cv="
286 cerr <<
"ii=" << c.ii() <<
", sz=" << e.gausspt.size() <<
"\n";
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);
294 ind_dof.at(e.inddof.at(i)) = i;
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)
300 if (ind_dof[mf.ind_basic_dof_of_element(cv)[i*rdim+j*mdim]]
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;
309 for (
size_type i = 0; i < pf->nb_dof(cv); ++i)
310 for (
size_type j = 0; j < target_dim(); ++j)
312 if (ind_dof[mf.ind_basic_dof_of_element(cv)[i*rdim+j*mdim]]
314 t(ind_dof[mf.ind_basic_dof_of_element(cv)[i*rdim+j*mdim]],j,k)
315 = taux(i,j*(1-mdim),k);
325 { GMM_ASSERT1(
false,
"Sorry, to be done."); }
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);
341 scalar_type &meang)
const {
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]++;
350 ming = 100000; maxg = 0; meang = 0;
352 !cv.finished(); ++cv) {
353 ming = std::min(ming, v[cv]);
354 maxg = std::max(maxg, v[cv]);
360 size_type interpolated_fem::memsize()
const {
362 sz += blocked_dof.memsize();
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);
375 interpolated_fem::interpolated_fem(
const mesh_fem &mef,
377 pinterpolated_func pif_,
378 dal::bit_vector blocked_dof_,
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");
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();
389 update_from_context();
392 DAL_SIMPLE_KEY(special_intfem_key,
pfem);
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);