20 #include "GmshConfig.h"
39 auto it =
all_.find(&l);
40 if(it ==
all_.end())
return nullptr;
49 int id1 =
box.getCellContainingPoint(x - rmax, y - rmax,
z - rmax);
50 int id2 =
box.getCellContainingPoint(x + rmax, y + rmax,
z + rmax);
52 box.getCellIJK(id1, i1, j1, k1);
54 box.getCellIJK(id2, i2, j2, k2);
55 for(
int i = i1; i <= i2; i++)
56 for(
int j = j1; j <= j2; j++)
57 for(
int k = k1; k <= k2; k++)
58 box.insertActiveCell(
box.getCellIndex(i, j, k));
64 double t_min = t_bounds.
low();
65 double t_max = t_bounds.
high();
68 for(
int i = 0; i < N; i++) {
69 double t = t_min + (double)i / (
double)(N - 1) * (t_max - t_min);
71 points.push_back(
SPoint3(p.
x(), p.
y(), p.
z()));
82 for(
int i = 0; i < I; i++)
83 for(
int j = 0; j < J; j++)
84 for(
int k = 0; k < K; k++) {
95 bool exist[8], atLeastOne =
false, butNotAll =
false;
96 for(
int ii = 0; ii < 8; ii++) {
98 if(exist[ii]) atLeastOne =
true;
99 if(!exist[ii]) butNotAll =
true;
122 if(!
box->getChildBox())
return;
123 for(
int i = 0; i <
box->getNxi(); i++)
124 for(
int j = 0; j <
box->getNeta(); j++)
125 for(
int k = 0; k <
box->getNzeta(); k++) {
126 if(
box->activeCellExists(
box->getCellIndex(i, j, k))) {
128 int ii = i, jj = j, kk = k;
133 if(child->activeCellExists(child->getCellIndex(ii, jj, kk))) {
134 box->eraseActiveCell(
box->getCellIndex(i, j, k));
146 std::vector<SPoint3> nodes;
147 std::vector<int> indices;
148 for(
auto it =
box.nodalValuesBegin(); it !=
box.nodalValuesEnd(); ++it) {
149 nodes.push_back(
box.getNodeCoordinates(it->first));
150 indices.push_back(it->first);
154 std::vector<double> dist, localdist;
155 std::vector<SPoint3> dummy;
158 for(std::size_t k = 0; k < (*fit)->getNumMeshElements(); k++) {
159 std::vector<double> iDistances;
160 std::vector<SPoint3> iClosePts;
161 std::vector<double> iDistancesE;
162 MElement *e = (*fit)->getMeshElement(k);
176 for(std::size_t j = 0; j < localdist.size(); j++) {
178 (fabs(dist[j]) < fabs(localdist[j])) ? dist[j] : localdist[j];
185 Msg::Info(
"Levelset computation on CAD surface not implemented");
189 for(std::size_t j = 0; j < dist.size(); j++)
190 box.setNodalValue(indices[j], dist[j]);
195 inline double det3(
double d11,
double d12,
double d13,
double d21,
double d22,
196 double d23,
double d31,
double d32,
double d33)
198 return d11 * (d22 * d33 - d23 * d32) - d21 * (d12 * d33 - d13 * d32) +
199 d31 * (d12 * d23 - d13 * d22);
204 double n = sqrt(vec[0] * vec[0] + vec[1] * vec[1] + vec[2] * vec[2]);
205 norm[0] = vec[0] / n;
206 norm[1] = vec[1] / n;
207 norm[2] = vec[2] / n;
210 inline void cross(
const double *pt0,
const double *pt1,
const double *pt2,
213 double v1[3] = {pt1[0] - pt0[0], pt1[1] - pt0[1], pt1[2] - pt0[2]};
214 double v2[3] = {pt2[0] - pt0[0], pt2[1] - pt0[1], pt2[2] - pt0[2]};
215 cross[0] = v1[1] * v2[2] - v2[1] * v1[2];
216 cross[1] = v2[0] * v1[2] - v1[0] * v2[2];
217 cross[2] = v1[0] * v2[1] - v2[0] * v1[1];
220 inline bool isPlanar(
const double *pt1,
const double *pt2,
const double *pt3,
226 cross(pt1, pt2, pt4, c2);
231 return (n1[0] == n2[0] && n1[1] == n2[1] && n1[2] == n2[2]);
237 double r2 = dx * dx + dy * dy + dz * dz;
241 case 0:
return exp(-ep * ep * r2);
242 case 1:
return -2 * ep * ep * dx * exp(-ep * ep * r2);
243 case 2:
return -2 * ep * ep * dy * exp(-ep * ep * r2);
244 case 3:
return -2 * ep * ep * dz * exp(-ep * ep * r2);
248 case 0:
return sqrt(1 + ep * ep * r2);
249 case 1:
return ep * ep * dx / sqrt(1 + ep * ep * r2);
250 case 2:
return ep * ep * dy / sqrt(1 + ep * ep * r2);
251 case 3:
return ep * ep * dz / sqrt(1 + ep * ep * r2);
259 FILE *xyz =
Fopen(
"myNodes.pos",
"w");
261 fprintf(xyz,
"View \"\"{\n");
262 for(
int itv = 1; itv != myNodes.
size1(); ++itv) {
263 fprintf(xyz,
"SP(%g,%g,%g){%g};\n", myNodes(itv, 0), myNodes(itv, 1),
264 myNodes(itv, 2), surf(itv, 0));
266 fprintf(xyz,
"};\n");
275 std::queue<gLevelset *> Q;
279 std::vector<gLevelset *> pp;
281 if(pp.empty()) gLsPrimitives.push_back(p);
284 for(
int i = 0; i < (int)pp.size(); i++) Q.push(pp[i]);
293 std::stack<gLevelset *>
S;
294 std::stack<gLevelset *> Sc;
298 std::vector<gLevelset *> pp;
301 gLsPrimitives.push_back(p);
305 if(Sc.empty() || p != Sc.top()) {
306 for(
int i = 1; i < (int)pp.size(); i++) Sc.push(p);
307 for(
int i = (
int)pp.size() - 1; i >= 0; i--) {
323 std::stack<gLevelset *>
S;
324 std::stack<gLevelset *> Sc;
328 std::vector<gLevelset *> pp;
335 if(Sc.empty() || p != Sc.top()) {
336 for(
int i = 1; i < (int)pp.size(); i++) Sc.push(p);
337 for(
int i = (
int)pp.size() - 1; i >= 0; i--) {
354 const double &
z,
const double &R,
int tag)
361 double &dfdy,
double &dfdz)
const
363 const double xx = x -
xc, yy = y -
yc, zz =
z -
zc,
364 dist = sqrt(xx * xx + yy * yy + zz * zz);
372 double &dfdxy,
double &dfdxz,
double &dfdyx,
373 double &dfdyy,
double &dfdyz,
double &dfdzx,
374 double &dfdzy,
double &dfdzz)
const
376 const double xx = x -
xc, yy = y -
yc, zz =
z -
zc;
377 const double distSq = xx * xx + yy * yy, fact = 1. / (distSq * sqrt(distSq));
379 dfdxx = (zz * zz + yy * yy) * fact;
380 dfdxy = -xx * yy * fact;
381 dfdxz = -xx * zz * fact;
383 dfdyy = (xx * xx + zz * zz) * fact;
384 dfdyz = -yy * zz * fact;
387 dfdzz = (xx * xx + yy * yy) * fact;
391 const std::vector<double> &
norm,
int tag)
397 d = -
a * pt[0] -
b * pt[1] -
c * pt[2];
406 d = -
a * pt[0] -
b * pt[1] -
c * pt[2];
410 const double *pt3,
int tag)
413 a =
det3(1., pt1[1], pt1[2], 1., pt2[1], pt2[2], 1., pt3[1], pt3[2]);
414 b =
det3(pt1[0], 1., pt1[2], pt2[0], 1., pt2[2], pt3[0], 1., pt3[2]);
415 c =
det3(pt1[0], pt1[1], 1., pt2[0], pt2[1], 1., pt3[0], pt3[1], 1.);
416 d = -
det3(pt1[0], pt1[1], pt1[2], pt2[0], pt2[1], pt2[2], pt3[0], pt3[1],
435 int m = nodes2.
size1();
436 int n = nodes1.
size1();
439 double eps = 0.5 /
delta;
440 for(
int i = 0; i < m; i++) {
441 for(
int j = 0; j < n; j++) {
442 double dx = nodes2(i, 0) - nodes1(j, 0);
443 double dy = nodes2(i, 1) - nodes1(j, 1);
444 double dz = nodes2(i, 2) - nodes1(j, 2);
467 D.gemm(rbfMatB, rbfInvA, 1.0, 0.0);
479 RbfOp(p, index, cntrs, nodes,
D, isLocal);
480 fApprox.
gemm(
D, fValues, 1.0, 0.0);
487 int numNodes = cntrs.
size1();
488 int nTot = 3 * numNodes;
490 level_set_nodes.
resize(nTot, 3);
491 level_set_funvals.
resize(nTot, 1);
494 cntrsPlus(numNodes + 1, 3);
497 double dist_min = 1.e6;
498 double dist_max = 1.e-6;
499 for(
int i = 0; i < numNodes; ++i) {
501 cntrsPlus(i, 0) = cntrs(i, 0);
502 cntrsPlus(i, 1) = cntrs(i, 1);
503 cntrsPlus(i, 2) = cntrs(i, 2);
504 for(
int j = i + 1; j < numNodes; j++) {
506 sqrt((cntrs(i, 0) - cntrs(j, 0)) * (cntrs(i, 0) - cntrs(j, 0)) +
507 (cntrs(i, 1) - cntrs(j, 1)) * (cntrs(i, 1) - cntrs(j, 1)) +
508 (cntrs(i, 2) - cntrs(j, 2)) * (cntrs(i, 2) - cntrs(j, 2)));
509 if(dist < dist_min) dist_min = dist;
510 if(dist > dist_max) dist_max = dist;
513 ONES(numNodes, 0) = -1.0;
514 cntrsPlus(numNodes, 0) = cntrs(0, 0) + dist_max;
515 cntrsPlus(numNodes, 1) = cntrs(0, 1) + dist_max;
516 cntrsPlus(numNodes, 2) = cntrs(0, 2) + dist_max;
518 delta = 0.23 * dist_min;
519 evalRbfDer(1, 1, cntrsPlus, cntrs, ONES, sx,
true);
520 evalRbfDer(2, 1, cntrsPlus, cntrs, ONES, sy,
true);
521 evalRbfDer(3, 1, cntrsPlus, cntrs, ONES, sz,
true);
522 for(
int i = 0; i < numNodes; ++i) {
524 sqrt(sx(i, 0) * sx(i, 0) + sy(i, 0) * sy(i, 0) + sz(i, 0) * sz(i, 0));
525 sx(i, 0) = sx(i, 0) / normFactor;
526 sy(i, 0) = sy(i, 0) / normFactor;
527 sz(i, 0) = sz(i, 0) / normFactor;
528 norms(i, 0) = sx(i, 0);
529 norms(i, 1) = sy(i, 0);
530 norms(i, 2) = sz(i, 0);
533 for(
int i = 0; i < numNodes; ++i) {
534 for(
int j = 0; j < 3; ++j) {
535 level_set_nodes(i, j) = cntrs(i, j);
536 level_set_nodes(i + numNodes, j) = cntrs(i, j) -
delta * norms(i, j);
537 level_set_nodes(i + 2 * numNodes, j) = cntrs(i, j) +
delta * norms(i, j);
539 level_set_funvals(i, 0) = 0.0;
540 level_set_funvals(i + numNodes, 0) = -1;
541 level_set_funvals(i + 2 * numNodes, 0) = 1;
548 int nbNodes = 3 * centers.
size1();
569 Msg::Info(
"Levelset Points : call computeLS() before calling operator()\n");
572 auto it =
mapP.find(sp);
573 if(it !=
mapP.end())
return it->second;
574 printf(
"Levelset Points : Point not found\n");
588 for(std::size_t i = 0; i < vert.size(); i++) {
589 xyz_eval(i, 0) = vert[i]->x();
590 xyz_eval(i, 1) = vert[i]->y();
591 xyz_eval(i, 2) = vert[i]->z();
594 for(std::size_t i = 0; i < vert.size(); i++) {
595 mapP[
SPoint3(vert[i]->x(), vert[i]->y(), vert[i]->
z())] = surf_eval(i, 0);
613 for(
int i = 0; i < 3; i++) {
615 for(
int j = 0; j < 3; j++)
A[i][j] = lv.
A[i][j];
622 for(
int i = 0; i < 3; i++) {
624 for(
int j = 0; j < 3; j++) { res[i] +=
A[i][j] * x[j] * fact; }
630 res = fact * (
A[0][0] * x[0] * x[0] +
A[1][1] * x[1] * x[1] +
631 A[2][2] * x[2] * x[2] +
A[1][0] * x[1] * x[0] * 2. +
632 A[2][0] * x[2] * x[0] * 2. +
A[1][2] * x[1] * x[2] * 2.);
639 C += (-
B[0] * transl[0] -
B[1] * transl[1] -
B[2] * transl[2] + b);
650 double a11 = 0., a12 = 0., a13 = 0., a22 = 0., a23 = 0., a33 = 0.;
651 double b1 = 0., b2 = 0., b3 = 0.;
652 for(
int i = 0; i < 3; i++) {
656 for(
int j = 0; j < 3; j++) {
666 A[0][1] =
A[1][0] = a12;
667 A[0][2] =
A[2][0] = a13;
669 A[1][2] =
A[2][1] = a23;
680 double norm = sqrt(dir[0] * dir[0] + dir[1] * dir[1]);
681 double n[3] = {1., 0., 0.};
682 double ct = 1., st = 0.;
684 double theta = atan(
norm / dir[2]);
685 n[0] = dir[1] /
norm;
686 n[1] = -dir[0] /
norm;
690 t[0][0] = n[0] * n[0] * (1. - ct) + ct;
691 t[0][1] = n[0] * n[1] * (1. - ct) - n[2] * st;
692 t[0][2] = n[0] * n[2] * (1. - ct) + n[1] * st;
693 t[1][0] = n[1] * n[0] * (1. - ct) + n[2] * st;
694 t[1][1] = n[1] * n[1] * (1. - ct) + ct;
695 t[1][2] = n[1] * n[2] * (1. - ct) - n[0] * st;
696 t[2][0] = n[2] * n[0] * (1. - ct) - n[1] * st;
697 t[2][1] = n[2] * n[1] * (1. - ct) + n[0] * st;
698 t[2][2] = n[2] * n[2] * (1. - ct) + ct;
702 const double dir2[3],
705 double norm = sqrt((dir1[1] * dir2[2] - dir1[2] * dir2[1]) *
706 (dir1[1] * dir2[2] - dir1[2] * dir2[1]) +
707 (dir1[2] * dir2[0] - dir1[0] * dir2[2]) *
708 (dir1[2] * dir2[0] - dir1[0] * dir2[2]) +
709 (dir1[0] * dir2[1] - dir1[1] * dir2[0]) *
710 (dir1[0] * dir2[1] - dir1[1] * dir2[0]));
711 double n[3] = {1., 0., 0.};
712 double ct = 1., st = 0.;
714 st =
norm / ((dir1[0] * dir1[0] + dir1[1] * dir1[1] + dir1[2] * dir1[2]) *
715 (dir2[0] * dir2[0] + dir2[1] * dir2[1] + dir2[2] * dir2[2]));
716 n[0] = (dir1[1] * dir2[2] - dir1[2] * dir2[1]) /
norm;
717 n[1] = (dir1[2] * dir2[0] - dir1[0] * dir2[2]) /
norm;
718 n[2] = (dir1[0] * dir2[1] - dir1[1] * dir2[0]) /
norm;
721 t[0][0] = n[0] * n[0] * (1. - ct) + ct;
722 t[0][1] = n[0] * n[1] * (1. - ct) - n[2] * st;
723 t[0][2] = n[0] * n[2] * (1. - ct) + n[1] * st;
724 t[1][0] = n[1] * n[0] * (1. - ct) + n[2] * st;
725 t[1][1] = n[1] * n[1] * (1. - ct) + ct;
726 t[1][2] = n[1] * n[2] * (1. - ct) - n[0] * st;
727 t[2][0] = n[2] * n[0] * (1. - ct) - n[1] * st;
728 t[2][1] = n[2] * n[1] * (1. - ct) + n[0] * st;
729 t[2][2] = n[2] * n[2] * (1. - ct) + ct;
734 for(
int i = 0; i < 3; i++) {
736 for(
int j = 0; j < 3; j++)
A[i][j] = 0.;
743 return (
A[0][0] * x * x + 2. *
A[0][1] * x * y + 2. *
A[0][2] * x *
z +
744 A[1][1] * y * y + 2. *
A[1][2] * y *
z +
A[2][2] *
z *
z +
B[0] * x +
745 B[1] * y +
B[2] *
z +
C);
749 double _a,
double _b,
int _c,
int tag)
755 while(
angle <= 2. * M_PI) {
759 angle += 2. * M_PI / 1000.;
766 double dx, dy, xi, yi, d;
769 double distance = sqrt(dx * dx + dy * dy);
770 auto itx =
iso_x.begin();
771 auto itxen =
iso_x.end();
772 auto ity =
iso_y.begin();
773 auto itminx =
iso_x.begin();
774 auto itminy =
iso_y.begin();
777 for(; itx != itxen; itx++, ity++) {
782 d = sqrt(dx * dx + dy * dy);
792 p(0) = (*itminx) - x;
793 p(1) = (*itminy) - y;
794 if(itminx == (itxen - 1)) {
795 t(0) =
iso_x[0] - (*itminx);
796 t(1) =
iso_y[0] - (*itminy);
799 t(0) = (*(itminx + 1)) - (*itminx);
800 t(1) = (*(itminy + 1)) - (*itminy);
802 double vectprod = (p(0) * t(1) - p(1) * t(0));
803 double sign = (vectprod < 0.) ? -1. : 1.;
809 double myr0,
double myA,
double mysigma,
825 sqrt((x -
xc) * (x -
xc) + (y -
yc) * (y -
yc) + (
z -
zc) * (
z -
zc));
828 for(
int k = 0; k < 5; k++) {
829 double xk =
r0 / sqrt(5.0) * (2. * cos(2 * k * M_PI / 5.0));
830 double yk =
r0 / sqrt(5.0) * (2. * sin(2 * k * M_PI / 5.0));
831 double zk =
r0 / sqrt(5.0);
833 A * exp(-((x -
xc - xk) * (x -
xc - xk) + (y -
yc - yk) * (y -
yc - yk) +
834 (
z -
zc - zk) * (
z -
zc - zk)) /
837 for(
int k = 5; k < 10; k++) {
838 double xk =
r0 / sqrt(5.0) * (2. * cos((2. * (k - 5.) - 1.) * M_PI / 5.0));
839 double yk =
r0 / sqrt(5.0) * (2. * sin((2. * (k - 5.) - 1.) * M_PI / 5.0));
840 double zk = -
r0 / sqrt(5.0);
842 A * exp(-((x -
xc - xk) * (x -
xc - xk) + (y -
yc - yk) * (y -
yc - yk) +
843 (
z -
zc - zk) * (
z -
zc - zk)) /
846 val -=
A * exp(-((x -
xc) * (x -
xc) + (y -
yc) * (y -
yc) +
849 val -=
A * exp(-((x -
xc) * (x -
xc) + (y -
yc) * (y -
yc) +
858 std::vector<std::string> expressions(1,
f);
859 std::vector<std::string> variables(3);
873 std::vector<double> values(3), res(1);
877 if(
_expr->
eval(values, res))
return res[0];
886 std::vector<std::string> variables(3);
900 std::vector<double> values(3), res(13);
904 if(
_expr->
eval(values, res))
return res[0];
909 double &dfdy,
double &dfdz)
const
911 std::vector<double> values(3), res(13);
923 double &dfdxy,
double &dfdxz,
double &dfdyx,
924 double &dfdyy,
double &dfdyz,
double &dfdzx,
925 double &dfdzy,
double &dfdzz)
const
927 std::vector<double> values(3), res(13);
944 #if defined(HAVE_ANN)
945 gLevelsetDistMesh::gLevelsetDistMesh(
GModel *gm,
const std::string &physical,
946 int nbClose,
int tag)
949 std::map<int, std::vector<GEntity *> > groups[4];
952 if(it->second == physical) {
953 _entities = groups[it->first.first][it->first.second];
956 if(_entities.size() == 0) {
958 "distanceToMesh: the physical name '%s' does not exist in the GModel",
963 std::set<MVertex *> _all;
964 for(std::size_t i = 0; i < _entities.size(); i++) {
965 for(std::size_t k = 0; k < _entities[i]->getNumMeshElements(); k++) {
966 MElement *e = _entities[i]->getMeshElement(k);
968 MVertex *v = _entities[i]->getMeshElement(k)->getVertex(j);
970 _v2e.insert(std::make_pair(v, e));
975 nodes = annAllocPts(_all.size(), 3);
976 auto itp = _all.begin();
978 while(itp != _all.end()) {
980 nodes[ind][0] = pt->
x();
981 nodes[ind][1] = pt->
y();
982 nodes[ind][2] = pt->
z();
983 _vertices.push_back(pt);
987 _kdtree =
new ANNkd_tree(nodes, _all.size(), 3);
990 gLevelsetDistMesh::~gLevelsetDistMesh()
993 ANNpointArray nodes = _kdtree->thePoints();
994 annDeallocPts(nodes);
999 double gLevelsetDistMesh::operator()(
double x,
double y,
double z)
const
1001 std::vector<ANNidx> index(_nbClose);
1002 std::vector<ANNdist> dist(_nbClose);
1003 double point[3] = {x, y,
z};
1005 _kdtree->annkSearch(point, _nbClose, &index[0],
1007 std::set<MElement *> elements;
1009 for(
int i = 0; i < _nbClose; i++) {
1010 int iVertex = index[i];
1011 MVertex *v = _vertices[iVertex];
1012 for(
auto itm = _v2e.lower_bound(v); itm != _v2e.upper_bound(v); ++itm) {
1013 elements.insert(itm->second);
1014 dimE = itm->second->getDim();
1017 double minDistance = 1.e22;
1019 std::vector<MElement *> closestElements;
1020 for(
auto it = elements.begin(); it != elements.end(); ++it) {
1022 MVertex *v1 = (*it)->getVertex(0);
1023 MVertex *v2 = (*it)->getVertex(1);
1030 else if(dimE == 2) {
1031 MVertex *v3 = (*it)->getVertex(2);
1033 if(p1 == p2 || p1 == p3 || p2 == p3)
1039 Msg::Error(
"Cannot compute a distance to an entity of dimension %d\n",
1042 if(dimE == 1 && fabs(
distance) < fabs(minDistance)) {
1045 else if(dimE == 2) {
1046 if(fabs(
distance) - fabs(minDistance) < 1.e-9) {
1047 closestElements.push_back(*it);
1049 else if(fabs(
distance) < fabs(minDistance)) {
1052 closestElements.clear();
1053 closestElements.push_back(*it);
1057 if(closestElements.size() > 1) {
1060 if(closestElements.size() == 2) {
1061 meanNorm = closestElements[0]->getFace(0).normal() +
1062 closestElements[1]->getFace(0).normal();
1065 for(std::size_t i = 0; i < closestElements.size(); i++) {
1069 for(
int j = 0; j < closestElements[i]->getNumEdges(); j++) {
1070 SPoint3 ep0 = closestElements[i]->getEdge(j).getVertex(0)->point();
1071 SPoint3 ep1 = closestElements[i]->getEdge(j).getVertex(1)->point();
1095 meanNorm = meanNorm + alpha * closestElements[i]->getFace(0).normal();
1098 double dotDN =
dot(vd, closestElements[0]->getFace(0).normal());
1099 double dotDE =
dot(vd, meanNorm);
1100 if(dotDN * dotDE < 0.) minDistance *= -1.;
1103 return -1.0 * minDistance;
1107 #if defined(HAVE_POST)
1108 gLevelsetPostView::gLevelsetPostView(
int index,
int tag)
1111 if(_viewIndex >= 0 && _viewIndex < (
int)
PView::list.size()) {
1116 Msg::Error(
"Unknown View[%d] in PostView levelset", _viewIndex);
1121 double gLevelsetPostView::operator()(
double x,
double y,
double z)
const
1123 if(!_octree)
return 1.;
1125 _octree->searchScalar(x, y,
z, &val, 0);
1131 const double &R,
int tag)
1149 const double &a,
const double &b,
1150 const double &
c,
int tag)
1153 A[0][0] = 1. / (a * a);
1154 A[1][1] = 1. / (b * b);
1155 A[2][2] = 1. / (
c *
c);
1169 const double &
angle,
int tag)
1184 const double *pt,
const double *dir,
const double &x2,
const double &y2,
1185 const double &z2,
const double &
z,
const double &
c,
int tag)
1207 std::vector<gLevelset *> _children = lv.
getChildren();
1208 unsigned siz = _children.size();
1210 for(
unsigned i = 0; i < siz; ++i)
children[i] = _children[i]->
clone();
1217 std::map<int, std::vector<GEntity *> > groups;
1221 printf(
"No physical %d found for levelset yarn!\n", phys);
1227 for(std::size_t i = 0; i <
entities.size(); i++) {
1257 const double *dir2,
const double *dir3,
1258 const double &a,
const double &b,
const double &
c,
1262 double dir1m[3] = {-dir1[0], -dir1[1], -dir1[2]};
1263 double dir2m[3] = {-dir2[0], -dir2[1], -dir2[2]};
1264 double dir3m[3] = {-dir3[0], -dir3[1], -dir3[2]};
1271 double pt2[3] = {pt[0] + a * n1[0] + b * n2[0] +
c * n3[0],
1272 pt[1] + a * n1[1] + b * n2[1] +
c * n3[1],
1273 pt[2] + a * n1[2] + b * n2[2] +
c * n3[2]};
1274 std::vector<gLevelset *> p;
1285 const double *pt3,
const double *pt4,
1286 const double *pt5,
const double *pt6,
1287 const double *pt7,
const double *pt8,
int tag)
1294 "WARNING : faces of the box are not planar! %d, %d, %d, %d, %d, %d\n",
1298 std::vector<gLevelset *> p;
1311 const std::vector<double> &dir,
1312 const double &R,
const double &H,
int tag)
1315 double pt1[3] = {pt[0], pt[1], pt[2]};
1316 double dir1[3] = {dir[0], dir[1], dir[2]};
1317 double dir2[3] = {-dir1[0], -dir1[1], -dir1[2]};
1320 double pt2[3] = {pt1[0] +
H * n[0], pt1[1] +
H * n[1], pt1[2] +
H * n[2]};
1321 std::vector<gLevelset *> p;
1329 const double &R,
const double &H,
int tag)
1332 double dir2[3] = {-dir[0], -dir[1], -dir[2]};
1335 double pt2[3] = {pt[0] +
H * n[0], pt[1] +
H * n[1], pt[2] +
H * n[2]};
1336 std::vector<gLevelset *> p;
1344 const double &R,
const double &r,
1345 const double &H,
int tag)
1348 double dir2[3] = {-dir[0], -dir[1], -dir[2]};
1351 double pt2[3] = {pt[0] +
H * n[0], pt[1] +
H * n[1], pt[2] +
H * n[2]};
1352 std::vector<gLevelset *> p1;
1356 std::vector<gLevelset *> p2;
1368 const double *dir2,
const double &H1,
1369 const double &H2,
const double &H3,
1370 const double &R1,
const double &r1,
1371 const double &R2,
const double &r2,
1372 const double &L1,
const double &L2,
1373 const double &E,
int tag)
1380 double pt1[3] = {pt[0] - n2[0] * H1 / 2., pt[1] - n2[1] * H1 / 2.,
1381 pt[2] - n2[2] * H1 / 2.};
1382 double pt2[3] = {pt[0] + n1[0] * E - n2[0] * H2 / 2.,
1383 pt[1] + n1[1] * E - n2[1] * H2 / 2.,
1384 pt[2] + n1[2] * E - n2[2] * H2 / 2.};
1386 cross(pt1, pt2, pt, dir3);
1389 double pt31[3] = {pt[0] - n2[0] * H3 / 2. + n3[0] * L1 / 2.,
1390 pt[1] - n2[1] * H3 / 2. + n3[1] * L1 / 2.,
1391 pt[2] - n2[2] * H3 / 2. + n3[2] * L1 / 2.};
1392 double pt32[3] = {pt31[0] - n3[0] * L1, pt31[1] - n3[1] * L1,
1393 pt31[2] - n3[2] * L1};
1394 double pt33[3] = {pt32[0] + n2[0] * H3, pt32[1] + n2[1] * H3,
1395 pt32[2] + n2[2] * H3};
1396 double pt34[3] = {pt33[0] + n3[0] * L1, pt33[1] + n3[1] * L1,
1397 pt33[2] + n3[2] * L1};
1398 double pt35[3] = {pt[0] + n1[0] * E - n2[0] * H3 / 2. + n3[0] * L2 / 2.,
1399 pt[1] + n1[1] * E - n2[1] * H3 / 2. + n3[1] * L2 / 2.,
1400 pt[2] + n1[2] * E - n2[2] * H3 / 2. + n3[2] * L2 / 2.};
1401 double pt36[3] = {pt35[0] - n3[0] * L2, pt35[1] - n3[1] * L2,
1402 pt35[2] - n3[2] * L2};
1403 double pt37[3] = {pt36[0] + n2[0] * H3, pt36[1] + n2[1] * H3,
1404 pt36[2] + n2[2] * H3};
1405 double pt38[3] = {pt37[0] + n3[0] * L2, pt37[1] + n3[1] * L2,
1406 pt37[2] + n3[2] * L2};
1407 std::vector<gLevelset *> p1;
1409 new gLevelsetBox(pt31, pt32, pt33, pt34, pt35, pt36, pt37, pt38, tag));
1412 std::vector<gLevelset *> p2;
1427 : _x0(x0), _y0(y0), _c(
c), _t(t)
1433 double &xb,
double &yb,
1434 double &curvRad,
bool &in)
const
1436 static const int maxIter = 100;
1437 static const double tol = 1.e-10;
1439 const double tolr = tol /
_c;
1443 const double xt = x -
_x0, yt = fabs(y -
_y0);
1445 if(xt -
_c > 1.21125 *
_t * yt) {
1453 const double fact = 5. *
_t *
_c;
1454 double xtb = std::max(xt, tolr), ytb;
1456 for(
int it = 0; it < maxIter; it++) {
1457 const double xbr = xtb /
_c, sxbr = sqrt(xbr), xbr32 = xbr * sxbr,
1458 xbr2 = xbr * xbr, xbr3 = xbr2 * xbr, xbr4 = xbr2 * xbr2;
1459 ytb = fact * (0.2969 * sxbr - 0.1260 * xbr - 0.3516 * xbr2 +
1460 0.2843 * xbr3 - 0.1036 * xbr4);
1462 (0.14845 / sxbr - 0.4144 * xbr3 + 0.8529 * xbr2 - 0.7032 * xbr -
1466 (-0.074225 / xbr32 - 1.2432 * xbr2 + 1.7058 * xbr - 0.7032) /
1468 const double xx = xt - xtb, yy = yt - ytb;
1469 in = (xt > 0) && (yy < 0);
1470 const double dDistSq = -2. * (xx + dyb * yy);
1471 const double ddDistSq = 2. * (1. - ddyb * yy + dyb * dyb);
1472 const double mIncr = dDistSq / ddDistSq;
1473 if(fabs(mIncr) < tolr)
1477 if(xtb < tolr) xtb = tolr;
1478 if(xtb >
_c - tolr) xtb =
_c - tolr;
1481 yb = (y >=
_y0) ?
_y0 + ytb :
_y0 - ytb;
1482 const double norm = sqrt(1. + dyb * dyb);
1489 double xb, yb, curvRadb;
1493 const double xx = x - xb, yy = y - yb, distSq = xx * xx + yy * yy;
1495 return in ? -sqrt(distSq) : sqrt(distSq);
1499 double &dfdy,
double &dfdz)
const
1501 double xb, yb, curvRadb;
1505 const double xx = x - xb, yy = y - yb, distSq = xx * xx + yy * yy;
1506 const double dist = in ? -sqrt(distSq) : sqrt(distSq);
1514 double &dfdxy,
double &dfdxz,
double &dfdyx,
1515 double &dfdyy,
double &dfdyz,
double &dfdzx,
1516 double &dfdzy,
double &dfdzz)
const
1518 double xb, yb, curvRadb;
1523 const double xx = x - xb, yy = y - yb, distSq = xx * xx + yy * yy;
1524 const double dist = in ? -sqrt(distSq) : sqrt(distSq);
1525 const double curvRad = curvRadb + dist,
1526 fact = 1. / (curvRad * curvRad * curvRad);
1528 dfdxx = yy * yy * fact;
1529 dfdxy = -xx * yy * fact;
1532 dfdyy = xx * xx * fact;