Belle II Software  release-08-01-10
BelleLathe.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 /* Own header. */
10 #include <ecl/geometry/BelleLathe.h>
11 
12 /* Geant4 headers. */
13 #include <G4AffineTransform.hh>
14 #include <G4VoxelLimits.hh>
15 #include <G4VGraphicsScene.hh>
16 #include <G4VPVParameterisation.hh>
17 
18 /* CLHEP headers. */
19 #include <CLHEP/Random/RandFlat.h>
20 
21 using namespace std;
22 using namespace Belle2;
23 using namespace ECL;
24 
25 #define COMPARE 0
26 #define PERFCOUNTER 0
27 #if PERFCOUNTER==1
28 typedef int counter_t[6];
29 map<string, counter_t> counterl;
30 #define COUNTER(x) counterl[GetName()][x]++
31 //#define MATCHOUT(x) //if(GetName().find("sv_crystalcontainersolid")==0) cout<<x<<endl;
32 #else
33 #define COUNTER(x)
34 //
35 #endif
36 
37 //#define MATCHOUT(x) G4cout<<GetName()<<" "<<x<<G4endl;
38 #define MATCHOUT(x)
39 
41 struct Plane_t {
42  G4ThreeVector n;
43  double d;
44  // => n.x*x + n.y*y + n.z*z + d = 0
45 };
46 
47 namespace Belle2 {
52  namespace ECL {
53  inline double dotxy(const G4ThreeVector& p, const G4ThreeVector& n)
54  {
55  return p.x() * n.x() + p.y() * n.y();
56  }
57  }
59 }
60 
61 ostream& operator <<(ostream& o, const zr_t& v)
62 {
63  return o << "{" << v.z << ", " << v.r << "}";
64 }
65 
67 struct curl_t {
68  G4ThreeVector v;
70  explicit curl_t(const G4ThreeVector& _v): v(_v) {}
71 };
72 
73 ostream& operator <<(ostream& o, const curl_t& c)
74 {
75  return o << "{" << c.v.x() << ", " << c.v.y() << ", " << c.v.z() << "}, ";
76 }
77 
78 
79 BelleLathe::BelleLathe(const G4String& pName, double phi0, double dphi, const vector<zr_t>& c)
80  : G4CSGSolid(pName)
81 {
82  Init(c, phi0, dphi);
83 }
84 
85 BelleLathe::BelleLathe(const G4String& pName, double phi0, double dphi, int n, double* z, double* rin, double* rout)
86  : G4CSGSolid(pName)
87 {
88  vector<zr_t> contour;
89  for (int i = 0; i < n; i++) {
90  zr_t t = {z[i], rin[i]};
91  contour.push_back(t);
92  }
93  for (int i = n - 1; i >= 0; i--) {
94  zr_t t = {z[i], rout[i]};
95  contour.push_back(t);
96  }
97 
98  Init(contour, phi0, dphi);
99 }
100 
101 void BelleLathe::Init(const vector<zr_t>& c, double phi0, double dphi)
102 {
103  vector<zr_t> contour = c;
104  // remove duplicated vertices
105  do {
106  vector<zr_t>::iterator it0 = contour.begin(), it1 = it0 + 1;
107  for (; it1 != contour.end();) {
108  const zr_t& s0 = *it0, &s1 = *it1;
109  if (abs(s0.z - s1.z) < kCarTolerance && abs(s0.r - s1.r) < kCarTolerance)
110  it1 = contour.erase(it1);
111  else {
112  ++it0; ++it1;
113  }
114  }
115  const zr_t& s0 = *it0, &s1 = contour[0]; // cppcheck-suppress invalidContainer ; contour should be valid here
116  if (abs(s0.z - s1.z) < kCarTolerance && abs(s0.r - s1.r) < kCarTolerance) contour.erase(it0);
117  } while (0);
118 
119  // remove vertices on the same line
120  do {
121  vector<zr_t>::iterator it0 = contour.begin(), it1 = it0 + 1, it2 = it1 + 1;
122  for (; it0 != contour.end();) {
123  const zr_t& s0 = *it0, &s1 = *it1, &s2 = *it2;
124  double dr2 = s2.r - s0.r, dz2 = s2.z - s0.z;
125  double d = (s1.z - s0.z) * dr2 - (s1.r - s0.r) * dz2;
126 
127  if (d * d < kCarTolerance * kCarTolerance * (dr2 * dr2 + dz2 * dz2)) {
128  it1 = contour.erase(it1);
129  it2 = it1;
130  if (++it2 >= contour.end()) it2 = contour.begin();
131  it0 = it1;
132  if (--it0 < contour.begin()) it0 = (++contour.rbegin()).base();
133 
134  } else {
135  ++it0;
136  if (++it1 >= contour.end()) it1 = contour.begin();
137  if (++it2 >= contour.end()) it2 = contour.begin();
138  }
139 
140  }
141  } while (0);
142 
143  double sum = 0;
144  zr_t p0 = contour[0];
145  for (int i = 1, imax = contour.size(); i < imax; i++) {
146  zr_t p1 = contour[i];
147  sum += (p1.z - p0.z) * (p1.r + p0.r);
148  p0 = p1;
149  }
150  zr_t p1 = contour[0];
151  sum += (p1.z - p0.z) * (p1.r + p0.r);
152 
153  // If contour is Clockwise: reverse contour
154  if (sum > 0)
155  std::reverse(contour.begin(), contour.end());
156 
157  fcontour = contour;
158 
159  auto convexside = [this](cachezr_t& s, double eps) -> void {
160  s.isconvex = false;
161  if (s.dz > 0) return;
162  vector<zr_t>::const_iterator it = fcontour.begin();
163  double a = s.dz * s.is, b = s.dr * s.is, cc = b * s.z - a * s.r;
164  bool dp = false, dm = false;
165  s.isconvex = true;
166  do
167  {
168  const zr_t& p = *it;
169  double d = a * p.r - b * p.z + cc; // distance to line
170  dm = dm || (d < -eps);
171  dp = dp || (d > eps);
172  if (dm && dp) {s.isconvex = false; return;}
173  } while (++it != fcontour.end());
174  };
175 
176  frmin = kInfinity;
177  frmax = -kInfinity;
178  fzmin = kInfinity;
179  fzmax = -kInfinity;
180  fcache.reserve(fcontour.size());
181  for (int i = 0, n = fcontour.size(); i < n; i++) {
182  const zr_t& s0 = fcontour[i], &s1 = fcontour[(i + 1) % n];
183  cachezr_t t;
184  t.z = s0.z;
185  t.r = s0.r;
186  t.dz = s1.z - s0.z;
187  t.dr = s1.r - s0.r;
188  t.s2 = t.dz * t.dz + t.dr * t.dr;
189  t.is2 = 1 / t.s2;
190  t.is = sqrt(t.is2);
191  t.zmin = min(s0.z, s1.z);
192  t.zmax = max(s0.z, s1.z);
193  t.r2min = pow(min(s0.r, s1.r), 2);
194  t.r2max = pow(max(s0.r, s1.r), 2);
195  t.ta = (s1.r - s0.r) / (s1.z - s0.z);
196  convexside(t, kCarTolerance);
197  fcache.push_back(t);
198 
199  frmax = max(frmax, s0.r);
200  frmin = min(frmin, s0.r);
201  fzmax = max(fzmax, s0.z);
202  fzmin = min(fzmin, s0.z);
203  }
204 
205  fphi = phi0;
206  fdphi = dphi;
207 
208  fdphi = std::min(2 * M_PI, fdphi);
209  fdphi = std::max(0.0, fdphi);
210  fc0 = cos(fphi);
211  fs0 = sin(fphi);
212  fc1 = cos(fphi + fdphi);
213  fs1 = sin(fphi + fdphi);
214 
215  fn0x = fs0;
216  fn0y = -fc0;
217  fn1x = -fs1;
218  fn1y = fc1;
219  fgtpi = fdphi > M_PI;
220  ftwopi = abs(fdphi - 2 * M_PI) < kCarTolerance;
221 
222  // cout << ftwopi << " " << fgtpi << " " << fn0y << " " << fn0x << " " << fn1y << " " << fn1x << endl;
223 
224  for (int i = 0, n = fcontour.size(); i < n; i++) {
225  const zr_t& s = fcontour[i];
226  fz.push_back(s.z);
227  }
228  sort(fz.begin(), fz.end());
229  fz.erase(std::unique(fz.begin(), fz.end()), fz.end());
230 
231  for (int i = 1, ni = fz.size(); i < ni; i++) {
232  double a = fz[i - 1], b = fz[i];
233  findx.push_back(fseg.size());
234  for (int j = 0, nj = fcache.size(); j < nj; j++) {
235  const cachezr_t& sj = fcache[j];
236  double cc = sj.zmin, d = sj.zmax;
237  if (cc != d and b > cc and d > a) { // overlap
238  fseg.push_back(j);
239  }
240  }
241  }
242  findx.push_back(fseg.size());
243 
244  getvolarea();
245 
246 #if COMPARE>0
247  auto getpolycone = [](const G4String & pName, double phi0, double dphi, const vector<zr_t>& c) -> G4GenericPolycone* {
248  vector<double> r, z;
249  r.reserve(c.size());
250  z.reserve(c.size());
251  for (int i = 0, imax = c.size(); i < imax; i++)
252  {
253  r.push_back(c[i].r);
254  z.push_back(c[i].z);
255  }
256  return new G4GenericPolycone(pName, phi0, dphi, c.size(), r.data(), z.data());
257  };
258  fshape = getpolycone(GetName(), phi0, dphi, fcontour);
259 #else
260  fshape = nullptr;
261 #endif
262 // StreamInfo(G4cout);
263 }
264 
265 // Nominal constructor for BelleLathe whose parameters are to be set by
266 // a G4VParamaterisation later. Check and set half-widths as well as
267 // angles: final check of coplanarity
268 BelleLathe::BelleLathe(const G4String& pName)
269  : G4CSGSolid(pName)
270 {
271  vector<zr_t> a;
272  Init(a, 0, 2 * M_PI);
273 }
274 
275 // Fake default constructor - sets only member data and allocates memory
276 // for usage restricted to object persistency.
278  : G4CSGSolid(a)
279 {
280  vector<zr_t> b;
281  Init(b, 0, 2 * M_PI);
282 }
283 
284 // Destructor
286 {
287 #if PERFCOUNTER==1
288  cout << GetName() << " ";
289  for (int i = 0; i < 6; i++) cout << counterl[GetName()][i] << " "; cout << endl;
290 #endif
291 }
292 
293 // Copy constructor
295  : G4CSGSolid(rhs), fcontour(rhs.fcontour), fcache(rhs.fcache), fz(rhs.fz),
296  findx(rhs.findx), fseg(rhs.fseg), farea(rhs.farea), ftlist(rhs.ftlist),
297  fphi(rhs.fphi), fdphi(rhs.fdphi), fs0(rhs.fs0), fc0(rhs.fc0), fs1(rhs.fs1),
298  fc1(rhs.fc1), fn0x(rhs.fn0x), fn0y(rhs.fn0y), fn1x(rhs.fn1x), fn1y(rhs.fn1y),
299  frmin(rhs.frmin), frmax(rhs.frmax), fzmin(rhs.fzmin), fzmax(rhs.fzmax),
300  fgtpi(rhs.fgtpi), ftwopi(rhs.ftwopi), fshape(rhs.fshape), fsurf(rhs.fsurf)
301 {
302 }
303 
304 // Assignment operator
306 {
307  // Check assignment to self
308  if (this == &rhs) { return *this; }
309 
310  // Copy base class data
311  G4CSGSolid::operator=(rhs);
312 
313  // Copy data
314  fcontour = rhs.fcontour;
315  fcache = rhs.fcache;
316  fz = rhs.fz;
317  findx = rhs.findx;
318  fseg = rhs.fseg;
319  farea = rhs.farea;
320  ftlist = rhs.ftlist;
321  fphi = rhs.fphi;
322  fdphi = rhs.fdphi;
323  fs0 = rhs.fs0;
324  fc0 = rhs.fc0;
325  fs1 = rhs.fs1;
326  fc1 = rhs.fc1;
327  fn0x = rhs.fn0x;
328  fn0y = rhs.fn0y;
329  fn1x = rhs.fn1x;
330  fn1y = rhs.fn1y;
331  frmin = rhs.frmin;
332  frmax = rhs.frmax;
333  fzmin = rhs.fzmin;
334  fzmax = rhs.fzmax;
335  fgtpi = rhs.fgtpi;
336  ftwopi = rhs.ftwopi;
337  fshape = rhs.fshape;
338  fsurf = rhs.fsurf;
339  return *this;
340 }
341 
342 
343 // Dispatch to parameterisation for replication mechanism dimension
344 // computation & modification.
345 void BelleLathe::ComputeDimensions(G4VPVParameterisation*,
346  const G4int,
347  const G4VPhysicalVolume*)
348 {
349  G4Exception("BelleLathe::ComputeDimensions()",
350  "GeomSolids0001", FatalException,
351  "BelleLathe does not support Parameterisation.");
352  // std::cout<<"ComputeDimensions"<<std::endl;
353  // p->ComputeDimensions(*this,n,pRep);
354 }
355 
356 vector<double> quadsolve(double a, double b, double c)
357 {
358  // solve equation a*t^2 + b*t + c = 0 taking care intermediate rounding errors
359  vector<double> t(2);
360  b *= 0.5;
361  double D = b * b - a * c;
362  if (D >= 0) {
363  double sD = sqrt(D);
364  double sum = b + ((b > 0) ? sD : -sD);
365  double t0 = -c / sum;
366  double t1 = -sum / a;
367  t[0] = t0;
368  t[1] = t1;
369  } else {
370  t.clear();
371  }
372 
373  return t;
374 }
375 
376 inline int quadsolve(double a, double b, double c, double& t0, double& t1)
377 {
378  // solve equation a*t^2 + b*t + c = 0 taking care intermediate rounding errors
379  b *= 0.5;
380  double D = b * b - a * c;
381  if (D >= 0) {
382  double sD = sqrt(D);
383  double sum = b + ((b > 0) ? sD : -sD);
384  t0 = -c / sum;
385  t1 = -sum / a;
386  return 2;
387  }
388  return 0;
389 }
390 
392 struct solution_t {
393  double t;
394  double s;
395 };
396 vector<solution_t> extremum(double A, double B, double C, double D, double E, double F)
397 {
398  // extremum of Fun(t,s) = A*t*t + B*t*s + C*s*s + D*t + E*s + F => dFun/ds = 0
399  vector<solution_t> res;
400  if (abs(B) < abs(A)) {
401  double a = 4 * A * C - B * B;
402  double b = 2 * (2 * A * E - B * D);
403  double c = 4 * A * F - D * D;
404  vector<double> ss = quadsolve(a, b, c);
405  for (auto s : ss) {
406  if (fpclassify(s) == FP_INFINITE) continue;
407  double t = -(s * B + D) / (2 * A);
408  solution_t r = {t, s};
409  res.push_back(r);
410  }
411  } else {
412  double B2 = B * B, CD = C * D, BE = B * E;
413  double a = A * (4 * A * C - B2);
414  double b = 2 * A * (2 * CD - BE);
415  double c = D * (CD - BE) + B2 * F;
416  vector<double> ts = quadsolve(a, b, c);
417  for (auto t : ts) {
418  if (fpclassify(t) == FP_INFINITE) continue;
419  double s = -(2 * t * A + D) / B;
420  solution_t r = {t, s};
421  res.push_back(r);
422  }
423  }
424  return res;
425 }
426 
427 // calculate all ray solid's surface intersection return ordered vector
428 vector<double> BelleLathe::linecross(const G4ThreeVector& p, const G4ThreeVector& n) const
429 {
430  auto hitside = [this, &p, &n](double t, double zmin, double zmax) -> bool {
431  double z = p.z() + n.z() * t;
432  bool k = zmin < z && z <= zmax;
433  if (k && !ftwopi)
434  {
435  double x = p.x() + n.x() * t;
436  double y = p.y() + n.y() * t;
437  k = k && insector(x, y);
438  }
439  return k;
440  };
441 
442  auto hitzside = [this, &p, &n](double t, double r2min, double r2max) -> bool {
443  double x = p.x() + n.x() * t;
444  double y = p.y() + n.y() * t;
445  double r2 = x * x + y * y;
446  bool k = r2min <= r2 && r2 < r2max;
447  if (k && !ftwopi)
448  {
449  k = k && insector(x, y);
450  }
451  return k;
452  };
453 
454  vector<double> tc;
455  double inz = 1 / n.z();
456  double nn = Belle2::ECL::dotxy(n, n), np = dotxy(n, p), pp = dotxy(p, p);
457  for (const cachezr_t& s : fcache) { // loop over sides
458  if (s.dz == 0.0) { // z-plane
459  double t = (s.z - p.z()) * inz;
460  if (hitzside(t, s.r2min, s.r2max)) { tc.push_back(t); }
461  } else {
462  double ta = s.ta;
463  double A, B, R2;
464  if (s.dr == 0.0) { // cylinder
465  double R = s.r;
466  R2 = R * R;
467 
468  A = -nn;
469  B = np;
470  } else { // cone
471  double taz = ta * (p.z() - s.z);
472  double R = taz + s.r;
473  R2 = R * R;
474 
475  double nzta = n.z() * ta;
476  A = nzta * nzta - nn;
477  B = np - nzta * R;
478  }
479  double D = B * B + (pp - R2) * A;
480  if (D > 0) {
481  double sD = sqrt(D), iA = 1 / A;
482  double t0 = (B + sD) * iA, t1 = (B - sD) * iA;
483  if (hitside(t0, s.zmin, s.zmax)) tc.push_back(t0);
484  if (hitside(t1, s.zmin, s.zmax)) tc.push_back(t1);
485  }
486  }
487  }
488 
489  if (!ftwopi) {
490  do { // side at phi0
491  double d = fn0x * p.x() + fn0y * p.y();
492  double vn = fn0x * n.x() + fn0y * n.y();
493  double t = -d / vn;
494  G4ThreeVector r = p + n * t;
495  zr_t zr = {r.z(), fc0 * r.x() + fs0 * r.y()};
496  if (vn != 0 && wn_poly(zr) == 2) tc.push_back(t);
497  } while (0);
498 
499  do { // side at phi0+dphi
500  double d = fn1x * p.x() + fn1y * p.y();
501  double vn = fn1x * n.x() + fn1y * n.y();
502  double t = -d / vn;
503  G4ThreeVector r = p + n * t;
504  zr_t zr = {r.z(), fc1 * r.x() + fs1 * r.y()};
505  if (vn != 0 && wn_poly(zr) == 2) tc.push_back(t);
506  } while (0);
507  }
508 
509  sort(tc.begin(), tc.end());
510  return tc;
511 }
512 
513 // Calculate extent under transform and specified limit
514 G4bool BelleLathe::CalculateExtent(const EAxis A,
515  const G4VoxelLimits& bb,
516  const G4AffineTransform& T,
517  G4double& pMin, G4double& pMax) const
518 {
519  auto maxdist = [this](const G4ThreeVector & n) -> G4ThreeVector {
520  G4ThreeVector r;
521  int i = 0, nsize = fcache.size();
522  if (ftwopi || insector(n.x(), n.y())) // n in sector
523  {
524  double nr = hypot(n.x(), n.y()), nz = n.z();
525  double dmax = -kInfinity, R = 0, Z = 0;
526  do {
527  const cachezr_t& s = fcache[i];
528  double d1 = nz * s.z + nr * s.r;
529  if (dmax < d1) { R = s.r; Z = s.z; dmax = d1;}
530  } while (++i < nsize);
531  if (nr > 0) {
532  R /= nr;
533  r.set(R * n.x(), R * n.y(), Z);
534  } else {
535  double phi = fphi + 0.5 * fdphi;
536  r.set(R * cos(phi), R * sin(phi), Z);
537  }
538  } else
539  {
540  double dmax = -kInfinity;
541  do {
542  const cachezr_t& s = fcache[i];
543  // check both sides
544  G4ThreeVector rf(-fn0y * s.r, fn0x * s.r, s.z), rl(fn1y * s.r, -fn1x * s.r, s.z);
545  double d0 = rf * n, d1 = rl * n;
546  // cout<<rf<<" "<<rl<<endl;
547  if (dmax < d0) { r = rf; dmax = d0;}
548  if (dmax < d1) { r = rl; dmax = d1;}
549  } while (++i < nsize);
550  }
551  return r;
552  };
553 
554  struct seg_t {int i0, i1;};
555  auto clip = [](vector<G4ThreeVector>& vlist, vector<seg_t>& slist, const G4ThreeVector & n, double dist) {
556  vector<seg_t> snew;
557  vector<int> lone;
558 
559  vector<double> d;
560  for (const G4ThreeVector& v : vlist) d.push_back(v * n + dist);
561 
562  for (seg_t s : slist) {
563  double prod = d[s.i0] * d[s.i1];
564  // cout<<d[s.i0]<<" "<<d[s.i1]<<endl;
565  if (prod < 0) { // segment crosses plane - break it
566  G4ThreeVector rn = (vlist[s.i0] * d[s.i1] - vlist[s.i1] * d[s.i0]) * (1 / (d[s.i1] - d[s.i0]));
567  lone.push_back(vlist.size()); vlist.push_back(rn);
568  if (d[s.i0] < 0) {
569  s = {lone.back(), s.i1};
570  } else {
571  s = {s.i0, lone.back()};
572  }
573  } else if (prod == 0) { // segment end on plane
574  if (d[s.i0] == 0 && d[s.i1] > 0) {
575  lone.push_back(s.i0);
576  } else if (d[s.i0] > 0 && d[s.i1] == 0) {
577  lone.push_back(s.i1);
578  } else continue;
579  } else {
580  if (d[s.i0] < 0) continue; // segment below plane
581  }
582  snew.push_back(s);
583  }
584 
585  double dmax = -1e99;
586  int imax = -1, jmax = -1;
587  // search for the most distant points on the clipping plane
588  for (unsigned int i = 0; i < lone.size(); i++) {
589  for (unsigned int j = i + 1; j < lone.size(); j++) {
590  double d2 = (vlist[lone[i]] - vlist[lone[j]]).mag2();
591  if (d2 > dmax) { imax = lone[i]; jmax = lone[j];}
592  }
593  }
594 
595  // close the new polygon by creating new segments
596  if (imax >= 0) {
597  G4ThreeVector k = vlist[jmax] - vlist[imax];
598  sort(lone.begin(), lone.end(), [&k, &vlist, &imax](int i, int j) {return k * (vlist[i] - vlist[imax]) < k * (vlist[j] - vlist[imax]);});
599 
600  for (unsigned int i = 0; i < lone.size(); i += 2) {
601  seg_t t = {lone[i], lone[i + 1]};
602  for (const seg_t& s : snew) {
603  if (t.i1 == s.i0) { snew.push_back(t); break;}
604  if (t.i0 == s.i0) { swap(t.i0, t.i1); snew.push_back(t); break;}
605  }
606  }
607  }
608  swap(slist, snew);
609  };
610 
611  auto PhiCrossN = [this, clip](const vector<Plane_t>& planes) {
612  // unordered clipped phi-sides vertices within
613  // limiting planes
614  vector<G4ThreeVector> vlist; // vertex list
615  vector<seg_t> slist; // segment list
616  vector<G4ThreeVector> res;
617 
618  int nsize = fcache.size();
619  vlist.reserve(nsize);
620  slist.reserve(nsize);
621  for (int iphi = 0; iphi < 2; iphi++) {
622  vlist.clear();
623  slist.clear();
624  // phi-side directional vector is (kx,ky)
625  double kx = iphi ? -fn0y : fn1y, ky = iphi ? fn0x : -fn1x;
626  do {
627  int i = 0;
628  do {
629  const cachezr_t& s = fcache[i];
630  G4ThreeVector r(kx * s.r, ky * s.r, s.z);
631  vlist.push_back(r);
632  seg_t t = {i, i + 1};
633  slist.push_back(t);
634  } while (++i < nsize - 1);
635  const cachezr_t& s = fcache[nsize - 1];
636  G4ThreeVector r(kx * s.r, ky * s.r, s.z);
637  vlist.push_back(r);
638  seg_t t = {nsize - 1, 0};
639  slist.push_back(t);
640  } while (0);
641 
642  // clip phi-side polygon by limiting planes
643  for (const Plane_t& p : planes) {
644  // cout<<p.n<<" "<<p.d<<endl;
645  clip(vlist, slist, p.n, p.d);
646  // for(auto t:vlist) cout<<t<<" "; cout<<endl;
647  }
648  vector<bool> bv(vlist.size(), false);
649 
650  for (vector<seg_t>::const_iterator it = slist.begin(); it != slist.end(); ++it) {
651  bv[(*it).i0] = true;
652  bv[(*it).i1] = true;
653  }
654 
655  for (unsigned int i = 0; i < vlist.size(); i++) {
656  if (!bv[i]) continue;
657  res.push_back(vlist[i]);
658  }
659  }
660  return res;
661  };
662 
663  auto RCross = [this](const G4ThreeVector & op, const G4ThreeVector & k, const G4ThreeVector & u) {
664  // plane with origin at op and normal vector n = [k x u], k and u are orthogonal k*u = 0
665  // plane equation r = t*k + s*u + op
666  vector<solution_t> ts;
667  int nsize = fcache.size();
668  int i = 0;
669  do {
670  const cachezr_t& seg = fcache[i];
671  // r0 -- cone radius at z0, c -- cone axis
672  // cone equation is (r0 + tg * ((r-c0)*c))^2 = (r-c0)^2 - ((r-c0)*c)^2
673  double r0 = seg.r, z0 = seg.z, tg = seg.ta, tg2 = tg * tg;
674  double rtg = r0 * tg;
675 
676  G4ThreeVector o(op.x(), op.y(), op.z() - z0);
677 
678  double ko = k * o, uo = u * o, ck = k.z(), cu = u.z(), co = o.z();
679  double k2 = 1, u2 = 1, o2 = o * o;
680  double ck2 = ck * ck, cu2 = cu * cu, co2 = co * co;
681  double dr2 = r0 * r0 - o2;
682  if (seg.dz != 0.0) {
683  double q0 = 1 + tg2;
684  double q1 = co * q0 + rtg;
685 
686  double F00 = co2 * q0 + 2 * co * rtg + dr2;
687  double F10 = 2 * (ck * q1 - ko);
688  double F20 = ck2 * q0 - k2;
689  double F01 = 2 * (cu * q1 - uo);
690  double F11 = 2 * ck * cu * q0;
691  double F02 = cu2 * q0 - u2;
692 
693  vector<solution_t> res = extremum(F02, F11, F20, F01, F10, F00);
694  for (const solution_t& r : res) {
695  double t = r.s, s = r.t;
696  G4ThreeVector p = t * k + s * u + op;
697  if (seg.zmin < p.z() && p.z() < seg.zmax) {
698  solution_t e = {t, s};
699  if (ftwopi || insector(p.x(), p.y()))
700  ts.push_back(e);
701  }
702  }
703  }
704  double a = -(ck2 * u2 + cu2 * k2);
705  if (a != 0) {
706  if (abs(cu) > abs(ck)) {
707  double b = 2 * (ck * (cu * uo - co * u2) - cu2 * ko);
708  double c = co * (2 * cu * uo - co * u2) + cu2 * dr2;
709  vector<double> tv = quadsolve(a, b, c);
710  for (double t : tv) {
711  double s = -(co + ck * t) / cu;
712  G4ThreeVector p = t * k + s * u + op;
713  if (ftwopi || insector(p.x(), p.y())) {
714  solution_t e = {t, s};
715  ts.push_back(e);
716  }
717  }
718  } else {
719  double b = 2 * (cu * (ck * ko - co * k2) - ck2 * uo);
720  double c = co * (2 * ck * ko - co * k2) + ck2 * dr2;
721  vector<double> sv = quadsolve(a, b, c);
722  for (double s : sv) {
723  double t = -(co + cu * s) / ck;
724  G4ThreeVector p = t * k + s * u + op;
725  if (ftwopi || insector(p.x(), p.y())) {
726  solution_t e = {t, s};
727  ts.push_back(e);
728  }
729  }
730  }
731  }
732  } while (++i < nsize);
733  return ts;
734  };
735 
736  bool b1 = false, b2 = false;
737  G4ThreeVector n0, n1, n2;
738  switch (A) {
739  case kXAxis: n0.set(1, 0, 0); n1.set(0, 1, 0); n2.set(0, 0, 1); b1 = bb.IsYLimited(); b2 = bb.IsZLimited(); break;
740  case kYAxis: n0.set(0, 1, 0); n1.set(1, 0, 0); n2.set(0, 0, 1); b1 = bb.IsXLimited(); b2 = bb.IsZLimited(); break;
741  case kZAxis: n0.set(0, 0, 1); n1.set(1, 0, 0); n2.set(0, 1, 0); b1 = bb.IsXLimited(); b2 = bb.IsYLimited(); break;
742  default: break;
743  }
744 
745  double dmin1 = -kInfinity, dmax1 = kInfinity;
746  if (b1) {
747  switch (A) {
748  case kXAxis: dmin1 = bb.GetMinYExtent(); dmax1 = bb.GetMaxYExtent(); break;
749  case kYAxis: dmin1 = bb.GetMinXExtent(); dmax1 = bb.GetMaxXExtent(); break;
750  case kZAxis: dmin1 = bb.GetMinXExtent(); dmax1 = bb.GetMaxXExtent(); break;
751  default: break;
752  }
753  }
754 
755  double dmin2 = -kInfinity, dmax2 = kInfinity;
756  if (b2) {
757  switch (A) {
758  case kXAxis: dmin2 = bb.GetMinZExtent(); dmax2 = bb.GetMaxZExtent(); break;
759  case kYAxis: dmin2 = bb.GetMinZExtent(); dmax2 = bb.GetMaxZExtent(); break;
760  case kZAxis: dmin2 = bb.GetMinYExtent(); dmax2 = bb.GetMaxYExtent(); break;
761  default: break;
762  }
763  }
764 
765  G4AffineTransform iT = T.Inverse();
766  // axis to solid coordinates
767  G4ThreeVector n0t = iT.TransformAxis(n0);
768  G4ThreeVector smin = n0t * kInfinity, smax = (-kInfinity) * n0t; // extremum points in solid coordinate system
769  double pmin = kInfinity, pmax = -pmin;
770  if (b1 && b2) {
771  G4ThreeVector corners[] = {n1* dmin1 + n2 * dmin2, n1* dmax1 + n2 * dmin2, n1* dmax1 + n2 * dmax2, n1* dmin1 + n2 * dmax2};
772  for (G4ThreeVector& c : corners) iT.ApplyPointTransform(c); // to solid coordinates
773 
774  vector<Plane_t> planes;
775  for (int i = 0; i < 4; i++) {
776  const G4ThreeVector& c0 = corners[i], &c1 = corners[(i + 1) % 4];
777  vector<double> dists = linecross(c0, n0t);
778  // cout<<"c0 "<<c0<<endl;
779  for (double t : dists) {
780  G4ThreeVector p = n0t * t + c0;
781  double tt = t + c0 * n0t;
782  // cout<<p<<" "<<tt<<endl;
783  if (pmax < tt) { pmax = tt; smax = p;}
784  if (pmin > tt) { pmin = tt; smin = p;}
785  }
786 
787  G4ThreeVector u = c1 - c0, un = u.unit();
788  vector<solution_t> ts = RCross(c0, n0t, un);
789  double umax = u.mag();
790  for (solution_t r : ts) {
791  if (0 < r.s && r.s < umax) {
792  double tt = r.t + c0 * n0t;
793  G4ThreeVector p = n0t * r.t + un * r.s + c0;
794  // cout<<r.t<<" "<<r.s<<" "<<smax<<endl;
795  if (pmax < tt) { pmax = tt; smax = p;}
796  if (pmin > tt) { pmin = tt; smin = p;}
797  }
798  }
799  planes.push_back({ -un, un * c1});
800  }
801 
802  vector<G4ThreeVector> vside = PhiCrossN(planes);
803  for (G4ThreeVector& p : vside) {
804  // cout<<p<<endl;
805  double tt = n0t * p;
806  if (pmax < tt) { pmax = tt; smax = p;}
807  if (pmin > tt) { pmin = tt; smin = p;}
808  }
809 
810  } else if (b1 || b2) {
811  G4ThreeVector limits[2], u;
812  if (b1) {
813  limits[0] = n1 * dmin1;
814  limits[1] = n1 * dmax1;
815  u = iT.TransformAxis(n2);
816  } else {
817  limits[0] = n2 * dmin2;
818  limits[1] = n2 * dmax2;
819  u = iT.TransformAxis(n1);
820  }
821 
822  for (G4ThreeVector& c : limits) iT.ApplyPointTransform(c); // to solid coordinates
823  for (int i = 0; i < 2; i++) {
824  vector<solution_t> ts = RCross(limits[i], n0t, u);
825  for (solution_t r : ts) {
826  double tt = r.t + limits[i] * n0t;
827  G4ThreeVector p = n0t * r.t + u * r.s + limits[i];
828  // cout<<r.t<<" "<<r.s<<" "<<endl;
829  if (pmax < tt) { pmax = tt; smax = p;}
830  if (pmin > tt) { pmin = tt; smin = p;}
831  }
832  }
833 
834  vector<Plane_t> planes(2);
835  G4ThreeVector n;
836  if (b1) {
837  n = iT.TransformAxis(n1);
838  } else {
839  n = iT.TransformAxis(n2);
840  }
841  planes[0] = { n, -limits[0]* n};
842  planes[1] = { -n, limits[1]* n};
843  vector<G4ThreeVector> vside = PhiCrossN(planes);
844 
845  for (G4ThreeVector& p : vside) {
846  // double t = n0t*(p-limits[0]);
847  double tt = n0t * p;
848  // cout<<tt<<" "<<p<<" "<<endl;
849  if (pmax < tt) { pmax = tt; smax = p;}
850  if (pmin > tt) { pmin = tt; smin = p;}
851  }
852  }
853  // maximal distance in +- directions
854  G4ThreeVector rp = maxdist(n0t), rm = maxdist(-n0t);
855  if (bb.Inside(T.TransformPoint(rm))) {
856  double tt = rm * n0t;
857  if (pmin > tt) {pmin = tt; smin = rm;}
858  }
859  if (bb.Inside(T.TransformPoint(rp))) {
860  double tt = rp * n0t;
861  if (pmax < tt) {pmax = tt; smax = rp;}
862  }
863 
864  // to mother volume coordinate system
865  T.ApplyPointTransform(smin);
866  T.ApplyPointTransform(smax);
867  pmin = n0 * smin;
868  pmax = n0 * smax;
869 
870  pmin -= kCarTolerance;
871  pmax += kCarTolerance;
872  // bool hit = pmin > -kInfinity && pmax < kInfinity;
873  bool hit = pmin < pmax;
874 
875 #if COMPARE==10
876  auto surfhit = [this, &bb, &T, &n0, &n0t](double & pmin, double & pmax, bool print = false)->bool {
877  const int N = 1000 * 1000;
878  if (fsurf.size() == 0) for (int i = 0; i < N; i++) fsurf.push_back(GetPointOnSurface());
879 
880  int umin = -1, umax = -1;
881  double wmin = 1e99, wmax = -1e99;
882  for (int i = 0; i < N; i++)
883  {
884  if (bb.Inside(T.TransformPoint(fsurf[i]))) {
885  double w = n0t * fsurf[i];
886  if (wmin > w) {wmin = w; umin = i;}
887  if (wmax < w) {wmax = w; umax = i;}
888  }
889  }
890  if (print)cout << umin << " " << umax << " " << wmin << " " << wmax << endl;
891  if (umin >= 0 && umax >= 0)
892  {
893  G4ThreeVector qmin = fsurf[umin], qmax = fsurf[umax];
894  T.ApplyPointTransform(qmin);
895  T.ApplyPointTransform(qmax);
896  pmin = n0 * qmin, pmax = n0 * qmax;
897  return true;
898  }
899  return false;
900  };
901 
902  bool res = fshape->CalculateExtent(A, bb, T, pMin, pMax);
903  double srfmin = kInfinity, srfmax = -srfmin;
904  bool sHit = surfhit(srfmin, srfmax);
905  double diff = kCarTolerance;
906  diff = 10;
907  // if (abs(pmin - pMin) > diff || abs(pmax - pMax) > diff || hit != res) {
908  if ((abs(pmin - srfmin) > diff || abs(pmax - srfmax) > diff) && sHit) {
909  cout << "===================================\n";
910  cout << GetName() << " " << fcache.size() << " " << fphi << " " << fdphi << " " << ftwopi << "\n";
911  cout << hit << " " << res << " " << b1 << " " << b2 << "\n";
912  if (sHit) {
913  cout << "ss " << srfmin << " " << srfmax << "\n";
914  } else {
915  cout << "ss : not in bounding box" << "\n";
916  }
917  cout << "my " << pmin << " " << pmax << "\n";
918  cout << "tc " << pMin << " " << pMax << "\n";
919  cout << "df " << pmin - pMin << " " << pmax - pMax << "\n";
920  G4ThreeVector bmin(bb.GetMinXExtent(), bb.GetMinYExtent(), bb.GetMinZExtent());
921  G4ThreeVector bmax(bb.GetMaxXExtent(), bb.GetMaxYExtent(), bb.GetMaxZExtent());
922  cout << "Axis=" << A << " " << bmin << " " << bmax << " " << T << "\n";
923  cout << rp << " " << rm << "\n";
924  cout << smin << " " << smax << "\n";
925  cout << flush;
926  // _exit(0);
927  }
928  // cout<<endl;
929 #endif
930  pMin = pmin;
931  pMax = pmax;
932 
933  return hit;
934 }
935 
936 // True if (x,y) is within the shape rotation
937 inline bool BelleLathe::insector(double x, double y) const
938 {
939  double d0 = fn0x * x + fn0y * y;
940  double d1 = fn1x * x + fn1y * y;
941  bool b0 = d0 < 0, b1 = d1 < 0;
942  return fgtpi ? b0 || b1 : b0 && b1;
943 }
944 
945 int BelleLathe::wn_poly(const zr_t& r) const
946 {
947  int wn = 0;
948  vector<double>::const_iterator it = upper_bound(fz.begin(), fz.end(), r.z);
949  // cout<<r<<" "<<fz.size()<<" "<<it-fz.begin()<<endl;
950  if (it != fz.begin() && it != fz.end()) {
951  int k = it - fz.begin();
952  for (int i = findx[k - 1]; i != findx[k]; i++) {
953  const cachezr_t& s = fcache[fseg[i]];
954  double dz = r.z - s.z, dr = r.r - s.r;
955  double crs = s.dr * dz - s.dz * dr;
956  wn -= (crs > 0) - (crs < 0);
957  }
958  }
959  return wn;
960 }
961 
962 double BelleLathe::mindist(const zr_t& r) const
963 {
964  double d = kInfinity;
965  int i = 0, n = fcache.size();
966  do {
967  const cachezr_t& s = fcache[i];
968  double dz = r.z - s.z, dr = r.r - s.r;
969  double dot = s.dz * dz + s.dr * dr; // projection of the point on the segement
970  if (dot < 0) {
971  d = min(d, dz * dz + dr * dr); // distance to the first point of the segement
972  } else if (dot <= s.s2) { // point should be within the segment
973  double crs = s.dr * dz - s.dz * dr;
974  d = min(d, crs * crs * s.is2);
975  }
976  } while (++i < n);
977  d = sqrt(d);
978  d = (wn_poly(r) == 2) ? -d : d;
979  return d;
980 }
981 
982 // Return whether point inside/outside/on surface, using tolerance
983 EInside BelleLathe::Inside(const G4ThreeVector& p) const
984 {
985  COUNTER(0);
986  const double delta = 0.5 * kCarTolerance;
987  EInside res = kInside;
988  if (!ftwopi) {
989  double d0 = fn0x * p.x() + fn0y * p.y();
990  double d1 = fn1x * p.x() + fn1y * p.y();
991  if (fgtpi) {
992  if (d0 > delta && d1 > delta) { res = kOutside;}
993  else if (d0 > -delta && d1 > -delta) { res = kSurface;}
994  } else {
995  if (d0 > delta || d1 > delta) { res = kOutside;}
996  else if (d0 > -delta || d1 > -delta) { res = kSurface;}
997  }
998  }
999  if (res != kOutside) {
1000  zr_t r = {p.z(), p.perp()};
1001  double d = mindist(r);
1002  if (res == kSurface && d < delta) res = kSurface;
1003  else if (d > delta) res = kOutside;
1004  else if (d > -delta) res = kSurface;
1005  else res = kInside;
1006  }
1007 
1008 #if COMPARE==1
1009  EInside dd = fshape->Inside(p);
1010  if (1 || dd != res) {
1011  double d0 = fn0x * p.x() + fn0y * p.y();
1012  double d1 = fn1x * p.x() + fn1y * p.y();
1013  // if (abs(d0) > kCarTolerance && abs(d1) > kCarTolerance) {
1014  int oldprec = cout.precision(16);
1015  zr_t r = {p.z(), p.perp()};
1016  cout << GetName() << " Inside(p) " << p << " " << r << " my=" << res << " tc=" << dd <<
1017  " dist=" << mindist(r) << " " << d0 << " " << d1 << endl;
1018  cout.precision(oldprec);
1019  // }
1020  }
1021 #endif
1022  MATCHOUT("BelleLathe::Inside(p) " << p << " res= " << res);
1023  return res;
1024 }
1025 
1026 zr_t BelleLathe::normal(const zr_t& r, double& d2) const
1027 {
1028  double d = std::numeric_limits<double>::infinity(), t = 0;
1029  int iseg = -1;
1030  for (int i = 0, imax = fcache.size(); i < imax; i++) {
1031  const cachezr_t& s = fcache[i];
1032  double dz = r.z - s.z, dr = r.r - s.r;
1033  double dot = s.dz * dz + s.dr * dr; // projection of the point on the segement
1034  if (dot < 0) {
1035  double dist = dz * dz + dr * dr; // distance to the first point of the segement
1036  if (dist < d) { d = dist; t = dot * s.is2; iseg = i;}
1037  } else if (dot <= s.s2) { // point should be within the segment
1038  double crs = s.dr * dz - s.dz * dr;
1039  double dist = crs * crs * s.is2;
1040  if (dist < d) { d = dist; t = dot * s.is2; iseg = i;}
1041  }
1042  }
1043  d2 = d;
1044 
1045  auto getn = [this](int i)->zr_t{
1046  int imax = fcache.size();
1047  int i0 = i;
1048  if (i == -1) i0 = imax;
1049  const cachezr_t& s = fcache[i0];
1050  double is = sqrt(s.is2);
1051  return {s.dr * is, -s.dz * is};
1052  };
1053  return getn(iseg);
1054 
1055  if (t < 0.0) {
1056  const cachezr_t& s = fcache[iseg];
1057  zr_t dist = {r.z - s.z, r.r - s.r};
1058  double dist2 = dist.z * dist.z + dist.r * dist.r;
1059  if (dist2 > 1e-18) {
1060  double q = 1 / sqrt(dist2);
1061  if (wn_poly(r) == 2) q = -q;
1062  return {dist.z * q, dist.r * q};
1063  } else {
1064  zr_t n = getn(iseg), np = getn(iseg - 1);
1065  n.z += np.z; n.r += np.r;
1066  double n2 = n.z * n.z + n.r * n.r;
1067  double q = 1 / sqrt(n2);
1068  n.z *= q, n.r *= q;
1069  return n;
1070  }
1071  }
1072  return getn(iseg);
1073 }
1074 
1075 // Calculate side nearest to p, and return normal
1076 // If 2+ sides equidistant, first side's normal returned (arbitrarily)
1077 G4ThreeVector BelleLathe::SurfaceNormal(const G4ThreeVector& p) const
1078 {
1079  COUNTER(1);
1080 
1081  auto side = [this](const zr_t & r, double d, int iside) {
1082  double nx = (iside) ? fn1x : fn0x, ny = (iside) ? fn1y : fn0y;
1083  if (wn_poly(r) == 2) return G4ThreeVector(nx, ny, 0);
1084  double cphi = (iside) ? fc1 : fc0, sphi = (iside) ? fs1 : fc0;
1085 
1086  double d2; zr_t n = normal(r, d2);
1087  double x = cphi * n.r, y = sphi * n.r;
1088  double u = sqrt(d2);
1089  d2 += d * d;
1090  G4ThreeVector res;
1091  if (d2 > 0) {
1092  double q = 1 / sqrt(d2);
1093  double cpsi = u * q, spsi = d * q;
1094  res.set(x * cpsi - y * spsi, x * spsi + y * cpsi, n.z);
1095  }
1096  res.set(x, y, n.z);
1097  return res;
1098  };
1099 
1100  G4ThreeVector res;
1101  zr_t r = {p.z(), p.perp()};
1102  double d2; zr_t n = normal(r, d2);
1103  double d = sqrt(d2);
1104  double pt = hypot(p.x(), p.y());
1105 
1106  if (pt > 0) {
1107  double ir = n.r / pt;
1108  res = G4ThreeVector(ir * p.x(), ir * p.y(), n.z);
1109  } else
1110  res = G4ThreeVector(n.r, 0, n.z);
1111 
1112  if (!ftwopi) {
1113  double d0 = fn0x * p.x() + fn0y * p.y();
1114  double d1 = fn1x * p.x() + fn1y * p.y();
1115  zr_t r0 = {p.z(), fc0 * p.x() + fs0 * p.y()}; // projection on plane phi
1116  zr_t r1 = {p.z(), fc1 * p.x() + fs1 * p.y()}; // projection on plane phi+dphi
1117  if (fgtpi) {
1118  if (d0 > 0 && d1 > 0) { // outside sector
1119  if (d0 < d1) {
1120  res = side(r0, d0, 0); goto exit;
1121  } else {
1122  res = side(r1, -d1, 1); goto exit;
1123  }
1124  } else {// inside sector
1125  if (wn_poly(r) == 2) { // point p inside the solid
1126  if (abs(d0) < d && abs(d0) < abs(d1)) { res = G4ThreeVector(fn0x, fn0y, 0); goto exit;}
1127  if (abs(d1) < d && abs(d1) < abs(d0)) { res = G4ThreeVector(fn1x, fn1y, 0); goto exit;}
1128  }
1129  }
1130  } else {
1131  if (d0 > 0 || d1 > 0) { // outside sector
1132  if (d0 < 0) {
1133  res = side(r1, -d1, 1); goto exit;
1134  } else {
1135  res = side(r0, d0, 0); goto exit;
1136  }
1137  } else {
1138  if (wn_poly(r) == 2) { // point p inside the solid
1139  if (abs(d0) < d && abs(d0) < abs(d1)) { res = G4ThreeVector(fn0x, fn0y, 0); goto exit;}
1140  if (abs(d1) < d && abs(d1) < abs(d0)) { res = G4ThreeVector(fn1x, fn1y, 0); goto exit;}
1141  }
1142  }
1143  }
1144  }
1145 exit:
1146 #if COMPARE==1
1147  G4ThreeVector dd = fshape->SurfaceNormal(p);
1148  if ((res - dd).mag() > 1e-11) {
1149  int oldprec = cout.precision(16);
1150  EInside inside = fshape->Inside(p);
1151  cout << GetName() << " SurfaceNormal(p) " << p << " " << res << " " << dd << " " << res - dd << " " << inside << endl;
1152  cout.precision(oldprec);
1153  // _exit(0);
1154  }
1155 #endif
1156  MATCHOUT("BelleLathe::SurfaceNormal(p,n) " << p << " res= " << res);
1157  return res;
1158 }
1159 
1160 // Calculate exact shortest distance to any boundary from outside
1161 // This is the best fast estimation of the shortest distance to trap
1162 // - Returns 0 is ThreeVector inside
1163 G4double BelleLathe::DistanceToIn(const G4ThreeVector& p) const
1164 {
1165  COUNTER(2);
1166  double d = 0;
1167  // int sector = 0, plane = 0;
1168  if (ftwopi) {
1169  zr_t r = {p.z(), p.perp()};
1170  d = max(mindist(r), 0.0);
1171  } else {
1172  double d0 = fn0x * p.x() + fn0y * p.y();
1173  double d1 = fn1x * p.x() + fn1y * p.y();
1174 
1175  if (fgtpi) {
1176  if (d0 > 0 && d1 > 0) { // outside sector
1177  if (d0 < d1) {
1178  zr_t r = {p.z(), -fn0y * p.x() + fn0x * p.y()}; // projection on plane
1179  d = sqrt(pow(max(mindist(r), 0.0), 2) + d0 * d0);
1180  } else {
1181  zr_t r = {p.z(), fn1y * p.x() - fn1x * p.y()}; // projection on plane
1182  d = sqrt(pow(max(mindist(r), 0.0), 2) + d1 * d1);
1183  }
1184  } else {
1185  zr_t r = {p.z(), p.perp()};
1186  d = max(mindist(r), 0.0);
1187  }
1188  } else {
1189  if (d0 > 0 || d1 > 0) { // outside sector
1190  if (d0 < 0) {
1191  zr_t r = {p.z(), fn1y * p.x() - fn1x * p.y()}; // projection on plane
1192  d = sqrt(pow(max(mindist(r), 0.0), 2) + d1 * d1);
1193  } else {
1194  zr_t r = {p.z(), -fn0y * p.x() + fn0x * p.y()}; // projection on plane
1195  d = sqrt(pow(max(mindist(r), 0.0), 2) + d0 * d0);
1196  }
1197  } else {
1198  zr_t r = {p.z(), p.perp()};
1199  d = max(mindist(r), 0.0);
1200  }
1201  }
1202  }
1203 #if COMPARE==1
1204  // double dd = fshape->Inside(p) == 2 ? 0.0 : fshape->DistanceToIn(p);
1205  double dd = fshape->DistanceToIn(p);
1206  // if (abs(d - dd) > kCarTolerance) {
1207  if (dd > d && abs(d - dd) > kCarTolerance) {
1208  int oldprec = cout.precision(16);
1209  EInside inside = fshape->Inside(p);
1210  zr_t r = {p.z(), p.perp()};
1211  cout << GetName() << " DistanceToIn(p) " << p << " " << r << " " << d << " " << dd << " " << d - dd << " " << inside << endl;
1212  cout.precision(oldprec);
1213  // exit(0);
1214  }
1215 #endif
1216  MATCHOUT("BelleLathe::DistanceToIn(p) " << p << " res= " << d);
1217  return d;
1218 }
1219 
1220 // Calculate exact shortest distance to any boundary from inside
1221 // - Returns 0 is ThreeVector outside
1222 G4double BelleLathe::DistanceToOut(const G4ThreeVector& p) const
1223 {
1224  // return ref->DistanceToOut(p);
1225  COUNTER(3);
1226  zr_t r = {p.z(), p.perp()};
1227  double d = mindist(r);
1228  if (!ftwopi) {
1229  double d0 = fn0x * p.x() + fn0y * p.y();
1230  double d1 = fn1x * p.x() + fn1y * p.y();
1231  if (fgtpi) {
1232  d = max(d, min(d0, d1));
1233  } else {
1234  d = max(d, max(d0, d1));
1235  }
1236  }
1237  d = max(-d, 0.0);
1238 #if COMPARE==1
1239  double dd = fshape->Inside(p) == 0 ? 0.0 : fshape->DistanceToOut(p);
1240  if (abs(d - dd) > kCarTolerance) {
1241  int oldprec = cout.precision(16);
1242  zr_t r = {p.z(), p.perp()};
1243  // cout<<r<<endl;
1244  cout << GetName() << " DistanceToOut(p) " << p << " " << r << " " << d << " " << dd << " " << d - dd << endl;
1245  cout.precision(oldprec);
1246  }
1247 #endif
1248  MATCHOUT("BelleLathe::DistanceToOut(p) " << p.x() << " " << p.y() << " " << p.z() << " res= " << d);
1249  return d;
1250 }
1251 
1252 // Calculate distance to shape from outside - return kInfinity if no intersection
1253 G4double BelleLathe::DistanceToIn(const G4ThreeVector& p, const G4ThreeVector& n) const
1254 {
1255  // return fshape->DistanceToIn(p, n);
1256  auto getnormal = [this, &p, &n](int i, double t) ->G4ThreeVector{
1257  const int imax = fcache.size();
1258  G4ThreeVector o;
1259  if (i < 0)
1260  {
1261  } else if (i < imax)
1262  {
1263  const cachezr_t& s = fcache[i];
1264  if (s.dz == 0.0) {
1265  o.setZ(copysign(1, s.dr));
1266  } else {
1267  double x = p.x() + n.x() * t;
1268  double y = p.y() + n.y() * t;
1269  double sth = s.dr * s.is, cth = -s.dz * s.is;
1270  double ir = cth / sqrt(x * x + y * y);
1271  o.set(x * ir, y * ir, sth);
1272  }
1273  } else if (i == imax)
1274  {
1275  o.setX(fn0x), o.setY(fn0y);
1276  } else
1277  {
1278  o.setX(fn1x), o.setY(fn1y);
1279  }
1280  return o;
1281  };
1282 
1283  auto hitside = [this, &p, &n](double t, const cachezr_t& s) -> bool {
1284  double z = p.z() + n.z() * t;
1285  // cout<<t<<" "<<x<<" "<<y<<" "<<z<<endl;
1286  bool k = s.zmin < z && z <= s.zmax;
1287  if (k && !ftwopi)
1288  {
1289  double x = p.x() + n.x() * t;
1290  double y = p.y() + n.y() * t;
1291  k = k && insector(x, y);
1292  }
1293  return k;
1294  };
1295 
1296  auto hitzside = [this, &p, &n](double t, const cachezr_t& s) -> bool {
1297  double x = p.x() + n.x() * t;
1298  double y = p.y() + n.y() * t;
1299  double r2 = x * x + y * y;
1300  bool k = s.r2min <= r2 && r2 < s.r2max;
1301  if (k && !ftwopi)
1302  {
1303  k = k && insector(x, y);
1304  }
1305  return k;
1306  };
1307 
1308  auto hitphi0side = [this, &p, &n](double t) -> bool {
1309  double x = p.x() + n.x() * t;
1310  double y = p.y() + n.y() * t;
1311  double r = x * fc0 + y * fs0;
1312  if (r >= frmin)
1313  {
1314  double z = p.z() + n.z() * t;
1315  zr_t zr = {z, r};
1316  return wn_poly(zr) == 2;
1317  }
1318  return false;
1319  };
1320 
1321  auto hitphi1side = [this, &p, &n](double t) -> bool {
1322  double x = p.x() + n.x() * t;
1323  double y = p.y() + n.y() * t;
1324  double r = x * fc1 + y * fs1;
1325  if (r >= frmin)
1326  {
1327  double z = p.z() + n.z() * t;
1328  zr_t zr = {z, r};
1329  return wn_poly(zr) == 2;
1330  }
1331  return false;
1332  };
1333 
1334  double tmin = kInfinity;
1335  const int imax = fcache.size();
1336  int iseg = -1, isurface = -1;
1337  const double delta = 0.5 * kCarTolerance;
1338  double inz = 1 / n.z();
1339  double nn = dotxy(n, n), np = dotxy(n, p), pp = dotxy(p, p);
1340  double pz = p.z(), pr = sqrt(pp);
1341  for (int i = 0; i < imax; i++) { // loop over sides
1342  const cachezr_t& s = fcache[i];
1343  double dz = pz - s.z, dr = pr - s.r;
1344  double d = dz * s.dr - dr * s.dz;
1345  bool surface = false;
1346  if (abs(d * s.is) < delta) {
1347  double dot = dz * s.dz + dr * s.dr;
1348  if (dot >= 0 && dot <= s.s2) {
1349  surface = true;
1350  isurface = i;
1351  }
1352  }
1353  if (s.dz == 0.0) { // z-plane
1354  if (!surface) {
1355  double t = -dz * inz;
1356  if (0 < t && t < tmin && hitzside(t, s)) { tmin = t; iseg = i;}
1357  }
1358  } else {
1359  double A, B, R2;
1360  if (s.dr == 0.0) { // cylinder
1361  double R = s.r;
1362  R2 = R * R;
1363 
1364  A = -nn;
1365  B = np;
1366  } else { // cone
1367  double taz = s.ta * dz;
1368  double R = taz + s.r;
1369  R2 = R * R;
1370 
1371  double nzta = n.z() * s.ta;
1372  A = nzta * nzta - nn;
1373  B = np - nzta * R;
1374  }
1375  double C = pp - R2;
1376  double D = B * B + C * A;
1377  if (D > 0) {
1378  // double sD = sqrt(D), iA = 1 / A;
1379  // double t0 = (B + sD) * iA, t1 = (B - sD) * iA;
1380  double sD = sqrt(D), sum = B + copysign(sD, B);
1381  double t0 = -C / sum, t1 = sum / A;
1382  if (surface) { // exclude solution on surface
1383  if (abs(t0) > abs(t1)) {
1384  if (t0 > 0 && t0 < tmin && hitside(t0, s)) { tmin = t0; iseg = i;}
1385  } else {
1386  if (t1 > 0 && t1 < tmin && hitside(t1, s)) { tmin = t1; iseg = i;}
1387  }
1388  } else {
1389  if (t0 > 0 && t0 < tmin && hitside(t0, s)) { tmin = t0; iseg = i;}
1390  if (t1 > 0 && t1 < tmin && hitside(t1, s)) { tmin = t1; iseg = i;}
1391  }
1392  }
1393  }
1394  }
1395 
1396  if (!ftwopi) {
1397  do { // side at phi0
1398  double vn = fn0x * n.x() + fn0y * n.y();
1399  if (vn < 0) {
1400  double d = fn0x * p.x() + fn0y * p.y();
1401  double t = -d / vn;
1402  if (hitphi0side(t)) {
1403  bool surface = std::abs(d) < delta;
1404  if (surface) {
1405  tmin = 0; iseg = imax + 0;
1406  } else {
1407  if (0 < t && t < tmin) {tmin = t; iseg = imax + 0;}
1408 
1409  }
1410  }
1411  }
1412  } while (0);
1413 
1414  do { // side at phi0+dphi
1415  double vn = fn1x * n.x() + fn1y * n.y();
1416  if (vn < 0) {
1417  double d = fn1x * p.x() + fn1y * p.y();
1418  double t = -d / vn;
1419  if (hitphi1side(t)) {
1420  bool surface = std::abs(d) < delta;
1421  if (surface) {
1422  tmin = 0; iseg = imax + 1;
1423  } else {
1424  if (0 < t && t < tmin) { tmin = t; iseg = imax + 1;}
1425  }
1426  }
1427  }
1428  } while (0);
1429  }
1430 
1431  if (iseg >= 0) {
1432  if (getnormal(iseg, tmin)*n > 0) tmin = 0;
1433  // if (getnormal(iseg, tmin)*n > 0) tmin = kInfinity; // mimic genericpolycone
1434  }
1435 
1436  auto convex = [this, imax](int i) -> bool{
1437  if (i < imax)
1438  return fcache[i].isconvex;
1439  else
1440  return !fgtpi;
1441  };
1442 
1443  if (tmin >= 0 && tmin < kInfinity) {
1444  if (isurface >= 0) if (convex(isurface) && getnormal(isurface, 0)*n >= 0) tmin = kInfinity;
1445  } else {
1446  if (Inside(p) == kSurface) {
1447  if (isurface >= 0) {
1448  tmin = (getnormal(isurface, 0) * n >= 0) ? kInfinity : 0; // mimic genericpolycone
1449  }
1450  }
1451  }
1452 
1453 #if COMPARE==1
1454  // double dd = fshape->Inside(p) == 2 ? 0.0 : fshape->DistanceToIn(p, n);
1455  double dd = fshape->DistanceToIn(p, n);
1456  if (abs(tmin - dd) > 1e-10) {
1457  int oldprec = cout.precision(16);
1458  EInside inside = fshape->Inside(p);
1459  cout << GetName() << " DistanceToIn(p,v) " << p << " " << n << " " << tmin << " " << dd << " " << tmin - dd << " " << inside << " "
1460  << Inside(p) << " iseg = " << iseg << " " << isurface << endl;
1461  if (isurface >= 0) cout << getnormal(isurface, 0) << endl;
1462  cout.precision(oldprec);
1463  }
1464 #endif
1465  tmin = max(0.0, tmin);
1466  MATCHOUT("BelleLathe::DistanceToIn(p,n) " << p << " " << n << " res= " << tmin);
1467  return tmin;
1468 }
1469 
1470 // Calculate distance to surface of shape from inside
1471 G4double BelleLathe::DistanceToOut(const G4ThreeVector& p, const G4ThreeVector& n,
1472  const G4bool calcNorm, G4bool* IsValid, G4ThreeVector* _n) const
1473 {
1474  // return fshape->DistanceToOut(p, n, calcNorm, IsValid, _n);
1475  auto getnormal = [this, &p, &n](int i, double t)->G4ThreeVector{
1476  const int imax = fcache.size();
1477  G4ThreeVector o;
1478  if (i < 0)
1479  {
1480  } else if (i < imax)
1481  {
1482  const cachezr_t& s = fcache[i];
1483  if (s.dz == 0.0) {
1484  o.setZ(copysign(1, s.dr));
1485  } else {
1486  double x = p.x() + n.x() * t;
1487  double y = p.y() + n.y() * t;
1488  double sth = s.dr * s.is, cth = -s.dz * s.is;
1489  double ir = cth / sqrt(x * x + y * y);
1490  o.set(x * ir, y * ir, sth);
1491  }
1492  } else if (i == imax)
1493  {
1494  o.setX(fn0x), o.setY(fn0y);
1495  } else
1496  {
1497  o.setX(fn1x), o.setY(fn1y);
1498  }
1499  return o;
1500  };
1501 
1502  double nn = dotxy(n, n), np = dotxy(n, p), pp = dotxy(p, p);
1503  auto hitside = [this, &p, &n, nn, np, pp](double t, const cachezr_t& s) -> bool {
1504  double z = p.z() + n.z() * t;
1505  // cout<<t<<" "<<x<<" "<<y<<" "<<z<<endl;
1506  double dot = n.z() * s.dr * sqrt(pp + ((np + np) + nn * t) * t) - s.dz * (np + nn * t);
1507  bool k = s.zmin < z && z <= s.zmax && dot > 0;
1508  if (k && !ftwopi)
1509  {
1510  double x = p.x() + n.x() * t;
1511  double y = p.y() + n.y() * t;
1512 
1513  k = k && insector(x, y);
1514  }
1515  return k;
1516  };
1517 
1518  auto hitzside = [this, &p, &n](double t, const cachezr_t& s) -> bool {
1519  double x = p.x() + n.x() * t;
1520  double y = p.y() + n.y() * t;
1521  double r2 = x * x + y * y;
1522  bool k = s.dr * n.z() > 0 && s.r2min <= r2 && r2 < s.r2max;
1523  if (k && !ftwopi)
1524  {
1525  k = k && insector(x, y);
1526  }
1527  return k;
1528  };
1529 
1530  auto hitphi0side = [this, &p, &n](double t) -> bool {
1531  double x = p.x() + n.x() * t;
1532  double y = p.y() + n.y() * t;
1533  double r = x * fc0 + y * fs0;
1534  if (r >= frmin)
1535  {
1536  double z = p.z() + n.z() * t;
1537  zr_t zr = {z, r};
1538  return wn_poly(zr) == 2;
1539  }
1540  return false;
1541  };
1542 
1543  auto hitphi1side = [this, &p, &n](double t) -> bool {
1544  double x = p.x() + n.x() * t;
1545  double y = p.y() + n.y() * t;
1546  double r = x * fc1 + y * fs1;
1547  if (r >= frmin)
1548  {
1549  double z = p.z() + n.z() * t;
1550  zr_t zr = {z, r};
1551  return wn_poly(zr) == 2;
1552  }
1553  return false;
1554  };
1555 
1556  COUNTER(5);
1557  double tmin = kInfinity;
1558 
1559  const int imax = fcache.size();
1560  int iseg = -1, isurface = -1;
1561 
1562  const double delta = 0.5 * kCarTolerance;
1563  double inz = 1 / n.z();
1564  double pz = p.z(), pr = sqrt(pp);
1565 
1566  for (int i = 0; i < imax; i++) {
1567  const cachezr_t& s = fcache[i];
1568 
1569  double d = (pz - s.z) * s.dr - (pr - s.r) * s.dz;
1570  bool surface = abs(d * s.is) < delta;
1571  if (surface) isurface = i;
1572  if (s.dz == 0.0) {
1573  double t = (s.z - p.z()) * inz;
1574  if (surface) {
1575  if (hitzside(t, s)) {tmin = 0; iseg = i; break;}
1576  } else {
1577  if (0 < t && t < tmin && hitzside(t, s)) {tmin = t; iseg = i;}
1578  }
1579  } else {
1580  double A, B, R2;
1581  if (s.dr == 0.0) { // cylinder
1582  double R = s.r;
1583  R2 = R * R;
1584 
1585  A = -nn;
1586  B = np;
1587  } else { // cone
1588  double taz = s.ta * (p.z() - s.z);
1589  double R = taz + s.r;
1590  R2 = R * R;
1591 
1592  double nzta = n.z() * s.ta;
1593  A = nzta * nzta - nn;
1594  B = np - nzta * R;
1595  }
1596  double C = pp - R2;
1597  double D = B * B + C * A;
1598  if (D > 0) {
1599  double sD = sqrt(D);
1600  double sum = B + copysign(sD, B);
1601  double t0 = -C / sum, t1 = sum / A;
1602  // cout<<t0<<" "<<t1<<" "<<endl;
1603  if (surface) { //exclude solution on surface
1604  if (abs(t0) < abs(t1)) {
1605  if (hitside(t0, s)) { tmin = 0; iseg = i; break;}
1606  if (0 < t1 && t1 < tmin && hitside(t1, s)) { tmin = t1; iseg = i;}
1607  } else {
1608  if (hitside(t1, s)) { tmin = 0; iseg = i; break;}
1609  if (0 < t0 && t0 < tmin && hitside(t0, s)) { tmin = t0; iseg = i;}
1610  }
1611  } else { // check both solutions
1612  if (0 < t0 && t0 < tmin && hitside(t0, s)) { tmin = t0; iseg = i;}
1613  if (0 < t1 && t1 < tmin && hitside(t1, s)) { tmin = t1; iseg = i;}
1614  }
1615  }
1616  }
1617  // cout<<i<<" "<<iseg<<" "<<sqrt(d*d*s.is2)<<" "<<tmin<<" "<<s.zmin<<" "<<s.zmax<<endl;
1618  }
1619 
1620  if (!ftwopi) {
1621  do { // side at phi0
1622  double vn = fn0x * n.x() + fn0y * n.y();
1623  if (vn > 0) {
1624  double d = fn0x * p.x() + fn0y * p.y();
1625  double t = -d / vn;
1626  if (hitphi0side(t)) {
1627  bool surface = std::abs(d) < delta;
1628  if (surface) {
1629  tmin = 0; iseg = imax + 0;
1630  } else {
1631  if (0 < t && t < tmin) {tmin = t; iseg = imax + 0;}
1632 
1633  }
1634  }
1635  }
1636  } while (0);
1637 
1638  do { // side at phi0+dphi
1639  double vn = fn1x * n.x() + fn1y * n.y();
1640  if (vn > 0) {
1641  double d = fn1x * p.x() + fn1y * p.y();
1642  double t = -d / vn;
1643  if (hitphi1side(t)) {
1644  bool surface = std::abs(d) < delta;
1645  if (surface) {
1646  tmin = 0; iseg = imax + 1;
1647  } else {
1648  if (0 < t && t < tmin) { tmin = t; iseg = imax + 1;}
1649  }
1650  }
1651  }
1652  } while (0);
1653  }
1654 
1655  auto convex = [this, imax](int i) -> bool{
1656  if (i < imax)
1657  return fcache[i].isconvex;
1658  else
1659  return !fgtpi;
1660  };
1661 
1662  if (calcNorm) {
1663  if (tmin >= 0 && tmin < kInfinity) {
1664  *_n = getnormal(iseg, tmin);
1665  *IsValid = convex(iseg);
1666  } else {
1667  if (Inside(p) == kSurface) {
1668  if (isurface >= 0) {
1669  *_n = getnormal(isurface, tmin);
1670  *IsValid = convex(isurface);
1671  tmin = 0;
1672  }
1673  } else {
1674  *IsValid = false;
1675  }
1676  }
1677  }
1678  // cout<<"tmin "<<tmin<<endl;
1679 #if COMPARE==1
1680  do {
1681  // G4ThreeVector p0(1210.555046, -292.9578965, -36.71671492);
1682  // if ((p - p0).mag() > 1e-2)continue;
1683  bool isvalid;
1684  G4ThreeVector norm;
1685  double dd = fshape->DistanceToOut(p, n, calcNorm, &isvalid, &norm);
1686  if (abs(tmin - dd) > 1e-10 || (calcNorm && *IsValid != isvalid)) {
1687  int oldprec = cout.precision(16);
1688  cout << GetName() << " DistanceToOut(p,v) p,n =" << curl_t(p) << curl_t(n) << " calcNorm=" << calcNorm
1689  << " myInside=" << Inside(p) << " tmin=" << tmin << " dd=" << dd << " d=" << tmin - dd << " iseg=" << iseg << " isurf=" << isurface
1690  << " ";
1691  if (calcNorm) cout << "myIsValid = " << *IsValid << " tIsValid=" << isvalid << " myn=" << (*_n) << " tn=" << (norm);
1692  cout << endl;
1693  cout.precision(oldprec);
1694  // _exit(0);
1695  }
1696  } while (0);
1697 #endif
1698  MATCHOUT("BelleLathe::DistanceToOut(p,n) " << p << " " << n << " res= " << tmin);
1699  return tmin;
1700 }
1701 
1703 {
1704  ftlist.clear();
1705  unsigned int n = fcontour.size();
1706  vector<int> indx;
1707  for (unsigned int i = 0; i < n; i++) indx.push_back(i);
1708  int count = 0;
1709  while (indx.size() > 3 && ++count < 200) {
1710  unsigned int ni = indx.size();
1711  for (unsigned int i = 0; i < ni; i++) {
1712  int i0 = indx[i], i1 = indx[(i + 1) % ni], i2 = indx[(i + 2) % ni];
1713  double nx = fcontour[i2].z - fcontour[i0].z;
1714  double ny = fcontour[i2].r - fcontour[i0].r;
1715  double d1 = nx * (fcontour[i1].r - fcontour[i0].r) - ny * (fcontour[i1].z - fcontour[i0].z);
1716  bool ear = true;
1717  for (unsigned int j = 0; j < ni - 3; j++) {
1718  int k = indx[(i + 3 + j) % ni];
1719  double d = nx * (fcontour[k].r - fcontour[i0].r) - ny * (fcontour[k].z - fcontour[i0].z);
1720  if (d * d1 > 0) {
1721  ear = false;
1722  break;
1723  }
1724  }
1725  if (ear) {
1726  triangle_t t = {i0, i1, i2};
1727  ftlist.push_back(t);
1728  indx.erase(indx.begin() + (i + 1) % ni);
1729  break;
1730  }
1731  }
1732  }
1733  if (indx.size() == 3) {
1734  triangle_t t = {indx[0], indx[1], indx[2]};
1735  ftlist.push_back(t);
1736  }
1737 }
1738 
1740 {
1741  double vol = 0;
1742  for (const cachezr_t& s : fcache) vol += s.dz * ((3 * s.r) * (s.r + s.dr) + s.dr * s.dr);
1743  fCubicVolume = -fdphi * vol / 6;
1744 
1745  double totalarea = 0;
1746  if (!ftwopi) {
1747  eartrim();
1748  for (const triangle_t& t : ftlist) {
1749  const zr_t& p0 = fcontour[t.i0], &p1 = fcontour[t.i1], &p2 = fcontour[t.i2];
1750  double area = (p1.z - p0.z) * (p2.r - p0.r) - (p1.r - p0.r) * (p2.z - p0.z);
1751  totalarea += abs(area);
1752  farea.push_back(totalarea);
1753  }
1754  }
1755 
1756  for (const cachezr_t& s : fcache) {
1757  double area = fdphi * (s.r + 0.5 * s.dr) * sqrt(s.dz * s.dz + s.dr * s.dr);
1758  totalarea += area;
1759  farea.push_back(totalarea);
1760  }
1761  fSurfaceArea = farea.back();
1762 }
1763 
1764 G4ThreeVector BelleLathe::GetPointOnSurface() const
1765 {
1766  auto GetPointOnTriangle = [this](const triangle_t& t)-> G4ThreeVector{
1767  // barycentric coordinates
1768  double a1 = CLHEP::RandFlat::shoot(0., 1.), a2 = CLHEP::RandFlat::shoot(0., 1.);
1769  if (a1 + a2 > 1) { a1 = 1 - a1; a2 = 1 - a2;}
1770  double a0 = 1 - (a1 + a2);
1771  const zr_t& p0 = fcontour[t.i0], &p1 = fcontour[t.i1], &p2 = fcontour[t.i2];
1772  zr_t p = {p0.z* a0 + p1.z* a1 + p2.z * a2, p0.r* a0 + p1.r* a1 + p2.r * a2};
1773  double c, s;
1774  if (CLHEP::RandFlat::shoot(0., 1.) > 0.5) // select phi side
1775  {
1776  c = -fn0y; s = fn0x;
1777  } else
1778  {
1779  c = fn1y; s = -fn1x;
1780  }
1781  G4ThreeVector r1(p.r * c, p.r * s, p.z);
1782  return r1;
1783  };
1784 
1785  double rnd = CLHEP::RandFlat::shoot(0., farea.back());
1786  std::vector<double>::const_iterator it = std::lower_bound(farea.begin(), farea.end(), rnd);
1787  unsigned int i = it - farea.begin();
1788 
1789  if (!ftwopi) {
1790  if (i < ftlist.size()) {
1791  return GetPointOnTriangle(ftlist[i]);
1792  } else {
1793  i -= ftlist.size();
1794  }
1795  }
1796 
1797  const cachezr_t& s = fcache[i];
1798  double I = 2 * s.r + s.dr;
1799  double Iw = CLHEP::RandFlat::shoot(0., I);
1800  double q = sqrt(Iw * s.dr + s.r * s.r);
1801  double t = Iw / (q + s.r);
1802  double z = s.z + s.dz * t;
1803  double r = s.r + s.dr * t;
1804  double phi = CLHEP::RandFlat::shoot(fphi, fphi + fdphi);
1805  double x = r * cos(phi), y = r * sin(phi);
1806  return G4ThreeVector(x, y, z);
1807 }
1808 
1809 // GetEntityType
1810 G4GeometryType BelleLathe::GetEntityType() const
1811 {
1812  return G4String("BelleLathe");
1813 }
1814 
1815 // Make a clone of the object
1816 G4VSolid* BelleLathe::Clone() const
1817 {
1818  return new BelleLathe(*this);
1819 }
1820 
1821 // Stream object contents to an output stream
1822 std::ostream& BelleLathe::StreamInfo(std::ostream& os) const
1823 {
1824  G4int oldprc = os.precision(16);
1825  os << "-----------------------------------------------------------\n"
1826  << " *** Dump for solid - " << GetName() << " ***\n"
1827  << " ===================================================\n"
1828  << " Solid type: BelleLathe\n"
1829  << " Contour: " << fcontour.size() << " sides, {z, r} points \n";
1830  for (int i = 0, imax = fcontour.size(); i < imax; i++) {
1831  os << fcontour[i] << ", ";
1832  }
1833  os << "\n";
1834  for (int i = 0, imax = fcontour.size(); i < imax; i++) {
1835  os << fcache[i].isconvex << ", ";
1836  }
1837  os << "\n";
1838  os << "phi0 = " << fphi << ", dphi = " << fdphi << ", Full Circle = " << (ftwopi ? "yes" : "no") << "\n";
1839  double xmin = fzmin - 0.05 * (fzmax - fzmin), xmax = fzmax + 0.05 * (fzmax - fzmin);
1840  double ymin = frmin - 0.05 * (frmax - frmin), ymax = frmax + 0.05 * (frmax - frmin);
1841  os << " BB: " << xmin << ", " << xmax << ", " << ymin << ", " << ymax << endl;
1842  os << "-----------------------------------------------------------\n";
1843  os.precision(oldprc);
1844 
1845  return os;
1846 }
1847 
1848 // Methods for visualisation
1849 void BelleLathe::DescribeYourselfTo(G4VGraphicsScene& scene) const
1850 {
1851  scene.AddSolid(*this);
1852 }
1853 
1854 // Function to define the bounding box
1855 void BelleLathe::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const
1856 {
1857  std::vector<vector_t> points;
1858  vector_t point;
1859 
1860  // Placeholder vectors
1861  const double inf = std::numeric_limits<double>::infinity();
1862  G4ThreeVector minimum(inf, inf, inf), maximum(-inf, -inf, -inf);
1863 
1864  // Outer vertices
1865  if (fdphi < 2 * M_PI) { // Only axis-crossings are relevant if the shape is a full circle
1866  point.x = frmax * cos(fphi); point.y = frmax * sin(fphi);
1867  points.push_back(point);
1868  point.x = frmax * cos(fphi + fdphi); point.y = frmax * sin(fphi + fdphi);
1869  points.push_back(point);
1870  }
1871 
1872  // Inner vertices
1873  point.x = frmin * cos(fphi); point.y = frmin * sin(fphi);
1874  points.push_back(point);
1875  if (frmin != 0) { // Avoid duplicate (0,0)
1876  point.x = frmin * cos(fphi + fdphi); point.y = frmin * sin(fphi + fdphi);
1877  points.push_back(point);
1878  }
1879 
1880  // Check if shape crosses an axis
1881  // Because the inside is concave, extremums will always be at vertices
1882  // Because the outside is convex, extremums may be at an axis-crossing
1883  if (insector(0, 1)) {
1884  point.x = 0; point.y = frmax;
1885  points.push_back(point);
1886  }
1887  if (insector(-1, 0)) {
1888  point.x = -frmax; point.y = 0;
1889  points.push_back(point);
1890  }
1891  if (insector(0, -1)) {
1892  point.x = 0; point.y = -frmax;
1893  points.push_back(point);
1894  }
1895  if (insector(1, 0)) {
1896  point.x = frmax; point.y = 0;
1897  points.push_back(point);
1898  }
1899 
1900  for (std::vector<vector_t>::size_type i = 0; i != points.size(); i++) {
1901  // global min in x
1902  if (points[i].x < minimum.x())
1903  minimum.setX(points[i].x);
1904 
1905  // global min in y
1906  if (points[i].y < minimum.y())
1907  minimum.setY(points[i].y);
1908 
1909  // global max in x
1910  if (points[i].x > maximum.x())
1911  maximum.setX(points[i].x);
1912 
1913  // global max in y
1914  if (points[i].y > maximum.y())
1915  maximum.setY(points[i].y);
1916  }
1917 
1918  // Set z extremum values
1919  minimum.setZ(fzmin);
1920  maximum.setZ(fzmax);
1921 
1922  // Assign
1923  pMin = minimum;
1924  pMax = maximum;
1925 }
1926 
1927 void takePolyhedron(const HepPolyhedron& p)
1928 {
1929  int i, nnode, iNodes[5], iVis[4], iFaces[4];
1930 
1931  for (int iface = 1; iface <= p.GetNoFacets(); iface++) {
1932  p.GetFacet(iface, nnode, iNodes, iVis, iFaces);
1933  for (i = 0; i < nnode; i++) {
1934  if (iNodes[i] < 1 || iNodes[i] > p.GetNoVertices()) { //G.Barrand
1935  // processor_error = 1;
1936  G4cerr
1937  << "BooleanProcessor::takePolyhedron : problem 1."
1938  << G4endl;
1939  }
1940  if (iFaces[i] < 1 || iFaces[i] > p.GetNoFacets()) { //G.Barrand
1941  // processor_error = 1;
1942  G4cerr
1943  << "BooleanProcessor::takePolyhedron : problem 2. "
1944  << i << " " << iFaces[i] << " " << p.GetNoFacets() << G4endl;
1945  }
1946  }
1947  }
1948 }
1949 
1950 PolyhedronBelleLathe::PolyhedronBelleLathe(const std::vector<zr_t>& v, const std::vector<triangle_t>& t, double phi, double dphi)
1951 {
1952  int nphi = GetNumberOfRotationSteps();
1953  bool twopi = abs(dphi - 2 * M_PI) < 1e-6;
1954  int n = v.size();
1955  if (twopi) {
1956  int nv = n * nphi;
1957  int nf = nv;
1958  AllocateMemory(nv, nf);
1959 
1960  auto vnum = [nphi, n](int iphi, int ip) {
1961  return (iphi % nphi) * n + (ip % n) + 1;
1962  };
1963 
1964  int fcount = 1;
1965  double dfi = dphi / nphi;
1966  for (int i = 0; i < nphi; i++) {
1967  double fi = phi + i * dfi;
1968  double cf = cos(fi), sf = sin(fi);
1969  for (int j = 0; j < n; j++) pV[vnum(i, j)].set(v[j].r * cf, v[j].r * sf, v[j].z);
1970  for (int j = 0; j < n; j++) pF[fcount++ ] = G4Facet(vnum(i, j), 0, vnum(i, j + 1), 0, vnum(i + 1, j + 1), 0, vnum(i + 1, j), 0);
1971  }
1972  } else {
1973  // cout<<"NPHI = "<<nphi<<" "<<phi<<" "<<dphi<<endl;
1974  nphi = int(nphi * (dphi / (2 * M_PI)) + 0.5);
1975  nphi = nphi > 3 ? nphi : 3;
1976 
1977  // cout<<"NPHI = "<<nphi<<endl;
1978 
1979  int nv = n * nphi;
1980  int nf = n * (nphi - 1) + 2 * t.size();
1981  AllocateMemory(nv, nf);
1982 
1983  auto vnum = [n](int iphi, int ip) {
1984  return iphi * n + (ip % n) + 1;
1985  };
1986 
1987  int fcount = 1;
1988  double dfi = dphi / (nphi - 1);
1989  for (int i = 0; i < nphi; i++) {
1990  double fi = phi + i * dfi;
1991  double cf = cos(fi), sf = sin(fi);
1992  for (int j = 0; j < n; j++) pV[vnum(i, j)].set(v[j].r * cf, v[j].r * sf, v[j].z);
1993  if (i == nphi - 1) break;
1994  for (int j = 0; j < n; j++) pF[fcount++] = G4Facet(vnum(i, j), 0, vnum(i, j + 1), 0, vnum(i + 1, j + 1), 0, vnum(i + 1, j), 0);
1995  }
1996 
1997  for (const triangle_t& k : t) pF[fcount++] = G4Facet(vnum(0, k.i0), 0, vnum(0, k.i2), 0, vnum(0, k.i1), 0, 0, 0);
1998  int i = nphi - 1;
1999  for (const triangle_t& k : t) pF[fcount++] = G4Facet(vnum(i, k.i0), 0, vnum(i, k.i1), 0, vnum(i, k.i2), 0, 0, 0);
2000 
2001  }
2002  SetReferences();
2003  // takePolyhedron(*this);
2004 }
2005 
2007 
2008 G4Polyhedron* BelleLathe::CreatePolyhedron() const
2009 {
2010  eartrim();
2012 }
2013 
2014 #if 0
2015 #include <immintrin.h>
2016 double mindistsimd(const zr_t& r, const vector<cachezr_t>& contour)
2017 {
2018  double d = kInfinity;
2019  int i = 0, n = contour.size();
2020  __m128d zero = _mm_set_sd(0);
2021  __m128d one = _mm_set_sd(1);
2022  double wn = 0;
2023  do {
2024  const cachezr_t& s = contour[i];
2025  double dz = r.z - s.z, dr = r.r - s.r;
2026  double crs = s.dr * dz - s.dz * dr;
2027  double dot = s.dz * dz + s.dr * dr; // projection of the point on the segement
2028  // if(s.zmin<=r.z&&r.z<s.zmax) wn -= (crs>0) - (crs<0);
2029  __m128d crssd = _mm_set_sd(crs);
2030  __m128d maskgt = _mm_cmpgt_sd(crssd, zero);
2031  __m128d masklt = _mm_cmplt_sd(crssd, zero);
2032  __m128d left = _mm_sub_sd(_mm_and_pd(maskgt, one), _mm_and_pd(masklt, one));
2033  __m128d z = _mm_set_sd(s.z);
2034  __m128d mask = _mm_and_pd(_mm_cmple_sd(_mm_set_sd(s.zmin), z), _mm_cmplt_sd(z, _mm_set_sd(s.zmax)));
2035  left = _mm_and_pd(mask, left);
2036  double du = dz * dz + dr * dr;
2037  double dv = crs * crs * s.is2;
2038 
2039  masklt = _mm_cmplt_sd(_mm_set_sd(dot), zero);
2040  maskgt = _mm_cmpgt_sd(_mm_set_sd(dot), _mm_set_sd(s.s2));
2041 
2042  __m128d uu = _mm_or_pd(_mm_and_pd(maskgt, _mm_set_sd(kInfinity)), _mm_andnot_pd(maskgt, _mm_set_sd(dv)));
2043  __m128d vv = _mm_or_pd(_mm_and_pd(masklt, _mm_set_sd(min(d, du))), _mm_andnot_pd(masklt, _mm_set_sd(min(d, uu[0]))));
2044  wn -= left[0];
2045  d = vv[0];
2046  } while (++i < n);
2047  d = sqrt(d);
2048  d = (wn == 2) ? -d : d;
2049  // cout<<wn<<" "<<d<<endl;
2050  // cout<<sqrt(dp)<<" "<<sqrt(dm)<<endl;
2051  return d;
2052 }
2053 inline int left(const zr_t& r0, const zr_t& r1, const zr_t& r)
2054 {
2055  double d = (r1.z - r0.z) * (r.r - r0.r) - (r.z - r0.z) * (r1.r - r0.r);
2056  return (d > 0) - (d < 0);
2057 }
2058 
2059 inline int checkside(const zr_t& s0, const zr_t& s1, const zr_t& r)
2060 {
2061  double zmin = min(s0.z, s1.z), zmax = max(s0.z, s1.z);
2062  if (zmin <= r.z && r.z < zmax) return left(s0, s1, r);
2063  return 0;
2064 }
2065 
2066 int wn_poly(const zr_t& r, const vector<zr_t>& contour)
2067 {
2068  int wn = 0; // the winding number counter
2069  int i = 0, n = contour.size() - 1;
2070  do {
2071  wn += checkside(contour[i], contour[i + 1], r);
2072  } while (++i < n);
2073  wn += checkside(contour[n], contour[0], r);
2074  return wn;
2075 }
2076 
2077 double mindist(const zr_t& r, const vector<zr_t>& contour)
2078 {
2079  int wn = 0;
2080  double d = kInfinity;
2081  auto dist = [&contour, &d, &r, &wn](int i0, int i1)->void {
2082  const zr_t& s0 = contour[i0], &s1 = contour[i1];
2083  double zmin = min(s0.z, s1.z), zmax = max(s0.z, s1.z);
2084  double sz = s1.z - s0.z, sr = s1.r - s0.r;
2085  double dz = r.z - s0.z, dr = r.r - s0.r;
2086  double crs = dz * sr - sz * dr;
2087  if (zmin <= r.z && r.z < zmax) wn -= (crs > 0) - (crs < 0);
2088  double dot = sz * dz + sr * dr; // projection of the point on the segement
2089  double s2 = sz * sz + sr * sr;
2090  if (dot > s2) return; // point should be within the segment
2091  if (dot < 0)
2092  {
2093  double d2 = dz * dz + dr * dr; // distance to the first point of the segement
2094  d = min(d, d2);
2095  } else
2096  {
2097  d = min(d, crs * crs / s2);
2098  }
2099  // cout<<i0<<" "<<s0.z<<" "<<s0.r<<" "<<d<<" "<<wn<<endl;
2100  };
2101  int i = 0, n = contour.size() - 1;
2102  do {dist(i, i + 1);} while (++i < n);
2103  dist(n, 0);
2104  d = sqrt(d);
2105  d = (wn == 2) ? -d : d;
2106  // cout<<wn<<" "<<d<<endl;
2107  // cout<<sqrt(dp)<<" "<<sqrt(dm)<<endl;
2108  return d;
2109 }
2110 
2111 double mindist(const zr_t& r, const vector<cachezr_t>& contour)
2112 {
2113  double d = kInfinity;
2114  int wn = 0, i = 0, n = contour.size();
2115  do {
2116  const cachezr_t& s = contour[i];
2117  double dz = r.z - s.z, dr = r.r - s.r;
2118  double crs = s.dr * dz - s.dz * dr;
2119  double dot = s.dz * dz + s.dr * dr; // projection of the point on the segement
2120  if (s.zmin <= r.z && r.z < s.zmax) wn -= (crs > 0) - (crs < 0);
2121  if (dot > s.s2) continue; // point should be within the segment
2122  if (dot < 0) {
2123  d = min(d, dz * dz + dr * dr); // distance to the first point of the segement
2124  } else {
2125  d = min(d, crs * crs * s.is2);
2126  }
2127  // cout<<i<<" "<<s.z<<" "<<s.r<<" "<<d<<" "<<wn<<endl;
2128  // cout<<i<<" "<<d<<" "<<wn<<endl;
2129  } while (++i < n);
2130  d = sqrt(d);
2131  d = (wn == 2) ? -d : d;
2132  // cout<<wn<<" "<<d<<endl;
2133  return d;
2134 }
2135 
2136 
2137 #endif
R E
internal precision of FFTW codelets
double R
typedef autogenerated by FFTW
BelleLathe class.
Definition: BelleLathe.h:67
double fzmax
maximal z value
Definition: BelleLathe.h:192
void getvolarea()
get volume area
Definition: BelleLathe.cc:1739
double fphi
starting angle
Definition: BelleLathe.h:179
std::vector< double > fz
vector of z values
Definition: BelleLathe.h:173
G4GeometryType GetEntityType() const
Get entity type.
Definition: BelleLathe.cc:1810
std::vector< int > findx
vector of indices
Definition: BelleLathe.h:174
std::vector< triangle_t > ftlist
vector of triangle structs
Definition: BelleLathe.h:177
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Visualisation function.
Definition: BelleLathe.cc:1849
void Init(const std::vector< zr_t > &, double, double)
initialize
Definition: BelleLathe.cc:101
bool ftwopi
bound within +- 2pi?
Definition: BelleLathe.h:194
std::vector< double > farea
vector of area values
Definition: BelleLathe.h:176
G4Polyhedron * CreatePolyhedron() const
create polyhedron
Definition: BelleLathe.cc:2008
std::vector< int > fseg
vector of segments
Definition: BelleLathe.h:175
G4VSolid * fshape
shape
Definition: BelleLathe.h:196
double fdphi
finishing angle
Definition: BelleLathe.h:180
std::vector< double > linecross(const G4ThreeVector &, const G4ThreeVector &) const
calculate all ray solid's surface intersection return ordered vector
Definition: BelleLathe.cc:428
BelleLathe & operator=(const BelleLathe &rhs)
assignment operator
Definition: BelleLathe.cc:305
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Two vectors define an axis-parallel bounding box for the shape.
Definition: BelleLathe.cc:1855
virtual ~BelleLathe()
Destructor.
Definition: BelleLathe.cc:285
std::vector< zr_t > fcontour
vector of zr structs
Definition: BelleLathe.h:171
bool fgtpi
greater than pi?
Definition: BelleLathe.h:193
std::vector< cachezr_t > fcache
vector of cached zr structs
Definition: BelleLathe.h:172
double mindist(const zr_t &) const
minimal distance
Definition: BelleLathe.cc:962
BelleLathe(const G4String &pName)
Constructor for "nominal" BelleLathe whose parameters are to be set by a G4VPVParamaterisation later.
Definition: BelleLathe.cc:268
double fzmin
minimal z value
Definition: BelleLathe.h:191
zr_t normal(const zr_t &, double &) const
return normal
Definition: BelleLathe.cc:1026
bool insector(double, double) const
True if (x,y) is within the shape rotation.
Definition: BelleLathe.cc:937
G4ThreeVector GetPointOnSurface() const
Get point on surface.
Definition: BelleLathe.cc:1764
double frmax
maximal r value
Definition: BelleLathe.h:190
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
calculate the extent of the volume
Definition: BelleLathe.cc:514
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Calculate side nearest to p, and return normal.
Definition: BelleLathe.cc:1077
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Calculate distance to shape from outside - return kInfinity if no intersection.
Definition: BelleLathe.cc:1253
double frmin
minimal r value
Definition: BelleLathe.h:189
void eartrim() const
ear trim
Definition: BelleLathe.cc:1702
std::vector< G4ThreeVector > fsurf
vector of surfaces
Definition: BelleLathe.h:197
G4VSolid * Clone() const
Make a clone of the object.
Definition: BelleLathe.cc:1816
std::ostream & StreamInfo(std::ostream &os) const
Stream object contents to an output stream.
Definition: BelleLathe.cc:1822
int wn_poly(const zr_t &) const
wn_poly
Definition: BelleLathe.cc:945
EInside Inside(const G4ThreeVector &p) const
Return whether point inside/outside/on surface, using tolerance.
Definition: BelleLathe.cc:983
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
compute the dimensions
Definition: BelleLathe.cc:345
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
Calculate distance to surface of shape from inside.
Definition: BelleLathe.cc:1471
Belle lathe polyhedron.
Definition: BelleLathe.h:201
virtual ~PolyhedronBelleLathe()
destructor
Definition: BelleLathe.cc:2006
PolyhedronBelleLathe(const std::vector< zr_t > &, const std::vector< triangle_t > &, double, double)
constructor
Definition: BelleLathe.cc:1950
std::ostream & operator<<(std::ostream &output, const IntervalOfValidity &iov)
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
T dot(GeneralVector< T > a, GeneralVector< T > b)
dot product of two general vectors
Definition: beamHelpers.h:163
TString rn()
Get random string.
Definition: tools.h:38
Abstract base class for different kinds of events.
struct for plane
Definition: BelleCrystal.h:24
cached z-r struct
Definition: BelleLathe.h:31
double dz
difference in z
Definition: BelleLathe.h:34
double r
r coordinate
Definition: BelleLathe.h:33
double zmin
minimal z value
Definition: BelleLathe.h:39
double z
z coordinate
Definition: BelleLathe.h:32
double ta
ratio of dr over dz
Definition: BelleLathe.h:43
double zmax
maximal z value
Definition: BelleLathe.h:40
struct for a triangle
Definition: BelleLathe.h:55
struct for a three vector
Definition: BelleLathe.h:48
double y
y coordinate
Definition: BelleLathe.h:50
double x
x coordinate
Definition: BelleLathe.h:49
simple struct with z and r coordinates
Definition: BelleLathe.h:25
double r
r coordinate
Definition: BelleLathe.h:27
double z
z coordinate
Definition: BelleLathe.h:26
G4ThreeVector n
Normal unit vector (x,y,z)
Definition: BelleLathe.cc:42
double d
offset (d)
Definition: BelleLathe.cc:43
curl struct
Definition: BelleLathe.cc:67
curl_t(const G4ThreeVector &_v)
constructor
Definition: BelleLathe.cc:70
G4ThreeVector v
vector
Definition: BelleLathe.cc:68
solution struct
Definition: BelleLathe.cc:392
double t
t
Definition: BelleLathe.cc:393
double s
s
Definition: BelleLathe.cc:394