38 #ifndef GETFEM_MESHER_H__ 39 #define GETFEM_MESHER_H__ 54 virtual scalar_type operator()(
const base_node &P)
const = 0;
55 virtual ~mesher_virtual_function() {}
58 class mvf_constant :
public mesher_virtual_function {
61 mvf_constant(scalar_type c_) : c(c_) {}
62 scalar_type operator()(
const base_node &)
const {
return c; }
71 class mesher_signed_distance :
public mesher_virtual_function {
75 mesher_signed_distance() : id(
size_type(-1)) {}
76 virtual ~mesher_signed_distance() {}
77 virtual bool bounding_box(base_node &bmin, base_node &bmax)
const = 0;
78 virtual scalar_type operator()(
const base_node &P,
79 dal::bit_vector &bv)
const = 0;
80 virtual scalar_type grad(
const base_node &P,
81 base_small_vector &G)
const = 0;
82 virtual void hess(
const base_node &P, base_matrix &H)
const = 0;
83 virtual void register_constraints(std::vector<
const 84 mesher_signed_distance*>& list)
const=0;
85 virtual scalar_type operator()(
const base_node &P)
const = 0;
88 typedef std::shared_ptr<const mesher_signed_distance> pmesher_signed_distance;
90 class mesher_half_space :
public mesher_signed_distance {
91 base_node x0; base_small_vector n; scalar_type xon;
93 mesher_half_space() =
default;
94 mesher_half_space(
const base_node &x0_,
const base_small_vector &n_)
97 bool bounding_box(base_node &, base_node &)
const 99 virtual scalar_type operator()(
const base_node &P)
const 100 {
return xon - gmm::vect_sp(P,n); }
101 virtual scalar_type operator()(
const base_node &P,
102 dal::bit_vector &bv)
const {
103 scalar_type d = xon - gmm::vect_sp(P,n);
104 bv[id] = (gmm::abs(d) < SEPS);
107 virtual void register_constraints(std::vector<
const 108 mesher_signed_distance*>& list)
const {
109 id = list.size(); list.push_back(
this);
111 scalar_type grad(
const base_node &P, base_small_vector &G)
const {
112 G = n; G *= scalar_type(-1);
113 return xon - gmm::vect_sp(P,n);
115 void hess(
const base_node &P, base_matrix &H)
const {
121 inline pmesher_signed_distance new_mesher_half_space
122 (
const base_node &x0,
const base_small_vector &n)
123 {
return std::make_shared<mesher_half_space>(x0, n); }
125 class mesher_level_set :
public mesher_signed_distance {
127 mutable std::vector<base_poly> gradient;
128 mutable std::vector<base_poly> hessian;
129 const fem<base_poly> *pf;
130 mutable int initialized;
131 scalar_type shift_ls;
133 bool is_initialized(
void)
const {
return initialized; }
134 mesher_level_set() : initialized(0) {}
135 template <
typename VECT>
136 mesher_level_set(
pfem pf_,
const VECT &coeff_,
137 scalar_type shift_ls_ = scalar_type(0)) {
138 init_base(pf_, coeff_);
139 set_shift(shift_ls_);
141 void set_shift(scalar_type shift_ls_) {
142 shift_ls = shift_ls_;
143 if (shift_ls != scalar_type(0)) {
144 base_node P(pf->dim()); base_small_vector G(pf->dim());
149 template <
typename VECT>
void init_base(
pfem pf_,
const VECT &coeff_);
150 void init_grad(
void)
const;
151 void init_hess(
void)
const;
153 bool bounding_box(base_node &, base_node &)
const 155 virtual scalar_type operator()(
const base_node &P)
const 156 {
return bgeot::to_scalar(base.eval(P.begin())) + shift_ls; }
157 virtual scalar_type operator()(
const base_node &P,
158 dal::bit_vector &bv)
const 159 { scalar_type d = (*this)(P); bv[id] = (gmm::abs(d) < SEPS);
return d; }
160 virtual void register_constraints(std::vector<
const 161 mesher_signed_distance*>& list)
const {
162 id = list.size(); list.push_back(
this);
164 scalar_type grad(
const base_node &P, base_small_vector &G)
const;
165 void hess(
const base_node &P, base_matrix &H)
const;
168 template <
typename VECT>
169 void mesher_level_set::init_base(
pfem pf_,
const VECT &coeff_) {
170 std::vector<scalar_type> coeff(coeff_.begin(), coeff_.end());
171 GMM_ASSERT1(gmm::vect_norm2(coeff) != 0,
"level is zero!");
172 pf =
dynamic_cast<const fem<base_poly>*
> (pf_.get());
173 GMM_ASSERT1(pf,
"A polynomial fem is required for level set (got " 174 <<
typeid(pf_).name() <<
")");
175 base = base_poly(pf->base()[0].dim(), pf->base()[0].degree());
176 for (
unsigned i=0; i < pf->nb_base(0); ++i) {
177 base += pf->base()[i] * coeff[i];
182 template <
typename VECT>
183 inline pmesher_signed_distance new_mesher_level_set
184 (
pfem pf,
const VECT &coeff, scalar_type shift_ls = scalar_type(0))
185 {
return std::make_shared<mesher_level_set>(pf, coeff, shift_ls); }
188 class mesher_ball :
public mesher_signed_distance {
189 base_node x0; scalar_type R;
191 mesher_ball(base_node x0_, scalar_type R_) : x0(x0_), R(R_) {}
192 bool bounding_box(base_node &bmin, base_node &bmax)
const {
194 for (
size_type i=0; i < x0.size(); ++i) { bmin[i] -= R; bmax[i] += R; }
197 virtual scalar_type operator()(
const base_node &P,
198 dal::bit_vector &bv)
const {
200 bv[id] = (gmm::abs(d) < SEPS);
203 virtual scalar_type operator()(
const base_node &P)
const 205 virtual void register_constraints(std::vector<
const 206 mesher_signed_distance*>& list)
const {
207 id = list.size(); list.push_back(
this);
209 scalar_type grad(
const base_node &P, base_small_vector &G)
const {
212 while (e == scalar_type(0))
217 void hess(
const base_node &, base_matrix &)
const {
218 GMM_ASSERT1(
false,
"Sorry, to be done");
222 inline pmesher_signed_distance new_mesher_ball(base_node x0, scalar_type R)
223 {
return std::make_shared<mesher_ball>(x0, R); }
225 class mesher_rectangle :
public mesher_signed_distance {
227 base_node rmin, rmax;
228 std::vector<mesher_half_space> hfs;
230 mesher_rectangle(base_node rmin_, base_node rmax_)
231 : rmin(rmin_), rmax(rmax_), hfs(rmin.size()*2) {
232 base_node n(rmin_.size());
233 for (
unsigned k = 0; k < rmin.size(); ++k) {
235 hfs[k*2] = mesher_half_space(rmin, n);
237 hfs[k*2+1] = mesher_half_space(rmax, n);
241 bool bounding_box(base_node &bmin, base_node &bmax)
const {
242 bmin = rmin; bmax = rmax;
245 virtual scalar_type operator()(
const base_node &P)
const {
247 scalar_type d = rmin[0] - P[0];
249 d = std::max(d, rmin[i] - P[i]);
250 d = std::max(d, P[i] - rmax[i]);
255 virtual scalar_type operator()(
const base_node &P, dal::bit_vector &bv)
257 scalar_type d = this->operator()(P);
258 if (gmm::abs(d) < SEPS)
259 for (
int k = 0; k < 2*rmin.size(); ++k) hfs[k](P, bv);
262 scalar_type grad(
const base_node &P, base_small_vector &G)
const {
263 unsigned i = 0; scalar_type di = hfs[i](P);
264 for (
int k = 1; k < 2*rmin.size(); ++k) {
265 scalar_type dk = hfs[k](P);
266 if (dk > di) { i = k; di = dk; }
268 return hfs[i].grad(P, G);
270 void hess(
const base_node &P, base_matrix &H)
const {
273 virtual void register_constraints(std::vector<
const 274 mesher_signed_distance*>& list)
const {
275 for (
int k = 0; k < 2*rmin.size(); ++k)
276 hfs[k].register_constraints(list);
280 inline pmesher_signed_distance new_mesher_rectangle(base_node rmin,
282 {
return std::make_shared<mesher_rectangle>(rmin, rmax); }
284 class mesher_simplex_ref :
public mesher_signed_distance {
286 std::vector<mesher_half_space> hfs;
290 mesher_simplex_ref(
unsigned N_) : hfs(N_+1), N(N_), org(N_) {
292 for (
unsigned k = 0; k < N; ++k) {
294 hfs[k] = mesher_half_space(org, no);
297 std::fill(org.begin(), org.end(), 1.0/N);
299 hfs[N] = mesher_half_space(org, no);
301 bool bounding_box(base_node &bmin, base_node &bmax)
const {
302 bmin.resize(N); bmax.resize(N);
303 std::fill(bmin.begin(), bmin.end(), scalar_type(0));
304 std::fill(bmax.begin(), bmax.end(), scalar_type(1));
307 virtual scalar_type operator()(
const base_node &P)
const {
308 scalar_type d = - P[0];
309 for (
size_type i=1; i < N; ++i) d = std::max(d, - P[i]);
310 d = std::max(d, gmm::vect_sp(P - org, org) / gmm::vect_norm2(org));
314 virtual scalar_type operator()(
const base_node &P, dal::bit_vector &bv)
316 scalar_type d = this->operator()(P);
317 if (gmm::abs(d) < SEPS)
for (
unsigned k = 0; k < N+1; ++k) hfs[k](P, bv);
320 scalar_type grad(
const base_node &P, base_small_vector &G)
const {
321 unsigned i = 0; scalar_type di = hfs[i](P);
322 for (
unsigned k = 1; k < N+1; ++k) {
323 scalar_type dk = hfs[k](P);
324 if (dk > di) { i = k; di = dk; }
326 return hfs[i].grad(P, G);
328 void hess(
const base_node &P, base_matrix &H)
const {
331 virtual void register_constraints(std::vector<
const 332 mesher_signed_distance*>& list)
const 333 {
for (
unsigned k = 0; k < N+1; ++k) hfs[k].register_constraints(list); }
336 inline pmesher_signed_distance new_mesher_simplex_ref(
unsigned N)
337 {
return std::make_shared<mesher_simplex_ref>(N); }
340 class mesher_prism_ref :
public mesher_signed_distance {
342 std::vector<mesher_half_space> hfs;
346 mesher_prism_ref(
unsigned N_) : hfs(N_+2), N(N_) {
349 for (
unsigned k = 0; k < N; ++k) {
351 hfs[k] = mesher_half_space(org, no);
356 hfs[N] = mesher_half_space(org, no);
357 std::fill(org.begin(), org.end(), 1.0/N);
360 hfs[N+1] = mesher_half_space(org, no);
362 bool bounding_box(base_node &bmin, base_node &bmax)
const {
363 bmin.resize(N); bmax.resize(N);
364 std::fill(bmin.begin(), bmin.end(), scalar_type(0));
365 std::fill(bmax.begin(), bmax.end(), scalar_type(1));
368 virtual scalar_type operator()(
const base_node &P)
const {
369 scalar_type d = - P[0];
370 for (
size_type i=1; i < N; ++i) d = std::max(d, - P[i]);
371 d = std::max(d, P[N-1] - scalar_type(1));
372 d = std::max(d, gmm::vect_sp(P - org, org) / gmm::vect_norm2(org));
376 virtual scalar_type operator()(
const base_node &P, dal::bit_vector &bv)
378 scalar_type d = this->operator()(P);
379 if (gmm::abs(d) < SEPS)
for (
unsigned k = 0; k < N+2; ++k) hfs[k](P, bv);
382 scalar_type grad(
const base_node &P, base_small_vector &G)
const {
383 unsigned i = 0; scalar_type di = hfs[i](P);
384 for (
unsigned k = 1; k < N+2; ++k) {
385 scalar_type dk = hfs[k](P);
386 if (dk > di) { i = k; di = dk; }
388 return hfs[i].grad(P, G);
390 void hess(
const base_node &P, base_matrix &H)
const {
393 virtual void register_constraints(std::vector<
const 394 mesher_signed_distance*>& list)
const 395 {
for (
unsigned k = 0; k < N+2; ++k) hfs[k].register_constraints(list); }
398 inline pmesher_signed_distance new_mesher_prism_ref(
unsigned N)
399 {
return std::make_shared<mesher_prism_ref>(N); }
408 class mesher_union :
public mesher_signed_distance {
409 std::vector<pmesher_signed_distance> dists;
410 mutable std::vector<scalar_type> vd;
414 mesher_union(
const std::vector<pmesher_signed_distance>
415 &dists_) : dists(dists_)
416 { vd.resize(dists.size()); with_min =
true; }
419 (
const pmesher_signed_distance &a,
420 const pmesher_signed_distance &b,
421 const pmesher_signed_distance &c = pmesher_signed_distance(),
422 const pmesher_signed_distance &d = pmesher_signed_distance(),
423 const pmesher_signed_distance &e = pmesher_signed_distance(),
424 const pmesher_signed_distance &f = pmesher_signed_distance(),
425 const pmesher_signed_distance &g = pmesher_signed_distance(),
426 const pmesher_signed_distance &h = pmesher_signed_distance(),
427 const pmesher_signed_distance &i = pmesher_signed_distance(),
428 const pmesher_signed_distance &j = pmesher_signed_distance(),
429 const pmesher_signed_distance &k = pmesher_signed_distance(),
430 const pmesher_signed_distance &l = pmesher_signed_distance(),
431 const pmesher_signed_distance &m = pmesher_signed_distance(),
432 const pmesher_signed_distance &n = pmesher_signed_distance(),
433 const pmesher_signed_distance &o = pmesher_signed_distance(),
434 const pmesher_signed_distance &p = pmesher_signed_distance(),
435 const pmesher_signed_distance &q = pmesher_signed_distance(),
436 const pmesher_signed_distance &r = pmesher_signed_distance(),
437 const pmesher_signed_distance &s = pmesher_signed_distance(),
438 const pmesher_signed_distance &t = pmesher_signed_distance()) {
439 dists.push_back(a); dists.push_back(b);
441 if (c) dists.push_back(c);
442 if (d) dists.push_back(d);
443 if (e) dists.push_back(e);
444 if (f) dists.push_back(f);
445 if (g) dists.push_back(g);
446 if (h) dists.push_back(h);
447 if (i) dists.push_back(i);
448 if (j) dists.push_back(j);
449 if (k) dists.push_back(k);
450 if (l) dists.push_back(l);
451 if (m) dists.push_back(m);
452 if (n) dists.push_back(n);
453 if (o) dists.push_back(o);
454 if (p) dists.push_back(p);
455 if (q) dists.push_back(q);
456 if (r) dists.push_back(r);
457 if (s) dists.push_back(s);
458 if (t) dists.push_back(t);
459 vd.resize(dists.size());
462 bool bounding_box(base_node &bmin, base_node &bmax)
const {
463 base_node bmin2, bmax2;
464 bool b = dists[0]->bounding_box(bmin, bmax);
465 if (!b)
return false;
466 for (
size_type k = 1; k < dists.size(); ++k) {
467 b = dists[k]->bounding_box(bmin2, bmax2);
468 if (!b)
return false;
469 for (
unsigned i=0; i < bmin.size(); ++i) {
470 bmin[i] = std::min(bmin[i],bmin2[i]);
471 bmax[i] = std::max(bmax[i],bmax2[i]);
476 virtual scalar_type operator()(
const base_node &P)
const {
477 scalar_type d, f(0), g(1);
479 d = (*(dists[0]))(P);
480 for (
size_type k = 1; k < dists.size(); ++k)
481 d = std::min(d, (*(dists[k]))(P));
485 for (
size_type k = 0; k < dists.size(); ++k) {
486 vd[k] = (*(dists[k]))(P);
487 if (vd[k] <= scalar_type(0)) isin =
true;
488 f += gmm::sqr(gmm::neg(vd[k]));
489 g *= gmm::pos(vd[k]);
491 d = isin ? -gmm::sqrt(f)
492 : pow(g, scalar_type(1) / scalar_type(dists.size()));
496 scalar_type operator()(
const base_node &P, dal::bit_vector &bv)
const {
498 scalar_type d = vd[0] = (*(dists[0]))(P);
499 bool ok = (d > -SEPS);
500 for (
size_type k = 1; k < dists.size(); ++k) {
501 vd[k] = (*(dists[k]))(P);
if (vd[k] <= -SEPS) ok =
false;
502 d = std::min(d,vd[k]);
504 for (
size_type k = 0; ok && k < dists.size(); ++k) {
505 if (vd[k] < SEPS) (*(dists[k]))(P, bv);
510 vd[0] = (*(dists[0]))(P);
511 bool ok = (vd[0] > -SEPS);
512 for (
size_type k = 1; k < dists.size(); ++k) {
513 vd[k] = (*(dists[k]))(P);
if (vd[k] <= -SEPS) ok =
false;
515 for (
size_type k = 0; ok && k < dists.size(); ++k) {
516 if (vd[k] < SEPS) (*(dists[k]))(P, bv);
518 return operator()(P);
521 virtual void register_constraints(std::vector<
const 522 mesher_signed_distance*>& list)
const {
523 for (
size_type k = 0; k < dists.size(); ++k)
524 dists[k]->register_constraints(list);
526 scalar_type grad(
const base_node &P, base_small_vector &G)
const {
529 d = (*(dists[0]))(P);
531 for (
size_type k = 1; k < dists.size(); ++k) {
532 scalar_type d2 = (*(dists[k]))(P);
533 if (d2 < d) { d = d2; i = k; }
535 return dists[i]->grad(P, G);
539 base_small_vector Gloc;
540 for (
size_type k = 0; k < dists.size(); ++k) {
541 dists[k]->grad(P, Gloc);
543 Gloc *= -gmm::neg(vd[k]);
545 Gloc *= pow(d, scalar_type(dists.size())) / vd[k];
546 if (!k) G = Gloc;
else G += Gloc;
549 G *= scalar_type(1) / d;
551 G /= pow(d, scalar_type(dists.size()-1)) * scalar_type(dists.size());
555 void hess(
const base_node &P, base_matrix &H)
const {
556 scalar_type d = (*(dists[0]))(P);
557 if (with_min || gmm::abs(d) < SEPS) {
559 for (
size_type k = 1; k < dists.size(); ++k) {
560 scalar_type d2 = (*(dists[k]))(P);
561 if (d2 < d) { d = d2; i = k; }
563 dists[i]->hess(P, H);
566 GMM_ASSERT1(
false,
"Sorry, to be done");
571 inline pmesher_signed_distance new_mesher_union
572 (
const std::vector<pmesher_signed_distance> &dists)
573 {
return std::make_shared<mesher_union>(dists); }
575 inline pmesher_signed_distance new_mesher_union
576 (
const pmesher_signed_distance &a,
577 const pmesher_signed_distance &b,
578 const pmesher_signed_distance &c = pmesher_signed_distance(),
579 const pmesher_signed_distance &d = pmesher_signed_distance(),
580 const pmesher_signed_distance &e = pmesher_signed_distance(),
581 const pmesher_signed_distance &f = pmesher_signed_distance(),
582 const pmesher_signed_distance &g = pmesher_signed_distance(),
583 const pmesher_signed_distance &h = pmesher_signed_distance(),
584 const pmesher_signed_distance &i = pmesher_signed_distance(),
585 const pmesher_signed_distance &j = pmesher_signed_distance(),
586 const pmesher_signed_distance &k = pmesher_signed_distance(),
587 const pmesher_signed_distance &l = pmesher_signed_distance(),
588 const pmesher_signed_distance &m = pmesher_signed_distance(),
589 const pmesher_signed_distance &n = pmesher_signed_distance(),
590 const pmesher_signed_distance &o = pmesher_signed_distance(),
591 const pmesher_signed_distance &p = pmesher_signed_distance(),
592 const pmesher_signed_distance &q = pmesher_signed_distance(),
593 const pmesher_signed_distance &r = pmesher_signed_distance(),
594 const pmesher_signed_distance &s = pmesher_signed_distance(),
595 const pmesher_signed_distance &t = pmesher_signed_distance()) {
596 return std::make_shared<mesher_union>
597 (a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s,t);
602 class mesher_intersection :
public mesher_signed_distance {
603 std::vector<pmesher_signed_distance> dists;
604 mutable std::vector<scalar_type> vd;
609 mesher_intersection(
const std::vector<pmesher_signed_distance>
610 &dists_) : dists(dists_)
611 { vd.resize(dists.size()); }
614 (
const pmesher_signed_distance &a,
615 const pmesher_signed_distance &b,
616 const pmesher_signed_distance &c = pmesher_signed_distance(),
617 const pmesher_signed_distance &d = pmesher_signed_distance(),
618 const pmesher_signed_distance &e = pmesher_signed_distance(),
619 const pmesher_signed_distance &f = pmesher_signed_distance(),
620 const pmesher_signed_distance &g = pmesher_signed_distance(),
621 const pmesher_signed_distance &h = pmesher_signed_distance(),
622 const pmesher_signed_distance &i = pmesher_signed_distance(),
623 const pmesher_signed_distance &j = pmesher_signed_distance(),
624 const pmesher_signed_distance &k = pmesher_signed_distance(),
625 const pmesher_signed_distance &l = pmesher_signed_distance(),
626 const pmesher_signed_distance &m = pmesher_signed_distance(),
627 const pmesher_signed_distance &n = pmesher_signed_distance(),
628 const pmesher_signed_distance &o = pmesher_signed_distance(),
629 const pmesher_signed_distance &p = pmesher_signed_distance(),
630 const pmesher_signed_distance &q = pmesher_signed_distance(),
631 const pmesher_signed_distance &r = pmesher_signed_distance(),
632 const pmesher_signed_distance &s = pmesher_signed_distance(),
633 const pmesher_signed_distance &t = pmesher_signed_distance()) {
634 dists.push_back(a); dists.push_back(b);
635 if (c) dists.push_back(c);
636 if (d) dists.push_back(d);
637 if (e) dists.push_back(e);
638 if (f) dists.push_back(f);
639 if (g) dists.push_back(g);
640 if (h) dists.push_back(h);
641 if (i) dists.push_back(i);
642 if (j) dists.push_back(j);
643 if (k) dists.push_back(k);
644 if (l) dists.push_back(l);
645 if (m) dists.push_back(m);
646 if (n) dists.push_back(n);
647 if (o) dists.push_back(o);
648 if (p) dists.push_back(p);
649 if (q) dists.push_back(q);
650 if (r) dists.push_back(r);
651 if (s) dists.push_back(s);
652 if (t) dists.push_back(t);
653 vd.resize(dists.size());
655 bool bounding_box(base_node &bmin, base_node &bmax)
const {
656 base_node bmin2, bmax2;
658 bool b = dists[0]->bounding_box(bmin, bmax); first = !b;
659 for (
size_type k = 1; k < dists.size(); ++k) {
660 bool bb = dists[k]->bounding_box(bmin2, bmax2);
661 for (
unsigned i=0; i < bmin.size() && bb && !first; ++i) {
662 bmin[i] = std::max(bmin[i],bmin2[i]);
663 bmax[i] = std::max(std::min(bmax[i],bmax2[i]), bmin[i]);
665 if (first && bb) { bmin = bmin2; bmax = bmax2; first =
false; }
670 virtual scalar_type operator()(
const base_node &P)
const {
671 scalar_type d = (*(dists[0]))(P);
672 for (
size_type k = 1; k < dists.size(); ++k)
673 d = std::max(d, (*(dists[k]))(P));
677 scalar_type operator()(
const base_node &P, dal::bit_vector &bv)
const {
678 scalar_type d = vd[0] = (*(dists[0]))(P);
679 bool ok = (d < SEPS);
680 for (
size_type k = 1; k < dists.size(); ++k) {
681 vd[k] = (*(dists[k]))(P);
if (vd[k] >= SEPS) ok =
false;
682 d = std::min(d,vd[k]);
684 for (
size_type k = 0; ok && k < dists.size(); ++k) {
685 if (vd[k] > -SEPS) (*(dists[k]))(P, bv);
689 virtual void register_constraints(std::vector<
const 690 mesher_signed_distance*>& list)
const {
691 for (
size_type k = 0; k < dists.size(); ++k) {
692 dists[k]->register_constraints(list);
695 scalar_type grad(
const base_node &P, base_small_vector &G)
const {
696 scalar_type d = (*(dists[0]))(P);
698 for (
size_type k = 1; k < dists.size(); ++k) {
699 scalar_type d2 = (*(dists[k]))(P);
700 if (d2 > d) { d = d2; i = k; }
702 return dists[i]->grad(P, G);
704 void hess(
const base_node &P, base_matrix &H)
const {
705 scalar_type d = (*(dists[0]))(P);
707 for (
size_type k = 1; k < dists.size(); ++k) {
708 scalar_type d2 = (*(dists[k]))(P);
709 if (d2 > d) { d = d2; i = k; }
711 dists[i]->hess(P, H);
715 inline pmesher_signed_distance new_mesher_intersection
716 (
const std::vector<pmesher_signed_distance> &dists)
717 {
return std::make_shared<mesher_intersection>(dists); }
719 inline pmesher_signed_distance new_mesher_intersection
720 (
const pmesher_signed_distance &a,
721 const pmesher_signed_distance &b,
722 const pmesher_signed_distance &c = pmesher_signed_distance(),
723 const pmesher_signed_distance &d = pmesher_signed_distance(),
724 const pmesher_signed_distance &e = pmesher_signed_distance(),
725 const pmesher_signed_distance &f = pmesher_signed_distance(),
726 const pmesher_signed_distance &g = pmesher_signed_distance(),
727 const pmesher_signed_distance &h = pmesher_signed_distance(),
728 const pmesher_signed_distance &i = pmesher_signed_distance(),
729 const pmesher_signed_distance &j = pmesher_signed_distance(),
730 const pmesher_signed_distance &k = pmesher_signed_distance(),
731 const pmesher_signed_distance &l = pmesher_signed_distance(),
732 const pmesher_signed_distance &m = pmesher_signed_distance(),
733 const pmesher_signed_distance &n = pmesher_signed_distance(),
734 const pmesher_signed_distance &o = pmesher_signed_distance(),
735 const pmesher_signed_distance &p = pmesher_signed_distance(),
736 const pmesher_signed_distance &q = pmesher_signed_distance(),
737 const pmesher_signed_distance &r = pmesher_signed_distance(),
738 const pmesher_signed_distance &s = pmesher_signed_distance(),
739 const pmesher_signed_distance &t = pmesher_signed_distance()) {
740 return std::make_shared<mesher_intersection>
741 (a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s,t);
746 class mesher_setminus :
public mesher_signed_distance {
747 pmesher_signed_distance a, b;
749 mesher_setminus(
const pmesher_signed_distance& a_,
750 const pmesher_signed_distance &b_) : a(a_), b(b_) {}
751 bool bounding_box(base_node &bmin, base_node &bmax)
const 752 {
return a->bounding_box(bmin,bmax); }
753 scalar_type operator()(
const base_node &P, dal::bit_vector &bv)
const {
754 scalar_type da = (*a)(P), db = -(*b)(P);
755 if (da < SEPS && db < SEPS) {
756 if (da > -SEPS) (*a)(P, bv);
757 if (db > -SEPS) (*b)(P, bv);
759 return std::max(da, db);
761 scalar_type operator()(
const base_node &P)
const 762 {
return std::max((*a)(P),-(*b)(P)); }
763 virtual void register_constraints(std::vector<
const 764 mesher_signed_distance*>& list)
const {
765 a->register_constraints(list); b->register_constraints(list);
767 scalar_type grad(
const base_node &P, base_small_vector &G)
const {
768 scalar_type da = (*a)(P), db = -(*b)(P);
769 if (da > db)
return a->grad(P, G);
770 else { b->grad(P, G); G *= scalar_type(-1);
return db; }
772 void hess(
const base_node &P, base_matrix &H)
const {
773 scalar_type da = (*a)(P), db = -(*b)(P);
774 if (da > db) a->hess(P, H);
775 else { b->hess(P, H); gmm::scale(H, scalar_type(-1)); }
780 inline pmesher_signed_distance new_mesher_setminus
781 (
const pmesher_signed_distance &a,
const pmesher_signed_distance &b)
782 {
return std::make_shared<mesher_setminus>(a, b); }
785 class mesher_tube :
public mesher_signed_distance {
786 base_node x0; base_node n; scalar_type R;
788 mesher_tube(base_node x0_, base_node n_, scalar_type R_)
789 : x0(x0_), n(n_), R(R_)
791 bool bounding_box(base_node &, base_node &)
const 793 virtual scalar_type operator()(
const base_node &P)
const {
794 base_node v(P); v -= x0;
795 gmm::add(gmm::scaled(n, -gmm::vect_sp(v, n)), v);
798 virtual scalar_type operator()(
const base_node &P,
799 dal::bit_vector &bv)
const {
800 scalar_type d = (*this)(P);
801 bv[id] = (gmm::abs(d) < SEPS);
804 virtual void register_constraints(std::vector<
const 805 mesher_signed_distance*>& list)
const {
806 id = list.size(); list.push_back(
this);
808 scalar_type grad(
const base_node &P, base_small_vector &G)
const {
810 gmm::add(gmm::scaled(n, -gmm::vect_sp(G, n)), G);
812 while (e == scalar_type(0)) {
814 gmm::add(gmm::scaled(n, -gmm::vect_sp(G, n)), G);
820 void hess(
const base_node &, base_matrix &)
const {
821 GMM_ASSERT1(
false,
"Sorry, to be done");
825 inline pmesher_signed_distance new_mesher_tube(base_node x0, base_node n,
827 {
return std::make_shared<mesher_tube>(x0,n,R); }
830 class mesher_cylinder :
public mesher_signed_distance {
831 base_node x0; base_small_vector n;
833 pmesher_signed_distance t, p1, p2, i1;
835 mesher_cylinder(
const base_node &c,
const base_small_vector &no,
836 scalar_type L_, scalar_type R_)
838 t(new_mesher_tube(x0, n, R)),
839 p1(new_mesher_half_space(x0, n)),
840 p2(new_mesher_half_space(x0+n*L, -1.0 * n)),
841 i1(new_mesher_intersection(p1, p2, t)) {}
842 bool bounding_box(base_node &bmin, base_node &bmax)
const {
843 base_node x1(x0+n*L);
845 for (
unsigned i = 0; i < gmm::vect_size(x0); ++i) {
846 bmin[i] = std::min(x0[i], x1[i]) - R;
847 bmax[i] = std::max(x0[i], x1[i]) + R;
851 virtual scalar_type operator()(
const base_node &P)
const {
return (*i1)(P); }
852 virtual scalar_type operator()(
const base_node &P,
853 dal::bit_vector& bv)
const 854 {
return (*i1)(P, bv); }
855 scalar_type grad(
const base_node &P, base_small_vector &G)
const 856 {
return i1->grad(P, G); }
857 void hess(
const base_node &, base_matrix &)
const {
858 GMM_ASSERT1(
false,
"Sorry, to be done");
860 virtual void register_constraints(std::vector<
const 861 mesher_signed_distance*>& list)
const 862 { i1->register_constraints(list); }
865 inline pmesher_signed_distance new_mesher_cylinder
866 (
const base_node &c,
const base_small_vector &no,
867 scalar_type L, scalar_type R)
868 {
return std::make_shared<mesher_cylinder>(c,no,L,R); }
871 class mesher_infinite_cone :
public mesher_signed_distance {
873 base_node x0; base_node n; scalar_type
alpha;
875 mesher_infinite_cone(base_node x0_, base_node n_, scalar_type alpha_)
876 : x0(x0_), n(n_), alpha(alpha_)
878 bool bounding_box(base_node &, base_node &)
const 880 virtual scalar_type operator()(
const base_node &P)
const {
881 base_node v(P); v -= x0;
882 scalar_type v_n = gmm::vect_sp(v, n);
883 gmm::add(gmm::scaled(n, -v_n), v);
886 virtual scalar_type operator()(
const base_node &P,
887 dal::bit_vector &bv)
const {
888 scalar_type d = (*this)(P);
889 bv[id] = (gmm::abs(d) < SEPS);
892 virtual void register_constraints(std::vector<
const 893 mesher_signed_distance*>& list)
const {
894 id = list.size(); list.push_back(
this);
896 scalar_type grad(
const base_node &P, base_small_vector &v)
const {
898 scalar_type v_n = gmm::vect_sp(v, n);
899 gmm::add(gmm::scaled(n, -v_n), v);
901 scalar_type d = no * cos(alpha) - gmm::abs(v_n) * sin(alpha);
902 while (no == scalar_type(0)) {
904 gmm::add(gmm::scaled(n, -gmm::vect_sp(v, n)), v);
907 v *= cos(alpha) / no;
908 v -= (sin(alpha) * gmm::sgn(v_n)) * n;
911 void hess(
const base_node &, base_matrix &)
const {
912 GMM_ASSERT1(
false,
"Sorry, to be done");
916 inline pmesher_signed_distance new_mesher_infinite_cone
917 (base_node x0, base_node n, scalar_type alpha)
918 {
return std::make_shared<mesher_infinite_cone>(x0,n,alpha); }
921 class mesher_cone :
public mesher_signed_distance {
922 base_node x0; base_small_vector n;
923 scalar_type L, alpha;
924 pmesher_signed_distance t, p1, p2, i1;
926 mesher_cone(
const base_node &c,
const base_small_vector &no,
927 scalar_type L_, scalar_type alpha_)
929 t(new_mesher_infinite_cone(x0, n, alpha)),
930 p1(new_mesher_half_space(x0, n)),
931 p2(new_mesher_half_space(x0+n*L, -1.0 * n)),
932 i1(new_mesher_intersection(p1, p2, t)) {}
933 bool bounding_box(base_node &bmin, base_node &bmax)
const {
934 base_node x1(x0+n*L);
935 scalar_type R = L * sin(alpha);
937 for (
unsigned i = 0; i < gmm::vect_size(x0); ++i) {
938 bmin[i] = std::min(x0[i], x1[i]) - R;
939 bmax[i] = std::max(x0[i], x1[i]) + R;
943 virtual scalar_type operator()(
const base_node &P)
const {
return (*i1)(P); }
944 virtual scalar_type operator()(
const base_node &P,
945 dal::bit_vector& bv)
const 946 {
return (*i1)(P, bv); }
947 scalar_type grad(
const base_node &P, base_small_vector &G)
const 948 {
return i1->grad(P, G); }
949 void hess(
const base_node &, base_matrix &)
const {
950 GMM_ASSERT1(
false,
"Sorry, to be done");
952 virtual void register_constraints(std::vector<
const 953 mesher_signed_distance*>& list)
const 954 { i1->register_constraints(list); }
957 inline pmesher_signed_distance new_mesher_cone
958 (
const base_node &c,
const base_small_vector &no,
959 scalar_type L, scalar_type alpha)
960 {
return std::make_shared<mesher_cone>(c,no,L,alpha); }
963 class mesher_ellipse :
public mesher_signed_distance {
964 base_node x0; base_small_vector n, t;
967 mesher_ellipse(
const base_node ¢er,
const base_small_vector &no,
968 scalar_type r_, scalar_type R_)
970 t[0] = -n[1]; t[1] = n[0];
971 if (R < r) { std::swap(r, R); std::swap(n, t); }
974 bool bounding_box(base_node &bmin, base_node &bmax)
const {
976 for (
unsigned i = 0; i < 2; ++i) { bmin[i] -= R; bmax[i] += R; }
979 virtual scalar_type operator()(
const base_node &P)
const {
980 base_small_vector v(P); v -= x0;
981 scalar_type vt = gmm::vect_sp(v, t);
982 vt = std::max(-a, std::min(a, vt));
983 base_node x1 = x0 + t*vt;
984 base_small_vector v1(P); v1 -= x1;
985 scalar_type v1n = gmm::vect_sp(v1, n), v1t = gmm::vect_sp(v1, t);
986 scalar_type x1n = gmm::vect_sp(x1, n), x1t = gmm::vect_sp(x1, t);
987 scalar_type ea = v1n*v1n / (r*r) + v1t * v1t / (R*R);
988 scalar_type eb = 2. * (x1n*v1n / (r*r) + x1t*v1t / (R*R));
989 scalar_type ec = x1n*x1n / (r*r) + x1t * x1t / (R*R);
991 scalar_type delta = eb*eb - 4 * ea * ec;
993 scalar_type lambda = (-eb + sqrt(delta)) / (2. * ea);
996 virtual scalar_type operator()(
const base_node &P,
997 dal::bit_vector& bv)
const {
998 scalar_type d = this->operator()(P);
999 bv[id] = (gmm::abs(d) < SEPS);
1002 scalar_type grad(
const base_node &, base_small_vector &)
const 1003 { GMM_ASSERT1(
false,
"Sorry, to be done");
return 0.; }
1004 void hess(
const base_node &, base_matrix &)
const {
1005 GMM_ASSERT1(
false,
"Sorry, to be done");
1007 virtual void register_constraints(std::vector<
const 1008 mesher_signed_distance*>& list)
const 1009 {
id = list.size(); list.push_back(
this); }
1012 inline pmesher_signed_distance new_mesher_ellipse
1013 (
const base_node ¢er,
const base_small_vector &no,
1014 scalar_type r, scalar_type R)
1015 {
return std::make_shared<mesher_ellipse>(center,no,r,R); }
1019 class mesher_torus :
public mesher_signed_distance {
1023 mesher_torus(scalar_type RR, scalar_type rr) : R(RR), r(rr) {}
1024 bool bounding_box(base_node &bmin, base_node &bmax)
const {
1025 bmin = base_node(3); bmax = base_node(3);
1026 bmin[0] = bmin[1] = -R-r; bmin[2] = -r;
1027 bmax[0] = bmax[1] = +R+r; bmax[2] = +r;
1030 virtual scalar_type operator()(
const base_node &P)
const {
1031 scalar_type x = P[0], y = P[1], z = P[2], c = sqrt(x*x + y*y);
1032 return (c == 0.) ? R - r : sqrt(gmm::sqr(c-R) + z*z) - r;
1034 virtual scalar_type operator()(
const base_node &P, dal::bit_vector&bv)
1036 scalar_type d = this->operator()(P);
1037 bv[id] = (gmm::abs(d) < SEPS);
1040 scalar_type grad(
const base_node &P, base_small_vector &G)
const {
1042 scalar_type x = P[0], y = P[1], z = P[2], c = sqrt(x*x + y*y), d(0);
1048 scalar_type w = 1. - R / c, e = sqrt(gmm::sqr(c-R) + z*z);
1054 G[0] = x * w / e; G[1] = y * w / e; G[2] = z / e;
1059 void hess(
const base_node &, base_matrix &)
const {
1060 GMM_ASSERT1(
false,
"Sorry, to be done");
1062 virtual void register_constraints(std::vector<
const 1063 mesher_signed_distance*>&list)
const 1064 {
id = list.size(); list.push_back(
this); }
1068 inline pmesher_signed_distance new_mesher_torus(scalar_type R, scalar_type r)
1069 {
return std::make_shared<mesher_torus>(R,r); }
1072 void build_mesh(mesh &m,
const pmesher_signed_distance& dist_,
1073 scalar_type h0,
const std::vector<base_node> &fixed_points
1074 = std::vector<base_node>(),
size_type K = 1,
int noise = -1,
1075 size_type iter_max = 500,
int prefind = 1,
1076 scalar_type dist_point_hull = 4,
1077 scalar_type boundary_threshold_flatness = 0.11);
1080 bool try_projection(
const mesher_signed_distance& dist, base_node &X,
1081 bool on_surface =
false);
1082 bool pure_multi_constraint_projection
1083 (
const std::vector<const mesher_signed_distance*> &list_constraints,
1084 base_node &X,
const dal::bit_vector &cts);
1085 scalar_type curvature_radius_estimate(
const mesher_signed_distance &dist,
1086 base_node X,
bool proj =
false);
1087 scalar_type min_curvature_radius_estimate
1088 (
const std::vector<const mesher_signed_distance*> &list_constraints,
1089 const base_node &X,
const dal::bit_vector &cts,
size_type hide_first = 0);
Define a getfem::getfem_mesh object.
base class for static stored objects
Export solutions to various formats.
computation of the condition number of dense matrices.
void fill_random(L &l, double cfill)
*/
number_traits< typename linalg_traits< V1 >::value_type >::magnitude_type vect_dist2(const V1 &v1, const V2 &v2)
Euclidean distance between two vectors.
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norm2(const V &v)
Euclidean norm of a vector.
size_t size_type
used as the common size type in the library
GEneric Tool for Finite Element Methods.
std::shared_ptr< const getfem::virtual_fem > pfem
type of pointer on a fem description
void clear(L &l)
clear (fill with zeros) a vector or matrix.
convenient initialization of vectors via overload of "operator,".
Implements BFGS (Broyden, Fletcher, Goldfarb, Shanno) algorithm.
Simple implementation of a KD-tree.
size_type alpha(short_type n, short_type d)
Return the value of which is the number of monomials of a polynomial of variables and degree ...
void resize(M &v, size_type m, size_type n)
*/