This file is indexed.

/usr/include/getfem/getfem_mesher.h is in libgetfem++-dev 4.2.1~beta1~svn4482~dfsg-2build1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
/* -*- c++ -*- (enables emacs c++ mode) */
/*===========================================================================
 
 Copyright (C) 2004-2012 Julien Pommier, Yves Renard
 
 This file is a part of GETFEM++
 
 Getfem++  is  free software;  you  can  redistribute  it  and/or modify it
 under  the  terms  of the  GNU  Lesser General Public License as published
 by  the  Free Software Foundation;  either version 3 of the License,  or
 (at your option) any later version along with the GCC Runtime Library
 Exception either version 3.1 or (at your option) any later version.
 This program  is  distributed  in  the  hope  that it will be useful,  but
 WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
 or  FITNESS  FOR  A PARTICULAR PURPOSE.  See the GNU Lesser General Public
 License and GCC Runtime Library Exception for more details.
 You  should  have received a copy of the GNU Lesser General Public License
 along  with  this program;  if not, write to the Free Software Foundation,
 Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301, USA.
 
 As a special exception, you  may use  this file  as it is a part of a free
 software  library  without  restriction.  Specifically,  if   other  files
 instantiate  templates  or  use macros or inline functions from this file,
 or  you compile this  file  and  link  it  with other files  to produce an
 executable, this file  does  not  by itself cause the resulting executable
 to be covered  by the GNU Lesser General Public License.  This   exception
 does not  however  invalidate  any  other  reasons why the executable file
 might be covered by the GNU Lesser General Public License.
 
===========================================================================*/

/**@file getfem_mesher.h
   @author  Julien Pommier <Julien.Pommier@insa-toulouse.fr>, Yves Renard <Yves.Renard@insa-lyon.fr>
   @date May 1, 2004.
   @brief a broken mesher.
*/

#ifndef GETFEM_MESHER_H__
#define GETFEM_MESHER_H__


#include "getfem_mesh.h"
#include "gmm/gmm_condition_number.h"
#include "bgeot_comma_init.h"
#include "gmm/gmm_solver_bfgs.h"
#include "getfem_export.h"
#include "bgeot_kdtree.h"
#include <typeinfo>

namespace getfem {

  class mesher_virtual_function {
  public:
    virtual scalar_type operator()(const base_node &P) const = 0;
    virtual ~mesher_virtual_function() {}
  };

  class mvf_constant : public mesher_virtual_function {
    scalar_type c;
  public: 
    mvf_constant(scalar_type c_) : c(c_) {}
    scalar_type operator()(const base_node &) const { return c; }
  };

  // Signed distance definition. Should not be a real signed distance but
  // should satisfy that dist(P) is less than the euclidean distance from
  // P to the boundary and more than 1/sqrt(N) time this distance.

#define SEPS 1e-8

  class mesher_signed_distance : public mesher_virtual_function {
  protected:
    mutable size_type id;
  public:
    mesher_signed_distance() : id(size_type(-1)) {}
    virtual bool bounding_box(base_node &bmin, base_node &bmax) const = 0;
    virtual scalar_type operator()(const base_node &P,
				   dal::bit_vector &bv) const = 0;
    virtual scalar_type grad(const base_node &P,
			     base_small_vector &G) const = 0;
    virtual void hess(const base_node &P, base_matrix &H) const = 0;
    virtual void register_constraints(std::vector<const
				      mesher_signed_distance*>& list) const=0;
    virtual scalar_type operator()(const base_node &P) const  = 0;
  };

  class mesher_half_space : public mesher_signed_distance {
    base_node x0; base_small_vector n; scalar_type xon;
  public:
    mesher_half_space(const base_node &x0_, const base_small_vector &n_)
      : x0(x0_), n(n_)
    { n /= gmm::vect_norm2(n); xon = gmm::vect_sp(x0, n); }
    bool bounding_box(base_node &, base_node &) const
    { return false; }
    virtual scalar_type operator()(const base_node &P) const
    { return xon - gmm::vect_sp(P,n); }
    virtual scalar_type operator()(const base_node &P,
				   dal::bit_vector &bv) const {
      scalar_type d = xon - gmm::vect_sp(P,n);
      bv[id] = (gmm::abs(d) < SEPS);
      return d;
    }
    virtual void register_constraints(std::vector<const
				      mesher_signed_distance*>& list) const {
      id = list.size(); list.push_back(this);
    }
    scalar_type grad(const base_node &P, base_small_vector &G) const {
      G = n; G *= scalar_type(-1); 
      return xon - gmm::vect_sp(P,n);
    }
    void hess(const base_node &P, base_matrix &H) const {
      gmm::resize(H, P.size(), P.size()); gmm::clear(H);
    }

  };

  class mesher_level_set : public mesher_signed_distance {
    bgeot::base_poly base;
    mutable std::vector<base_poly> gradient;
    mutable std::vector<base_poly> hessian;
    const fem<base_poly> *pf;
    mutable int initialized;
    scalar_type shift_ls;     // for the computation of a gap on a level_set.
  public:
    bool is_initialized(void) const { return initialized; }
    mesher_level_set() : initialized(0) {}
    template <typename VECT>
    mesher_level_set(pfem pf_, const VECT &coeff_,
		     scalar_type shift_ls_ = scalar_type(0)) {
      shift_ls = shift_ls_; init_base(pf_, coeff_);
      if (shift_ls != scalar_type(0)) {
	base_node P(pf->dim()); base_small_vector G(pf->dim());
	grad(P, G);
	shift_ls *= gmm::vect_norm2(G);
      }
    }
    template <typename VECT> void init_base(pfem pf_, const VECT &coeff_);
    void init_grad(void) const;
    void init_hess(void) const;

    bool bounding_box(base_node &, base_node &) const
    { return false; }
    virtual scalar_type operator()(const base_node &P) const
    {  return bgeot::to_scalar(base.eval(P.begin())) + shift_ls; }
    virtual scalar_type operator()(const base_node &P,
				   dal::bit_vector &bv) const
    { scalar_type d = (*this)(P); bv[id] = (gmm::abs(d) < SEPS); return d; }
    virtual void register_constraints(std::vector<const
				      mesher_signed_distance*>& list) const {
      id = list.size(); list.push_back(this);
    }
    scalar_type grad(const base_node &P, base_small_vector &G) const;
    void hess(const base_node &P, base_matrix &H) const;
  };

  template <typename VECT> 
  void mesher_level_set::init_base(pfem pf_, const VECT &coeff_) {
    std::vector<scalar_type> coeff(coeff_.begin(), coeff_.end());
    GMM_ASSERT1(gmm::vect_norm2(coeff) != 0, "level is zero!");
    pf = dynamic_cast<const fem<base_poly>* > (pf_.get());
    GMM_ASSERT1(pf, "PK fem are required for level set (got "
		<< typeid(pf_).name() << ")");
    base = base_poly(pf->base()[0].dim(), pf->base()[0].degree());
    for (unsigned i=0; i < pf->nb_base(0); ++i) {
      base += pf->base()[i] * coeff[i];
    }
    initialized = 0; 
  }




  class mesher_ball : public mesher_signed_distance {
    base_node x0; scalar_type R;
  public:
    mesher_ball(base_node x0_, scalar_type R_) : x0(x0_), R(R_) {}
    bool bounding_box(base_node &bmin, base_node &bmax) const { 
      bmin = bmax = x0; 
      for (size_type i=0; i < x0.size(); ++i) { bmin[i] -= R; bmax[i] += R; }
      return true;
    }
    virtual scalar_type operator()(const base_node &P,
				   dal::bit_vector &bv) const {
      scalar_type d = gmm::vect_dist2(P,x0)-R;
      bv[id] = (gmm::abs(d) < SEPS);
      return d;
    }
    virtual scalar_type operator()(const base_node &P) const
    { return gmm::vect_dist2(P,x0)-R; }
    virtual void register_constraints(std::vector<const
				      mesher_signed_distance*>& list) const {
      id = list.size(); list.push_back(this);
    }
    scalar_type grad(const base_node &P, base_small_vector &G) const {
      G = P; G -= x0;
      scalar_type e= gmm::vect_norm2(G), d = e - R;
      while (e == scalar_type(0))
	{ gmm::fill_random(G); e = gmm::vect_norm2(G); }
      G /= e;
      return d;
    }
    void hess(const base_node &, base_matrix &) const {
      GMM_ASSERT1(false, "Sorry, to be done");
    }
  };

  class mesher_rectangle : public mesher_signed_distance {
    // ajouter une rotation rigide et translation ..
    base_node rmin, rmax;
    std::vector<mesher_half_space> hfs;
  public:
    mesher_rectangle(base_node rmin_, base_node rmax_)
      : rmin(rmin_), rmax(rmax_) {
      base_node n(rmin_.size());
      for (unsigned k = 0; k < rmin.size(); ++k) {
	n[k] = 1.0;
	hfs.push_back(mesher_half_space(rmin, n));
	n[k] = -1.0;
	hfs.push_back(mesher_half_space(rmax, n));
	n[k] = 0.0;
      }
    }
    bool bounding_box(base_node &bmin, base_node &bmax) const {
      bmin = rmin; bmax = rmax;
      return true;
    }
    virtual scalar_type operator()(const base_node &P) const {
      size_type N = rmin.size();
      scalar_type d = rmin[0] - P[0];
      for (size_type i=0; i < N; ++i) {
	d = std::max(d, rmin[i] - P[i]);
	d = std::max(d, P[i] - rmax[i]);
      }
      return d;
    }

    virtual scalar_type operator()(const base_node &P, dal::bit_vector &bv)
      const {
      scalar_type d = this->operator()(P);
      if (gmm::abs(d) < SEPS)
	for (int k = 0; k < 2*rmin.size(); ++k) hfs[k](P, bv);
      return d;
    }
    scalar_type grad(const base_node &P, base_small_vector &G) const {
      unsigned i = 0; scalar_type di = hfs[i](P);
      for (int k = 1; k < 2*rmin.size(); ++k) {
	scalar_type dk = hfs[k](P);
	if (dk > di) { i = k; di = dk; }
      }
      return hfs[i].grad(P, G);
    }
    void hess(const base_node &P, base_matrix &H) const {
      gmm::resize(H, P.size(), P.size()); gmm::clear(H);
    }
    virtual void register_constraints(std::vector<const
				      mesher_signed_distance*>& list) const {
      for (int k = 0; k < 2*rmin.size(); ++k)
	hfs[k].register_constraints(list);
    }
  };


  class mesher_simplex_ref : public mesher_signed_distance {
    // ajouter une rotation rigide, homothetie,  et translation ..
    std::vector<mesher_half_space> hfs;
    unsigned N;
    base_node org;
  public:
    mesher_simplex_ref(unsigned N_) : N(N_) {
      base_node no(N);
      org = no;
      for (unsigned k = 0; k < N; ++k) {
	no[k] = scalar_type(1);
	hfs.push_back(mesher_half_space(org, no));
	no[k] = scalar_type(0);
      }
      std::fill(org.begin(), org.end(), scalar_type(1)/scalar_type(N));
      no = -org;
      hfs.push_back(mesher_half_space(org, no));
    }
    bool bounding_box(base_node &bmin, base_node &bmax) const {
      bmin.resize(N); bmax.resize(N);
      std::fill(bmin.begin(), bmin.end(), scalar_type(0));
      std::fill(bmax.begin(), bmax.end(), scalar_type(1));
      return true;
    }
    virtual scalar_type operator()(const base_node &P) const {
      scalar_type d = - P[0];
      for (size_type i=1; i < N; ++i) d = std::max(d, - P[i]);
      d = std::max(d, gmm::vect_sp(P - org, org) / gmm::vect_norm2(org));
      return d;
    }

    virtual scalar_type operator()(const base_node &P, dal::bit_vector &bv)
      const {
      scalar_type d = this->operator()(P);
      if (gmm::abs(d) < SEPS) for (unsigned k = 0; k < N+1; ++k) hfs[k](P, bv);
      return d;
    }
    scalar_type grad(const base_node &P, base_small_vector &G) const {
      unsigned i = 0; scalar_type di = hfs[i](P);
      for (unsigned k = 1; k < N+1; ++k) {
	scalar_type dk = hfs[k](P);
	if (dk > di) { i = k; di = dk; }
      }
      return hfs[i].grad(P, G);
    }
    void hess(const base_node &P, base_matrix &H) const {
      gmm::resize(H, P.size(), P.size()); gmm::clear(H);
    }
    virtual void register_constraints(std::vector<const
				      mesher_signed_distance*>& list) const
    { for (unsigned k = 0; k < N+1; ++k) hfs[k].register_constraints(list); }
  };

  class mesher_prism_ref : public mesher_signed_distance {
    // ajouter une rotation rigide, homothetie,  et translation ..
    std::vector<mesher_half_space> hfs;
    unsigned N;
    base_node org;
  public:
    mesher_prism_ref(unsigned N_) : N(N_) {
      base_node no(N);
      org = no;
      for (unsigned k = 0; k < N; ++k) {
	no[k] = scalar_type(1);
	hfs.push_back(mesher_half_space(org, no));
	no[k] = scalar_type(0);
      }
      no[N-1] = -scalar_type(1);
      org[N-1] = scalar_type(1);
      hfs.push_back(mesher_half_space(org, no));
      std::fill(org.begin(), org.end(), scalar_type(1)/scalar_type(N));
      org[N-1] = scalar_type(0);
      no = -org;
      hfs.push_back(mesher_half_space(org, no));
      
    }
    bool bounding_box(base_node &bmin, base_node &bmax) const {
      bmin.resize(N); bmax.resize(N);
      std::fill(bmin.begin(), bmin.end(), scalar_type(0));
      std::fill(bmax.begin(), bmax.end(), scalar_type(1));
      return true;
    }
    virtual scalar_type operator()(const base_node &P) const {
      scalar_type d = - P[0];
      for (size_type i=1; i < N; ++i) d = std::max(d, - P[i]);
      d = std::max(d, P[N-1] - scalar_type(1));
      d = std::max(d, gmm::vect_sp(P - org, org) / gmm::vect_norm2(org));
      return d;
    }

    virtual scalar_type operator()(const base_node &P, dal::bit_vector &bv)
      const {
      scalar_type d = this->operator()(P);
      if (gmm::abs(d) < SEPS) for (unsigned k = 0; k < N+2; ++k) hfs[k](P, bv);
      return d;
    }
    scalar_type grad(const base_node &P, base_small_vector &G) const {
      unsigned i = 0; scalar_type di = hfs[i](P);
      for (unsigned k = 1; k < N+2; ++k) {
	scalar_type dk = hfs[k](P);
	if (dk > di) { i = k; di = dk; }
      }
      return hfs[i].grad(P, G);
    }
    void hess(const base_node &P, base_matrix &H) const {
      gmm::resize(H, P.size(), P.size()); gmm::clear(H);
    }
    virtual void register_constraints(std::vector<const
				      mesher_signed_distance*>& list) const
    { for (unsigned k = 0; k < N+2; ++k) hfs[k].register_constraints(list); }
  };


  // Rem : It would be better for the convergence of Newton's methods to
  //       take somthing more smooth than min. for instance tthe square root
  //       of positive parts or negative parts depending on the position of
  //       point (in or out the domain).

  extern const mesher_half_space void_signed_distance;

  class mesher_union : public mesher_signed_distance {
    std::vector<const mesher_signed_distance *> dists;
    mutable std::vector<scalar_type> vd;
    mutable bool isin;
    bool with_min;
  public:
    mesher_union(const std::vector<const mesher_signed_distance *>
			&dists_) : dists(dists_) 
    { vd.resize(dists.size()); }

    mesher_union(const mesher_signed_distance &a_,
		 const mesher_signed_distance &b_,
		 const mesher_signed_distance &c_ = void_signed_distance,
		 const mesher_signed_distance &d_ = void_signed_distance,
		 const mesher_signed_distance &e_ = void_signed_distance,
		 const mesher_signed_distance &f_ = void_signed_distance,
		 const mesher_signed_distance &g_ = void_signed_distance,
		 const mesher_signed_distance &h_ = void_signed_distance,
		 const mesher_signed_distance &i_ = void_signed_distance,
		 const mesher_signed_distance &j_ = void_signed_distance,
		 const mesher_signed_distance &k_ = void_signed_distance,
		 const mesher_signed_distance &l_ = void_signed_distance,
		 const mesher_signed_distance &m_ = void_signed_distance,
		 const mesher_signed_distance &n_ = void_signed_distance,
		 const mesher_signed_distance &o_ = void_signed_distance,
		 const mesher_signed_distance &p_ = void_signed_distance,
		 const mesher_signed_distance &q_ = void_signed_distance,
		 const mesher_signed_distance &r_ = void_signed_distance,
		 const mesher_signed_distance &s_ = void_signed_distance,
		 const mesher_signed_distance &t_ = void_signed_distance) {
      dists.push_back(&a_); dists.push_back(&b_);
      size_type nb = 2; with_min = true;
      if (&c_ != &void_signed_distance) { dists.push_back(&c_); ++nb; }
      if (&d_ != &void_signed_distance) { dists.push_back(&d_); ++nb; }
      if (&e_ != &void_signed_distance) { dists.push_back(&e_); ++nb; }
      if (&f_ != &void_signed_distance) { dists.push_back(&f_); ++nb; }
      if (&g_ != &void_signed_distance) { dists.push_back(&g_); ++nb; }
      if (&h_ != &void_signed_distance) { dists.push_back(&h_); ++nb; }
      if (&i_ != &void_signed_distance) { dists.push_back(&i_); ++nb; }
      if (&j_ != &void_signed_distance) { dists.push_back(&j_); ++nb; }
      if (&k_ != &void_signed_distance) { dists.push_back(&k_); ++nb; }
      if (&l_ != &void_signed_distance) { dists.push_back(&l_); ++nb; }
      if (&m_ != &void_signed_distance) { dists.push_back(&m_); ++nb; }
      if (&n_ != &void_signed_distance) { dists.push_back(&n_); ++nb; }
      if (&o_ != &void_signed_distance) { dists.push_back(&o_); ++nb; }
      if (&p_ != &void_signed_distance) { dists.push_back(&p_); ++nb; }
      if (&q_ != &void_signed_distance) { dists.push_back(&q_); ++nb; }
      if (&r_ != &void_signed_distance) { dists.push_back(&r_); ++nb; }
      if (&s_ != &void_signed_distance) { dists.push_back(&s_); ++nb; }
      if (&t_ != &void_signed_distance) { dists.push_back(&t_); ++nb; }
      vd.resize(nb);
    }
    
    bool bounding_box(base_node &bmin, base_node &bmax) const {
      base_node bmin2, bmax2;
      bool b = dists[0]->bounding_box(bmin, bmax);
      if (!b) return false;
      for (size_type k = 1; k < dists.size(); ++k) {
	b = dists[k]->bounding_box(bmin2, bmax2);
	if (!b) return false;
	for (unsigned i=0; i < bmin.size(); ++i) { 
	  bmin[i] = std::min(bmin[i],bmin2[i]);
	  bmax[i] = std::max(bmax[i],bmax2[i]);
	}
      }
      return true;
    }
    virtual scalar_type operator()(const base_node &P) const {
      scalar_type d, f(0), g(1);
      if (with_min) {
	d = (*(dists[0]))(P);
	for (size_type k = 1; k < dists.size(); ++k)
	  d = std::min(d, (*(dists[k]))(P));
      }
      else { // essai raté ...
	isin = false;
	for (size_type k = 0; k < dists.size(); ++k) {
	  vd[k] = (*(dists[k]))(P);
	  if (vd[k] <= scalar_type(0)) isin = true;
	  f += gmm::sqr(gmm::neg(vd[k]));
	  g *= gmm::pos(vd[k]);
	}
	d = isin ? -gmm::sqrt(f)
	  : pow(g, scalar_type(1) / scalar_type(dists.size()));
      }
      return d;
    }
    scalar_type operator()(const base_node &P, dal::bit_vector &bv) const {
      if (with_min) {
	scalar_type d = vd[0] = (*(dists[0]))(P);
	bool ok = (d > -SEPS);
	for (size_type k = 1; k < dists.size(); ++k) {
	  vd[k] = (*(dists[k]))(P); if (vd[k] <= -SEPS) ok = false;
	  d = std::min(d,vd[k]);
	}
	for (size_type k = 0; ok && k < dists.size(); ++k) {
	  if (vd[k] < SEPS) (*(dists[k]))(P, bv);
	}
	return d;
      }
      else { // essai raté ...
	vd[0] = (*(dists[0]))(P);
	bool ok = (vd[0] > -SEPS);
	for (size_type k = 1; k < dists.size(); ++k) {
	  vd[k] = (*(dists[k]))(P); if (vd[k] <= -SEPS) ok = false;
	}
	for (size_type k = 0; ok && k < dists.size(); ++k) {
	  if (vd[k] < SEPS) (*(dists[k]))(P, bv);
	}
	return operator()(P);
      }
    }
    virtual void register_constraints(std::vector<const
				      mesher_signed_distance*>& list) const {
      for (size_type k = 0; k < dists.size(); ++k)
	dists[k]->register_constraints(list); 
    }
    scalar_type grad(const base_node &P, base_small_vector &G) const {
      scalar_type d;
      if (with_min /* || gmm::abs(d) < SEPS*/) {
	d = (*(dists[0]))(P);
	size_type i = 0;
	for (size_type k = 1; k < dists.size(); ++k) {
	  scalar_type d2 = (*(dists[k]))(P);
	  if (d2 < d) { d = d2; i = k; }
	}
	return dists[i]->grad(P, G);
      }
      else { // essai raté ...
	d = operator()(P);
	base_small_vector Gloc;
	for (size_type k = 0; k < dists.size(); ++k) {
	  dists[k]->grad(P, Gloc);
	  if (isin) 
	    Gloc *= -gmm::neg(vd[k]);
	  else 
	    Gloc *= pow(d, scalar_type(dists.size())) / vd[k];
	  if (!k) G = Gloc; else G += Gloc;
	}
	if (isin)
	  G *= scalar_type(1) / d;
	else
	  G /= pow(d, scalar_type(dists.size()-1)) * scalar_type(dists.size());
	return d;
      }
    }
    void hess(const base_node &P, base_matrix &H) const {
      scalar_type d = (*(dists[0]))(P);
      if (with_min || gmm::abs(d) < SEPS) {
	size_type i = 0;
	for (size_type k = 1; k < dists.size(); ++k) {
	  scalar_type d2 = (*(dists[k]))(P);
	  if (d2 < d) { d = d2; i = k; }
	}
	dists[i]->hess(P, H);
      }
      else {
	GMM_ASSERT1(false, "Sorry, to e done");
      }
    }
  };

  class mesher_intersection : public mesher_signed_distance {
    std::vector<const mesher_signed_distance *> dists;
    mutable std::vector<scalar_type> vd;

    // const mesher_signed_distance &a, &b;
  public:
    
    mesher_intersection(const std::vector<const mesher_signed_distance *>
			&dists_) : dists(dists_) 
    { vd.resize(dists.size()); }
    
    mesher_intersection(const mesher_signed_distance &a_,
		 const mesher_signed_distance &b_,
		 const mesher_signed_distance &c_ = void_signed_distance,
		 const mesher_signed_distance &d_ = void_signed_distance,
		 const mesher_signed_distance &e_ = void_signed_distance,
		 const mesher_signed_distance &f_ = void_signed_distance,
		 const mesher_signed_distance &g_ = void_signed_distance,
		 const mesher_signed_distance &h_ = void_signed_distance,
		 const mesher_signed_distance &i_ = void_signed_distance,
		 const mesher_signed_distance &j_ = void_signed_distance,
		 const mesher_signed_distance &k_ = void_signed_distance,
		 const mesher_signed_distance &l_ = void_signed_distance,
		 const mesher_signed_distance &m_ = void_signed_distance,
		 const mesher_signed_distance &n_ = void_signed_distance,
		 const mesher_signed_distance &o_ = void_signed_distance,
		 const mesher_signed_distance &p_ = void_signed_distance,
		 const mesher_signed_distance &q_ = void_signed_distance,
		 const mesher_signed_distance &r_ = void_signed_distance,
		 const mesher_signed_distance &s_ = void_signed_distance,
		 const mesher_signed_distance &t_ = void_signed_distance) {
      dists.push_back(&a_); dists.push_back(&b_);
      size_type nb = 2;
      if (&c_ != &void_signed_distance) { dists.push_back(&c_); ++nb; }
      if (&d_ != &void_signed_distance) { dists.push_back(&d_); ++nb; }
      if (&e_ != &void_signed_distance) { dists.push_back(&e_); ++nb; }
      if (&f_ != &void_signed_distance) { dists.push_back(&f_); ++nb; }
      if (&g_ != &void_signed_distance) { dists.push_back(&g_); ++nb; }
      if (&h_ != &void_signed_distance) { dists.push_back(&h_); ++nb; }
      if (&i_ != &void_signed_distance) { dists.push_back(&i_); ++nb; }
      if (&j_ != &void_signed_distance) { dists.push_back(&j_); ++nb; }
      if (&k_ != &void_signed_distance) { dists.push_back(&k_); ++nb; }
      if (&l_ != &void_signed_distance) { dists.push_back(&l_); ++nb; }
      if (&m_ != &void_signed_distance) { dists.push_back(&m_); ++nb; }
      if (&n_ != &void_signed_distance) { dists.push_back(&n_); ++nb; }
      if (&o_ != &void_signed_distance) { dists.push_back(&o_); ++nb; }
      if (&p_ != &void_signed_distance) { dists.push_back(&p_); ++nb; }
      if (&q_ != &void_signed_distance) { dists.push_back(&q_); ++nb; }
      if (&r_ != &void_signed_distance) { dists.push_back(&r_); ++nb; }
      if (&s_ != &void_signed_distance) { dists.push_back(&s_); ++nb; }
      if (&t_ != &void_signed_distance) { dists.push_back(&t_); ++nb; }
      vd.resize(nb);
    }
    bool bounding_box(base_node &bmin, base_node &bmax) const {
      base_node bmin2, bmax2;
      bool first;
      bool b = dists[0]->bounding_box(bmin, bmax); first = !b;
      for (size_type k = 1; k < dists.size(); ++k) {
	bool bb = dists[k]->bounding_box(bmin2, bmax2);
	for (unsigned i=0; i < bmin.size() && bb && !first; ++i) { 
	  bmin[i] = std::max(bmin[i],bmin2[i]);
	  bmax[i] = std::max(std::min(bmax[i],bmax2[i]), bmin[i]);
	}
	if (first && bb) { bmin = bmin2; bmax = bmax2; first = false; }
	b = b || bb;
      }
      return b;
    }
    virtual scalar_type operator()(const base_node &P) const {
      scalar_type d = (*(dists[0]))(P);
      for (size_type k = 1; k < dists.size(); ++k)
	d = std::max(d, (*(dists[k]))(P));
      return d;

    }
    scalar_type operator()(const base_node &P, dal::bit_vector &bv) const {
      scalar_type d = vd[0] = (*(dists[0]))(P);
      bool ok = (d < SEPS);
      for (size_type k = 1; k < dists.size(); ++k) {
	vd[k] = (*(dists[k]))(P); if (vd[k] >= SEPS) ok = false;
	d = std::min(d,vd[k]);
      }
      for (size_type k = 0; ok && k < dists.size(); ++k) {
	if (vd[k] > -SEPS) (*(dists[k]))(P, bv);
      }
      return d;
    }
    virtual void register_constraints(std::vector<const
				      mesher_signed_distance*>& list) const {
      for (size_type k = 0; k < dists.size(); ++k) {
	dists[k]->register_constraints(list); 
      }
    }
    scalar_type grad(const base_node &P, base_small_vector &G) const {
      scalar_type d = (*(dists[0]))(P);
      size_type i = 0;
      for (size_type k = 1; k < dists.size(); ++k) {
	scalar_type d2 = (*(dists[k]))(P);
	if (d2 > d) { d = d2; i = k; }
      }
      return dists[i]->grad(P, G);
    }
    void hess(const base_node &P, base_matrix &H) const {
      scalar_type d = (*(dists[0]))(P);
      size_type i = 0;
      for (size_type k = 1; k < dists.size(); ++k) {
	scalar_type d2 = (*(dists[k]))(P);
	if (d2 > d) { d = d2; i = k; }
      }
      dists[i]->hess(P, H);
    }
  };

  class mesher_setminus : public mesher_signed_distance {
    const mesher_signed_distance &a, &b;
  public:
    mesher_setminus(const mesher_signed_distance& a_,
		    const mesher_signed_distance &b_) : a(a_), b(b_) {}
    bool bounding_box(base_node &bmin, base_node &bmax) const
    { return a.bounding_box(bmin,bmax); }
    scalar_type operator()(const base_node &P, dal::bit_vector &bv) const {
      scalar_type da = a(P), db = -b(P);
      if (da < SEPS && db < SEPS) {
	if (da > -SEPS) a(P, bv);
	if (db > -SEPS) b(P, bv);
      }
      return std::max(da, db);
    }
    scalar_type operator()(const base_node &P) const
    { return std::max(a(P),-b(P)); }
    virtual void register_constraints(std::vector<const
				      mesher_signed_distance*>& list) const {
      a.register_constraints(list); b.register_constraints(list);
    }
    scalar_type grad(const base_node &P, base_small_vector &G) const {
      scalar_type da = a(P), db = -b(P);
      if (da > db) return a.grad(P, G);
      else { b.grad(P, G); G *= scalar_type(-1); return db; }
    }
    void hess(const base_node &P, base_matrix &H) const {
      scalar_type da = a(P), db = -b(P);
      if (da > db) a.hess(P, H);
      else { b.hess(P, H); gmm::scale(H, scalar_type(-1)); }
    }

  };
  

  class mesher_tube : public mesher_signed_distance {
    base_node x0; base_node n; scalar_type R;
  public:
      mesher_tube(base_node x0_, base_node n_, scalar_type R_)
	: x0(x0_), n(n_), R(R_)
    { n /= gmm::vect_norm2(n); }
    bool bounding_box(base_node &, base_node &) const
    { return false; }
    virtual scalar_type operator()(const base_node &P) const {
      base_node v(P); v -= x0;
      gmm::add(gmm::scaled(n, -gmm::vect_sp(v, n)), v);
      return gmm::vect_norm2(v) - R;
    }
    virtual scalar_type operator()(const base_node &P,
				   dal::bit_vector &bv) const {
      scalar_type d = (*this)(P);
      bv[id] = (gmm::abs(d) < SEPS);
      return d;
    }
    virtual void register_constraints(std::vector<const
				      mesher_signed_distance*>& list) const {
      id = list.size(); list.push_back(this);
    }
    scalar_type grad(const base_node &P, base_small_vector &G) const {
      G = P; G -= x0;
      gmm::add(gmm::scaled(n, -gmm::vect_sp(G, n)), G);
      scalar_type e = gmm::vect_norm2(G), d = e - R;
      while (e == scalar_type(0)) {
	gmm::fill_random(G);
	gmm::add(gmm::scaled(n, -gmm::vect_sp(G, n)), G);
	e = gmm::vect_norm2(G);
      }
      G /= e;
      return d;
    }
    void hess(const base_node &, base_matrix &) const {
      GMM_ASSERT1(false, "Sorry, to be done");
    }
  };

  class mesher_cylinder : public mesher_signed_distance {
    base_node x0; base_small_vector n;
    scalar_type L, R;
    mesher_tube t;
    mesher_half_space p1, p2;
    mesher_intersection i1;
  public:
    mesher_cylinder(const base_node &c, const base_small_vector &no,
		    scalar_type L_, scalar_type R_)
      : x0(c), n(no/gmm::vect_norm2(no)), L(L_), R(R_), t(x0, n, R),
	p1(x0, n), p2(x0+n*L, -1.0 * n), i1(p1, p2, t) {}
    bool bounding_box(base_node &bmin, base_node &bmax) const {
      base_node x1(x0+n*L);
      bmin = bmax = x0;
      for (unsigned i = 0; i < gmm::vect_size(x0); ++i) {
	bmin[i] = std::min(x0[i], x1[i]) - R;
	bmax[i] = std::max(x0[i], x1[i]) + R;
      }
      return true;
    }
    virtual scalar_type operator()(const base_node &P) const { return i1(P); }
    virtual scalar_type operator()(const base_node &P,
				   dal::bit_vector& bv) const
    { return i1(P, bv); }
    scalar_type grad(const base_node &P, base_small_vector &G) const
      { return i1.grad(P, G); }
    void hess(const base_node &, base_matrix &) const {
      GMM_ASSERT1(false, "Sorry, to be done");
    }
    virtual void register_constraints(std::vector<const
				      mesher_signed_distance*>& list) const
    { i1.register_constraints(list); }
  };



  class mesher_infinite_cone : public mesher_signed_distance {
    // uses the true distance to a cone.
    base_node x0; base_node n; scalar_type alpha;
  public:
      mesher_infinite_cone(base_node x0_, base_node n_, scalar_type alpha_)
	: x0(x0_), n(n_), alpha(alpha_)
    { n /= gmm::vect_norm2(n); }
    bool bounding_box(base_node &, base_node &) const
    { return false; }
    virtual scalar_type operator()(const base_node &P) const {
      base_node v(P); v -= x0;
      scalar_type v_n = gmm::vect_sp(v, n);
      gmm::add(gmm::scaled(n, -v_n), v);
      return gmm::vect_norm2(v) * cos(alpha) - gmm::abs(v_n) * sin(alpha);
    }
    virtual scalar_type operator()(const base_node &P,
				   dal::bit_vector &bv) const {
      scalar_type d = (*this)(P);
      bv[id] = (gmm::abs(d) < SEPS);
      return d;
    }
    virtual void register_constraints(std::vector<const
				      mesher_signed_distance*>& list) const {
      id = list.size(); list.push_back(this);
    }
    scalar_type grad(const base_node &P, base_small_vector &v) const {
      v = P; v -= x0;
      scalar_type v_n = gmm::vect_sp(v, n);
      gmm::add(gmm::scaled(n, -v_n), v);
      scalar_type no = gmm::vect_norm2(v);
      scalar_type d = no * cos(alpha) - gmm::abs(v_n) * sin(alpha);
      while (no ==  scalar_type(0)) {
	gmm::fill_random(v);
	gmm::add(gmm::scaled(n, -gmm::vect_sp(v, n)), v);
	no = gmm::vect_norm2(v);
      }
      v *= cos(alpha) / no;
      v -= (sin(alpha) * gmm::sgn(v_n)) * n;
      return d;
    }
    void hess(const base_node &, base_matrix &) const {
      GMM_ASSERT1(false, "Sorry, to be done");
    }
  };

  class mesher_cone : public mesher_signed_distance {
    base_node x0; base_small_vector n;
    scalar_type L, alpha;
    mesher_infinite_cone t;
    mesher_half_space p1, p2;
    mesher_intersection i1;
  public:
    mesher_cone(const base_node &c, const base_small_vector &no,
		    scalar_type L_, scalar_type alpha_)
      : x0(c), n(no/gmm::vect_norm2(no)), L(L_), alpha(alpha_),
	t(x0, n, alpha), p1(x0, n), p2(x0+n*L, -1.0 * n), i1(p1, p2, t) {}
    bool bounding_box(base_node &bmin, base_node &bmax) const {
      base_node x1(x0+n*L);
      scalar_type R = L * sin(alpha);
      bmin = bmax = x0;
      for (unsigned i = 0; i < gmm::vect_size(x0); ++i) {
	bmin[i] = std::min(x0[i], x1[i]) - R;
	bmax[i] = std::max(x0[i], x1[i]) + R;
      }
      return true;
    }
    virtual scalar_type operator()(const base_node &P) const { return i1(P); }
    virtual scalar_type operator()(const base_node &P,
				   dal::bit_vector& bv) const
    { return i1(P, bv); }
    scalar_type grad(const base_node &P, base_small_vector &G) const
      { return i1.grad(P, G); }
    void hess(const base_node &, base_matrix &) const {
      GMM_ASSERT1(false, "Sorry, to be done");
    }
    virtual void register_constraints(std::vector<const
				      mesher_signed_distance*>& list) const
    { i1.register_constraints(list); }
  };

  class mesher_ellipse : public mesher_signed_distance {
    base_node x0; base_small_vector n, t;
    scalar_type r, R, a;
  public:
    mesher_ellipse(const base_node &center, const base_small_vector &no,
		    scalar_type r_, scalar_type R_)
      : x0(center), n(no/gmm::vect_norm2(no)), r(r_), R(R_) {
      t[0] = -n[1]; t[1] = n[0];
      if (R < r) { std::swap(r, R); std::swap(n, t); }
      a = sqrt(R*R - r*r);
    }
    bool bounding_box(base_node &bmin, base_node &bmax) const {
      bmin = bmax = x0;
      for (unsigned i = 0; i < 2; ++i) { bmin[i] -= R; bmax[i] += R; }
      return true;
    }
    virtual scalar_type operator()(const base_node &P) const { 
      base_small_vector v(P); v -= x0;
      scalar_type vt = gmm::vect_sp(v, t);
      vt = std::max(-a, std::min(a, vt));
      base_node x1 = x0 + t*vt;
      base_small_vector v1(P); v1 -= x1;
      scalar_type v1n = gmm::vect_sp(v1, n), v1t = gmm::vect_sp(v1, t);
      scalar_type x1n = gmm::vect_sp(x1, n), x1t = gmm::vect_sp(x1, t);
      scalar_type ea = v1n*v1n / (r*r) + v1t * v1t / (R*R);
      scalar_type eb = 2. * (x1n*v1n / (r*r) + x1t*v1t / (R*R));
      scalar_type ec = x1n*x1n / (r*r) + x1t * x1t / (R*R);

      scalar_type delta = eb*eb - 4 * ea * ec;
      assert(delta >= 0);
      scalar_type lambda = (-eb + sqrt(delta)) / (2. * ea);
      return (1.-lambda)*gmm::vect_norm2(v1);
    }
    virtual scalar_type operator()(const base_node &P,
				   dal::bit_vector& bv) const {
      scalar_type d = this->operator()(P);
      bv[id] = (gmm::abs(d) < SEPS);
      return d;
    }
    scalar_type grad(const base_node &, base_small_vector &) const
    { GMM_ASSERT1(false, "Sorry, to be done"); return 0.; }
    void hess(const base_node &, base_matrix &) const {
      GMM_ASSERT1(false, "Sorry, to be done");
    }
    virtual void register_constraints(std::vector<const
				      mesher_signed_distance*>& list) const
    { id = list.size(); list.push_back(this); }
  };


  class mesher_torus : public mesher_signed_distance { // to be done
    // ajouter une rotation rigide et translation ..
    scalar_type R, r;
  public:
    mesher_torus(scalar_type RR, scalar_type rr) : R(RR), r(rr) {}
    bool bounding_box(base_node &bmin, base_node &bmax) const {
      bmin = base_node(3); bmax = base_node(3);
      bmin[0] = bmin[1] = -R-r; bmin[2] = -r;
      bmax[0] = bmax[1] = +R+r; bmax[2] = +r;
      return true;
    }
    virtual scalar_type operator()(const base_node &P) const {
      scalar_type x = P[0], y = P[1], z = P[2], c = sqrt(x*x + y*y);
      return (c == 0.) ? R - r : sqrt(gmm::sqr(c-R) + z*z) - r;
    }
    virtual scalar_type operator()(const base_node &P, dal::bit_vector&bv)
      const {
      scalar_type d = this->operator()(P);
      bv[id] = (gmm::abs(d) < SEPS);
      return d;
    }
    scalar_type grad(const base_node &P, base_small_vector &G) const {
      G.resize(3);
      scalar_type x = P[0], y = P[1], z = P[2], c = sqrt(x*x + y*y), d(0);
      if (c == 0.) {
	d = R - r;
	gmm::fill_random(G); G[2] = 0.0; G /= gmm::vect_norm2(G);
      }
      else {
	scalar_type w = 1. - R / c, e = sqrt(gmm::sqr(c-R) + z*z);
	d = e - r;
	if (e == 0.) {
	  gmm::fill_random(G); G[0] = x; G[1] = y; G /= gmm::vect_norm2(G);
	}
	else {
	  G[0] = x * w / e; G[1] = y * w / e; G[2] = z / e;
	}
      }
      return d;
    }
    void hess(const base_node &, base_matrix &) const {
      GMM_ASSERT1(false, "Sorry, to be done");
    }
    virtual void register_constraints(std::vector<const
				      mesher_signed_distance*>&list) const
    { id = list.size(); list.push_back(this); }
  };


  

  // interface with qhull
  void delaunay(const std::vector<base_node> &pts,
		gmm::dense_matrix<size_type>& simplexes);

  
  // mesher
  void build_mesh(mesh &m, const mesher_signed_distance& dist_,
		  scalar_type h0, const std::vector<base_node> &fixed_points
		  = std::vector<base_node>(), size_type K = 1, int noise = -1,
		  size_type iter_max = 500, int prefind = 1,
		  scalar_type dist_point_hull = 4,
		  scalar_type boundary_threshold_flatness = 0.11);

  // exported functions
  bool try_projection(const mesher_signed_distance& dist, base_node &X,
		      bool on_surface = false);
  bool pure_multi_constraint_projection
  (const std::vector<const mesher_signed_distance*> &list_constraints,
   base_node &X, const dal::bit_vector &cts);
  scalar_type curvature_radius_estimate(const mesher_signed_distance &dist,
					base_node X, bool proj = false);
  scalar_type min_curvature_radius_estimate
  (const std::vector<const mesher_signed_distance*> &list_constraints,
   const base_node &X, const dal::bit_vector &cts, size_type hide_first = 0);

}

#endif /* GETFEM_MESHER_H__ */