rsvs3D  0.0.0
Codes for the c++ implementation of the 3D RSVS
RSVSmath.cpp
1 #include "RSVSmath.hpp"
2 
3 #include <iostream>
4 
5 #include "RSVSmath_automatic.hpp"
6 #include "arraystructures.hpp" // for use of DisplayVector
7 #include "parameters.hpp"
8 #include "warning.hpp"
9 
10 using namespace std;
11 
12 bool TriFunc::MakeValidField(vector<double> *TriFunc::*mp)
13 {
14  // Tries to make a valid array
15  // Warning : This operates directly on the data stored in the field of the
16  // pointer
17 
18  int ii, n;
19  bool fieldReady = true;
20  if ((this->*mp) == NULL)
21  {
22  fieldReady = false;
23  }
24  else
25  {
26  n = int((this->*mp)->size());
27  if (n < nTarg)
28  {
29  for (ii = 0; ii < (nTarg - n); ++ii)
30  {
31  (this->*mp)->push_back(0);
32  }
33  }
34  else if (n > nTarg)
35  {
36  for (ii = 0; ii < (n - nTarg); ++ii)
37  {
38  (this->*mp)->pop_back();
39  }
40  }
41  }
42  return (fieldReady);
43 }
44 
45 bool TriFunc::CheckValid()
46 {
47  isReady = true;
48  if (p0 == NULL)
49  {
50  isReady = false;
51  }
52  else if (int(p0->size()) != nTarg)
53  {
54  isReady = false;
55  }
56  if (p1 == NULL)
57  {
58  isReady = false;
59  }
60  else if (int(p1->size()) != nTarg)
61  {
62  isReady = false;
63  }
64  if (p2 == NULL)
65  {
66  isReady = false;
67  }
68  else if (int(p2->size()) != nTarg)
69  {
70  isReady = false;
71  }
72 
73  return (isReady);
74 }
75 
76 bool TriFunc::MakeValid()
77 {
78  // Tries to make a valid array
79  // Warning : This operates directly on the data of the pointer
80 
81  isReady = true;
82  // bool fieldReady;
83  CheckValid();
84  // fieldReady=MakeValidField(&TriFunc::p0);isReady=isReady & fieldReady;
85  // fieldReady=MakeValidField(&TriFunc::p1);isReady=isReady & fieldReady;
86  // fieldReady=MakeValidField(&TriFunc::p2);isReady=isReady & fieldReady;
87 
88  return (isReady);
89 }
90 
91 void TriFunc::assign(const vector<double> &in0, const vector<double> &in1, const vector<double> &in2)
92 {
93  p0 = &in0;
94  p1 = &in1;
95  p2 = &in2;
96 
97  isCalc = false;
98  isReady = false;
99 }
100 void TriFunc::assign(const vector<double> *in0, const vector<double> *in1, const vector<double> *in2)
101 {
102  p0 = in0;
103  p1 = in1;
104  p2 = in2;
105 
106  isCalc = false;
107  isReady = false;
108 }
109 void TriFunc::assign(int pRepI, const vector<double> &pRep)
110 {
111  switch (pRepI)
112  {
113  case 0:
114  p0 = &pRep;
115  break;
116  case 1:
117  p1 = &pRep;
118  break;
119  case 2:
120  p2 = &pRep;
121  break;
122  }
123 
124  isCalc = false;
125  isReady = false;
126 }
127 
128 void TriFunc::PreCalc()
129 {
130  if (!isReady)
131  {
132  CheckValid();
133  }
134  if (!isReady)
135  {
136  MakeValid();
137  }
138  if (!isReady)
139  {
140  cerr << "TriFunc cannot be made to match the required paremeters " << endl;
141  }
142 }
143 
144 void TriFunc::ReturnDatPoint(double **a, ArrayVec<double> **b, ArrayVec<double> **c)
145 {
146  *a = &fun;
147  *b = &jac;
148  *c = &hes;
149 }
150 
153 
154 bool CoordFunc::MakeValidField(vector<double> const *mp)
155 {
156  // Tries to make a valid array
157  // Warning : This operates directly on the data stored in the field of the
158  // pointer
159 
160  // int ii,n;
161  bool fieldReady;
162  if ((mp) == NULL)
163  {
164  fieldReady = false;
165  }
166  else
167  {
168  fieldReady = false;
169  // n=int((mp)->size());
170  // if (n<nDim){
171  // for(ii=0;ii<(nDim-n);++ii){
172  // (mp)->push_back(0);
173  // }
174  // } else if (n>nDim){
175  // for(ii=0;ii<(n-nDim);++ii){
176  // (mp)->pop_back();
177  // }
178  // }
179  }
180  return (fieldReady);
181 }
182 
183 bool CoordFunc::CheckValid()
184 {
185  isReady = true;
186 
187  for (int ii = 0; ii < nCoord; ++ii)
188  {
189  if (coords[ii] == NULL)
190  {
191  isReady = false;
192  }
193  else if (int(coords[ii]->size()) != nDim)
194  {
195  isReady = false;
196  }
197  }
198  return (isReady);
199 }
200 
201 bool CoordFunc::MakeValid()
202 {
203  // Tries to make a valid array
204  // Warning : This operates directly on the data of the pointer
205 
206  isReady = true;
207  bool fieldReady;
208  for (int ii = 0; ii < nCoord; ++ii)
209  {
210  fieldReady = MakeValidField(coords[ii]);
211  isReady = isReady & fieldReady;
212  }
213 
214  return (isReady);
215 }
216 
217 void CoordFunc::assign(grid::coordlist &pRep)
218 {
219  if (int(pRep.size()) == nCoord)
220  {
221  coords = pRep;
222  }
223  else
224  {
225  ResetNCoord(int(pRep.size()));
226  coords = pRep;
227  }
228  isCalc = false;
229  isReady = false;
230 }
231 
232 void CoordFunc::assign(int pRepI, const vector<double> &pRep)
233 {
234  if ((pRepI < nCoord) & (pRepI >= 0))
235  {
236  coords[pRepI] = &pRep;
237  }
238  else
239  {
240  RSVS3D_ERROR_ARGUMENT("pointer index out off bound in CoordFunc::assign");
241  }
242 
243  isCalc = false;
244  isReady = false;
245 }
246 
247 void CoordFunc::ReturnDat(double &a, ArrayVec<double> &b, ArrayVec<double> &c)
248 {
249  a = fun;
250  b = jac;
251  c = hes;
252 }
253 void CoordFunc::ReturnDat(ArrayVec<double> &a, ArrayVec<double> &b, ArrayVec<double> &c)
254 {
255  a = funA;
256  b = jac;
257  c = hes;
258 }
259 void CoordFunc::ReturnDatPoint(double **a, ArrayVec<double> **b, ArrayVec<double> **c)
260 {
261  *a = &fun;
262  *b = &jac;
263  *c = &hes;
264 }
265 void CoordFunc::ReturnDatPoint(ArrayVec<double> **a, ArrayVec<double> **b, ArrayVec<double> **c)
266 {
267  *a = &funA;
268  *b = &jac;
269  *c = &hes;
270 }
271 
272 void CoordFunc::PreCalc()
273 {
274  if (!isReady)
275  {
276  CheckValid();
277  }
278  if (!isReady)
279  {
280  MakeValid();
281  }
282  if (!isReady)
283  {
284  cerr << "CoordFunc cannot be made to match the required paremeters " << endl;
285  }
286 }
287 
288 void CoordFunc::InitialiseArrays()
289 {
290  fun = 0;
291  if (nFun > 1)
292  {
293  funA.assign(nFun, 1, fun);
294  }
295  jac.assign(nFun, nCoord * nDim, fun);
296  hes.assign(nCoord * nDim, nCoord * nDim * nFun, fun);
297  coords.assign(nCoord, NULL);
298  isReady = false;
299  isCalc = false;
300 }
301 
302 // void CoordFunc::ResetDim(int nDim){ // implicitely inlined
303 // InitialiseArrays();
304 // }
305 // void CoordFunc::ResetNCoord(int nCoord){
306 // InitialiseArrays();
307 // }
308 
309 // Derived classes
310 
311 void Volume::Calc()
312 {
313 /* This function calculates Volume Jacobian and Hessian using
314 matlab automatically generated code
315 
316 the jacobian is arranged :
317  x0 y0 z0 x1 y1 z1 x2 y2 z2
318  V [ ]
319 
320  The Hessian is arranged:
321  V
322  x0 y0 z0 x1 y1 z1 x2 y2 z2
323  x0 [ ]
324  y0 [ ]
325  z0 [ ]
326  x1 [ ]
327  y1 [ ]
328  z1 [ ]
329  x2 [ ]
330  y2 [ ]
331  z2 [ ]
332 
333 */
334 #ifdef SAFE_ALGO
335  PreCalc();
336 #endif
337 
338  if (!isCalc)
339  {
340  Volume_f(*p0, *p1, *p2, fun);
341  Volume_df(*p0, *p1, *p2, jac);
342  Volume_ddf(*p0, *p1, *p2, hes);
343  }
344 }
345 
346 void Volume::CalcFD()
347 {
348 /* This function calculates Volume Jacobian and Hessian using
349 matlab automatically generated code
350 
351 the jacobian is arranged :
352  x0 y0 z0 x1 y1 z1 x2 y2 z2
353  V [ ]
354 
355  The Hessian is arranged:
356  V
357  x0 y0 z0 x1 y1 z1 x2 y2 z2
358  x0 [ ]
359  y0 [ ]
360  z0 [ ]
361  x1 [ ]
362  y1 [ ]
363  z1 [ ]
364  x2 [ ]
365  y2 [ ]
366  z2 [ ]
367 
368 */
369 #ifdef SAFE_ALGO
370  PreCalc();
371 #endif
372  std::vector<double> v;
373  double fdStep = 1e-6;
374  double tempVal;
375 
376  vector<vector<double>> vecs;
377 
378  v.reserve(3);
379 
380  if (!isCalc)
381  {
382  this->Calc();
383 
384  vecs.clear();
385  vecs.reserve(3);
386  vecs.push_back(*p0);
387  vecs.push_back(*p1);
388  vecs.push_back(*p2);
389  for (int ii = 0; ii < 3; ++ii)
390  {
391  for (int jj = 0; jj < 3; ++jj)
392  {
393  vecs[ii][jj] = vecs[ii][jj] + fdStep;
394  Volume_f(vecs[0], vecs[1], vecs[2], tempVal);
395  jac[0][ii * 3 + jj] = (tempVal - this->fun) / fdStep;
396  vecs[ii][jj] = vecs[ii][jj] - fdStep;
397  }
398  }
399  }
400 }
401 
402 void Area::Calc()
403 {
404 /* This function calculates Area Jacobian and Hessian using
405 matlab automatically generated code
406 
407 the jacobian is arranged :
408  x0 y0 z0 x1 y1 z1 x2 y2 z2
409  A [ ]
410 
411  The Hessian is arranged:
412  A
413  x0 y0 z0 x1 y1 z1 x2 y2 z2
414  x0 [ ]
415  y0 [ ]
416  z0 [ ]
417  x1 [ ]
418  y1 [ ]
419  z1 [ ]
420  x2 [ ]
421  y2 [ ]
422  z2 [ ]
423 
424 */
425 #ifdef SAFE_ALGO
426  PreCalc();
427 #endif
428 
429  if (!isCalc)
430  {
431  Area_f(*p0, *p1, *p2, fun);
432  Area_df(*p0, *p1, *p2, jac);
433  Area_ddf(*p0, *p1, *p2, hes);
434  }
435 }
436 
437 void LengthEdge::Calc()
438 {
439 #ifdef SAFE_ALGO
440  PreCalc();
441 #endif
442  if (!isCalc)
443  {
444  LengthEdge_f(*coords[0], *coords[1], fun);
445  LengthEdge_df(*coords[0], *coords[1], jac);
446  LengthEdge_ddf(*coords[0], *coords[1], hes);
447  }
448 }
449 
450 // void SurfCentroid::Calcddf(){
451 
452 // }
453 
454 void SurfCentroid::CalcFD()
455 {
456  this->Calc();
457  auto baseVal = this->funA;
458  std::vector<double> vtemp;
459  double fdStep = 1e-6;
460 
461  this->calcDeriv = false;
462  this->jac.assign(3, nCoord * 3, 0.0);
463 
464  for (int i = 0; i < nCoord; ++i)
465  {
466  auto saveVec = this->coords[i];
467  this->coords[i] = &vtemp;
468  for (int j = 0; j < 3; ++j)
469  {
470  this->isCalc = false;
471  vtemp = *saveVec;
472  vtemp[j] = vtemp[j] + fdStep;
473  this->Calc();
474  for (int k = 0; k < 3; ++k)
475  {
476  this->jac[k][i + nCoord * j] = (this->funA[k][0] - baseVal[k][0]) / fdStep;
477  }
478  }
479  this->coords[i] = saveVec;
480  }
481 }
482 
483 void SurfCentroid::Calc()
484 {
485  /* This function calculates centroid Jacobian and Hessian using
486  matlab automatically generated code
487 
488  the jacobian is arranged :
489  x0 x1 xn y0 y1 yn z0 z1 zn
490  xc [ ]
491  yc [ ]
492  zc [ ]
493 
494  The Hessian is arranged:
495  xc
496  yc zc x0 x1 xn y0
497  y1 yn z0 z1 zn x0 x1 xn y0 y1 yn z0 z1 zn x0 x1 xn y0 y1 yn
498  z0 z1 zn x0 [ ] [ ] [ ] x1 [ ] [ ] [ ] xn [
499  ] [ ] [ ] y0 [
500  ] [ ] [ ] y1 [
501  ] [ ] [ ] yn [
502  ] [ ] [ ] z0 [
503  ] [ ] [ ] z1 [
504  ] [ ] [ ] zn [
505  ] [ ] [ ]
506 
507  */
508 
509  vector<double> x, y, z, centroidTLength;
510  ArrayVec<double> temp;
511  double currEdge;
512  int ii, ii1, ii2, jj, jj2, nR, nC, ind, ind2;
513 
514  this->PreCalc();
515 
516  if (!isCalc)
517  {
518  // Calculate centroid and
519  centroid.assign(nDim, 0);
520  edgeLength = 0;
521  for (jj = 0; jj < nCoord; ++jj)
522  {
523  currEdge = 0;
524  LengthEdge_f(*coords[jj], *coords[(jj + 1) % nCoord], currEdge);
525  edgeLength += currEdge;
526  for (ii = 0; ii < nDim; ++ii)
527  {
528  centroid[ii] += currEdge * ((*(coords[jj]))[ii] + (*(coords[(jj + 1) % nCoord]))[ii]) / 2;
529  }
530  }
531  if (edgeLength < 0.000000000000001)
532  {
533  cout << "Warning edgeLength is 0" << endl;
534  for (ii = 0; ii < nDim; ++ii)
535  {
536  funA[ii][0] = (*(coords[0]))[ii];
537  }
538  edgeLength = 0.0000001;
539  }
540  else
541  {
542  for (ii = 0; ii < nDim; ++ii)
543  {
544  funA[ii][0] = centroid[ii] / edgeLength;
545  }
546  }
547  if (!this->calcDeriv)
548  {
549  return;
550  }
551  if (nCoord < 4)
552  {
553  RSVS3D_ERROR_ARGUMENT("nCoord <= 3 surfCentroid Not implemented");
554  }
555  else if (nCoord < 7)
556  {
557  x.reserve(nCoord);
558  y.reserve(nCoord);
559  z.reserve(nCoord);
560  // Build Full x,y,z vectors
561  for (ii = 0; ii < nCoord; ++ii)
562  {
563  x.push_back((*coords[ii])[0]);
564  y.push_back((*coords[ii])[1]);
565  z.push_back((*coords[ii])[2]);
566  }
567 
568  switch (nCoord)
569  {
570  case 4:
571  // SurfCentroid4_f(x,y,z,edgeLength,centroid[0],centroid[1],centroid[2],funA);
572  SurfCentroid4_df(x, y, z, edgeLength, centroid[0], centroid[1], centroid[2], jac);
573  break;
574  case 5:
575  // SurfCentroid5_f(x,y,z,edgeLength,centroid[0],centroid[1],centroid[2],funA);
576  SurfCentroid5_df(x, y, z, edgeLength, centroid[0], centroid[1], centroid[2], jac);
577  break;
578  case 6:
579  // SurfCentroid6_f(x,y,z,edgeLength,centroid[0],centroid[1],centroid[2],funA);
580  SurfCentroid6_df(x, y, z, edgeLength, centroid[0], centroid[1], centroid[2], jac);
581  break;
582  }
583  }
584  else
585  {
586  /*the jacobian is arranged :
587  x0 x1 xn y0 y1 yn z0 z1 zn
588  xc [ ]
589  yc [ ]
590  zc [ ]*/
591  x.assign(3, 0);
592  y.assign(3, 0);
593  z.assign(3, 0);
594  temp.assign(nDim, nDim, 0);
595  for (jj = 0; jj < nCoord; jj++)
596  {
597  for (ii = 0; ii < 3; ++ii)
598  {
599  x[ii] = ((*coords[(jj + ii) % nCoord])[0]);
600  y[ii] = ((*coords[(jj + ii) % nCoord])[1]);
601  z[ii] = ((*coords[(jj + ii) % nCoord])[2]);
602  }
603  SurfCentroidSelf_df(x, y, z, edgeLength, centroid[0], centroid[1], centroid[2], temp);
604  ind = (jj + 1) % nCoord;
605  // WARNING : This needs to be checked to make sure it matches the
606  // automaticly generated ones
607  for (ii1 = 0; ii1 < nDim; ii1++)
608  {
609  for (ii2 = 0; ii2 < nDim; ii2++)
610  {
611  jac[ii1][ind + (ii2 * nCoord)] = temp[ii1][ii2];
612  }
613  }
614  }
615  }
616 
617  // Calculate Hessian
618 
619  /* The Hessian is arranged:
620  xc yc
621  zc x0 x1 xn y0 y1 yn z0 z1 zn x0 x1 xn y0 y1 yn z0 z1 zn x0
622  x1 xn y0 y1 yn z0 z1 zn x0 [ ] [ ] [ ] x1 [
623  ] [ ] [ ] xn [ ]
624  [ ] [ ] y0 [ ] [
625  ] [ ] y1 [ ] [ ]
626  [ ] yn [ ] [ ] [ ]
627  z0 [ ] [ ] [ ] z1
628  [ ] [ ] [ ] zn [
629  ] [ ] [ ] */
630 
631  // Squared terms of the Hessian
632  hes.size(nR, nC);
633  x.assign(3, 0);
634  y.assign(3, 0);
635  z.assign(3, 0);
636  temp.assign(nDim, nDim * nFun, 0);
637  for (jj = 0; jj < nCoord; jj++)
638  {
639  for (ii = 0; ii < 3; ++ii)
640  {
641  x[ii] = ((*coords[(jj + ii) % nCoord])[0]);
642  y[ii] = ((*coords[(jj + ii) % nCoord])[1]);
643  z[ii] = ((*coords[(jj + ii) % nCoord])[2]);
644  }
645  SurfCentroidSelf_ddf(x, y, z, edgeLength, centroid[0], centroid[1], centroid[2], temp);
646  ind = (jj + 1) % nCoord;
647  // WARNING : This needs to be checked to make sure it matches the
648  // automaticly generated ones
649  for (ii1 = 0; ii1 < nDim; ii1++)
650  {
651  for (ii2 = 0; ii2 < nDim; ii2++)
652  {
653  for (ii = 0; ii < nFun; ii++)
654  {
655  hes[ind + (ii1 * nCoord)][ind + (ii2 * nCoord) + (ii * nDim * nCoord)] =
656  temp[ii1][ii2 + (nDim * ii)];
657  }
658  }
659  }
660  }
661 
662  // Connected terms of the Hessian
663  x.assign(6, 0);
664  y.assign(6, 0);
665  z.assign(6, 0);
666  temp.assign(nDim * 2, 2 * nDim * nFun, 0);
667  for (jj = 0; jj < nCoord; jj++)
668  {
669  jj2 = (jj + 1) % nCoord;
670  for (ii = 0; ii < 4; ++ii)
671  {
672  x[ii] = ((*coords[(jj + ii) % nCoord])[0]);
673  y[ii] = ((*coords[(jj + ii) % nCoord])[1]);
674  z[ii] = ((*coords[(jj + ii) % nCoord])[2]);
675  }
676  SurfCentroidConnec_ddf(x, y, z, edgeLength, centroid[0], centroid[1], centroid[2], temp);
677  ind = (jj + 1) % nCoord;
678  ind2 = (jj2 + 1) % nCoord;
679  // WARNING : This needs to be checked to make sure it matches the
680  // automaticly generated ones
681  for (ii1 = 0; ii1 < nDim; ii1++)
682  {
683  for (ii2 = 0; ii2 < nDim; ii2++)
684  {
685  for (ii = 0; ii < nFun; ii++)
686  {
687  hes[ind + (ii1 * nCoord)][ind2 + (ii2 * nCoord) + (ii * nDim * nCoord)] =
688  temp[ii1][nDim + ii2 + (nDim * ii)];
689  hes[ind2 + (ii1 * nCoord)][ind + (ii2 * nCoord) + (ii * nDim * nCoord)] =
690  temp[ii1 + nDim][ii2 + (nDim * ii)];
691  }
692  }
693  }
694  }
695 
696  // Not Connected terms of the Hessian
697  x.assign(8, 0);
698  y.assign(8, 0);
699  z.assign(8, 0);
700  temp.assign(nDim * 2, 2 * nDim * nFun, 0);
701  for (jj = 0; jj < nCoord; jj++)
702  {
703  for (jj2 = jj + 2; jj2 < nCoord; ++jj2)
704  {
705  for (ii = 0; ii < 3; ++ii)
706  {
707  x[ii] = ((*coords[(jj + ii) % nCoord])[0]);
708  y[ii] = ((*coords[(jj + ii) % nCoord])[1]);
709  z[ii] = ((*coords[(jj + ii) % nCoord])[2]);
710  }
711 
712  for (ii = 0; ii < 3; ++ii)
713  {
714  x[ii + 4] = ((*coords[(jj2 + ii) % nCoord])[0]);
715  y[ii + 4] = ((*coords[(jj2 + ii) % nCoord])[1]);
716  z[ii + 4] = ((*coords[(jj2 + ii) % nCoord])[2]);
717  }
718  SurfCentroidNoConnec_ddf(x, y, z, edgeLength, centroid[0], centroid[1], centroid[2], temp);
719  ind = (jj + 1) % nCoord;
720  ind2 = (jj2 + 1) % nCoord;
721  // WARNING : This needs to be checked to make sure it matches the
722  // automaticly generated ones
723  for (ii1 = 0; ii1 < nDim; ii1++)
724  {
725  for (ii2 = 0; ii2 < nDim; ii2++)
726  {
727  for (ii = 0; ii < nFun; ii++)
728  {
729  hes[ind + (ii1 * nCoord)][ind2 + (ii2 * nCoord) + (ii * nDim * nCoord)] =
730  temp[ii1][nDim + ii2 + (nDim * ii)];
731  hes[ind2 + (ii1 * nCoord)][ind + (ii2 * nCoord) + (ii * nDim * nCoord)] =
732  temp[ii1 + nDim][ii2 + (nDim * ii)];
733  }
734  }
735  }
736  }
737  }
738  }
739 }
740 
741 void SurfCentroid::Disp()
742 {
743  for (int ii = 0; ii < nCoord; ii++)
744  {
745  cout << "c " << ii << " ";
746  DisplayVector(*coords[ii]);
747  cout << endl;
748  }
749  cout << "edgeLength" << edgeLength << endl << "centroid ";
750  DisplayVector(centroid);
751  cout << endl;
752 }
753 
754 void Volume2::Calc()
755 {
756 #ifdef SAFE_ALGO
757  PreCalc();
758 #endif
759  if (!isCalc)
760  {
761  Volume2_f((*coords[0])[0], (*coords[0])[1], (*coords[0])[2], *coords[1], *coords[2], *coords[3], *coords[4],
762  *coords[5], *coords[6], fun);
763  Volume2_df((*coords[0])[0], (*coords[0])[1], (*coords[0])[2], *coords[1], *coords[2], *coords[3], *coords[4],
764  *coords[5], *coords[6], jac);
765  Volume2_ddf((*coords[0])[0], (*coords[0])[1], (*coords[0])[2], *coords[1], *coords[2], *coords[3], *coords[4],
766  *coords[5], *coords[6], hes);
767  }
768 }
769 
770 void SetEnvironmentEpsilon(const param::dev::rsvseps &rsvsEpsilons)
771 {
772  rsvsmath_automatic_eps_edge = rsvsEpsilons.rsvsmath_automatic_eps_edge;
773  rsvsmath_automatic_eps_surf = rsvsEpsilons.rsvsmath_automatic_eps_surf;
774  rsvsmath_automatic_eps_volu = rsvsEpsilons.rsvsmath_automatic_eps_volu;
775  rsvsmath_automatic_eps_centre = rsvsEpsilons.rsvsmath_automatic_eps_centre;
776  rsvsmath_automatic_eps_centre2 = rsvsEpsilons.rsvsmath_automatic_eps_centre2;
777 }
Performs Volume and Area calculations for the RSVS process.
Provide std::vector container with hashed index mapping.
bool MakeValidField(std::vector< double > const *mp)
CoordFunc supports the same stuff as tri func but can have any number of points.
Definition: RSVSmath.cpp:154
\TODO refactor mesh.hpp into mesh and coord headers to allow smaller includes
Definition: RSVSmath.hpp:38
Class for control of rsvs epsilon.
Definition: parameters.hpp:257
std::vector< const std::vector< double > * > coordlist
Defines a list of coordinates.
Definition: mesh.hpp:87
Parameters for the integrated 3DRSVS.
Provides the error and warning system used by the RSVS3D project.
#define RSVS3D_ERROR_ARGUMENT(M)
Throw a invalid_argument.
Definition: warning.hpp:148