Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TwistTrapParallelSide.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4TwistTrapParallelSide implementation
27//
28// Author: Oliver Link (CERN), 27.10.2004 - Created
29// --------------------------------------------------------------------
30
31#include <cmath>
32
36
37//=====================================================================
38//* constructors ------------------------------------------------------
39
41 G4double PhiTwist, // twist angle
42 G4double pDz, // half z lenght
43 G4double pTheta, // direction between end planes
44 G4double pPhi, // by polar and azimutal angles
45 G4double pDy1, // half y length at -pDz
46 G4double pDx1, // half x length at -pDz,-pDy
47 G4double pDx2, // half x length at -pDz,+pDy
48 G4double pDy2, // half y length at +pDz
49 G4double pDx3, // half x length at +pDz,-pDy
50 G4double pDx4, // half x length at +pDz,+pDy
51 G4double pAlph, // tilt angle at +pDz
52 G4double AngleSide // parity
53 )
54 : G4VTwistSurface(name)
55{
56
57 fAxis[0] = kXAxis; // in local coordinate system
58 fAxis[1] = kZAxis;
59 fAxisMin[0] = -kInfinity ; // X Axis boundary
60 fAxisMax[0] = kInfinity ; // depends on z !!
61 fAxisMin[1] = -pDz ; // Z Axis boundary
62 fAxisMax[1] = pDz ;
63
64 fDx1 = pDx1 ;
65 fDx2 = pDx2 ;
66 fDx3 = pDx3 ;
67 fDx4 = pDx4 ;
68
69 fDy1 = pDy1 ;
70 fDy2 = pDy2 ;
71
72 fDz = pDz ;
73
74 fAlph = pAlph ;
75 fTAlph = std::tan(fAlph) ;
76
77 fTheta = pTheta ;
78 fPhi = pPhi ;
79
80 // precalculate frequently used parameters
81 //
82 fDx4plus2 = fDx4 + fDx2 ;
83 fDx4minus2 = fDx4 - fDx2 ;
84 fDx3plus1 = fDx3 + fDx1 ;
85 fDx3minus1 = fDx3 - fDx1 ;
86 fDy2plus1 = fDy2 + fDy1 ;
87 fDy2minus1 = fDy2 - fDy1 ;
88
89 fa1md1 = 2*fDx2 - 2*fDx1 ;
90 fa2md2 = 2*fDx4 - 2*fDx3 ;
91
92 fPhiTwist = PhiTwist ; // dphi
93 fAngleSide = AngleSide ; // 0,90,180,270 deg
94
95 fdeltaX = 2*fDz*std::tan(fTheta)*std::cos(fPhi); // dx in surface equation
96 fdeltaY = 2*fDz*std::tan(fTheta)*std::sin(fPhi); // dy in surface equation
97
98 fRot.rotateZ( AngleSide ) ;
99
100 fTrans.set(0, 0, 0); // No Translation
101 fIsValidNorm = false;
102
103 SetCorners() ;
104 SetBoundaries() ;
105}
106
107//=====================================================================
108//* Fake default constructor ------------------------------------------
109
114
115//=====================================================================
116//* GetNormal ---------------------------------------------------------
117
119 G4bool isGlobal)
120{
121 // GetNormal returns a normal vector at a surface (or very close
122 // to surface) point at tmpxx.
123 // If isGlobal=true, it returns the normal in global coordinate.
124 //
125
126 G4ThreeVector xx;
127 if (isGlobal)
128 {
129 xx = ComputeLocalPoint(tmpxx);
130 if ((xx - fCurrentNormal.p).mag() < 0.5 * kCarTolerance)
131 {
133 }
134 }
135 else
136 {
137 xx = tmpxx;
138 if (xx == fCurrentNormal.p)
139 {
140 return fCurrentNormal.normal;
141 }
142 }
143
144 G4double phi ;
145 G4double u ;
146
147 GetPhiUAtX(xx,phi,u) ; // phi,u for point xx close to surface
148
149 G4ThreeVector normal = NormAng(phi,u) ; // the normal vector at phi,u
150
151#ifdef G4TWISTDEBUG
152 G4cout << "normal vector = " << normal << G4endl ;
153 G4cout << "phi = " << phi << " , u = " << u << G4endl ;
154#endif
155
156 // normal = normal/normal.mag() ;
157
158 if (isGlobal)
159 {
160 fCurrentNormal.normal = ComputeGlobalDirection(normal.unit());
161 }
162 else
163 {
164 fCurrentNormal.normal = normal.unit();
165 }
166 return fCurrentNormal.normal;
167}
168
169//=====================================================================
170//* DistanceToSurface -------------------------------------------------
171
173 const G4ThreeVector& gv,
174 G4ThreeVector gxx[],
175 G4double distance[],
176 G4int areacode[],
177 G4bool isvalid[],
178 EValidate validate)
179{
180 static const G4double pihalf = pi/2 ;
181 const G4double ctol = 0.5 * kCarTolerance;
182
183 G4bool IsParallel = false ;
184 G4bool IsConverged = false ;
185
186 G4int nxx = 0 ; // number of physical solutions
187
188 fCurStatWithV.ResetfDone(validate, &gp, &gv);
189
190 if (fCurStatWithV.IsDone())
191 {
192 for (G4int i=0; i<fCurStatWithV.GetNXX(); ++i)
193 {
194 gxx[i] = fCurStatWithV.GetXX(i);
195 distance[i] = fCurStatWithV.GetDistance(i);
196 areacode[i] = fCurStatWithV.GetAreacode(i);
197 isvalid[i] = fCurStatWithV.IsValid(i);
198 }
199 return fCurStatWithV.GetNXX();
200 }
201
202 // initialize
203 for (G4int i=0; i<G4VSURFACENXX ; ++i)
204 {
205 distance[i] = kInfinity;
206 areacode[i] = sOutside;
207 isvalid[i] = false;
208 gxx[i].set(kInfinity, kInfinity, kInfinity);
209 }
210
213
214#ifdef G4TWISTDEBUG
215 G4cout << "Local point p = " << p << G4endl ;
216 G4cout << "Local direction v = " << v << G4endl ;
217#endif
218
219 G4double phi,u ; // parameters
220
221 // temporary variables
222
223 G4double tmpdist = kInfinity ;
224 G4ThreeVector tmpxx;
225 G4int tmpareacode = sOutside ;
226 G4bool tmpisvalid = false ;
227
228 std::vector<Intersection> xbuf ;
229 Intersection xbuftmp ;
230
231 // prepare some variables for the intersection finder
232
233 G4double L = 2*fDz ;
234
235 G4double phixz = fPhiTwist * ( p.x() * v.z() - p.z() * v.x() ) ;
236 G4double phiyz = fPhiTwist * ( p.y() * v.z() - p.z() * v.y() ) ;
237
238 // special case vz = 0
239
240 if ( v.z() == 0. )
241 {
242 if ( std::fabs(p.z()) <= L ) // intersection possible in z
243 {
244 phi = p.z() * fPhiTwist / L ; // phi is determined by the z-position
245
246 u = (2*(fdeltaY*phi*v.x() - fPhiTwist*p.y()*v.x() - fdeltaX*phi*v.y()
247 + fPhiTwist*p.x()*v.y()) + (fDy2plus1*fPhiTwist
248 + 2*fDy2minus1*phi)*(v.x()*std::cos(phi) + v.y()*std::sin(phi)))
249 / (2.* fPhiTwist*(v.y()*std::cos(phi) - v.x()*std::sin(phi)));
250
251 xbuftmp.phi = phi ;
252 xbuftmp.u = u ;
253 xbuftmp.areacode = sOutside ;
254 xbuftmp.distance = kInfinity ;
255 xbuftmp.isvalid = false ;
256
257 xbuf.push_back(xbuftmp) ; // store it to xbuf
258 }
259 else // no intersection possible
260 {
261 distance[0] = kInfinity;
262 gxx[0].set(kInfinity,kInfinity,kInfinity);
263 isvalid[0] = false ;
264 areacode[0] = sOutside ;
265 fCurStatWithV.SetCurrentStatus(0, gxx[0], distance[0],
266 areacode[0], isvalid[0],
267 0, validate, &gp, &gv);
268
269 return 0;
270 } // end std::fabs(p.z() <= L
271 } // end v.z() == 0
272 else // general solution for non-zero vz
273 {
274 G4double c[9],srd[8],si[8] ;
275
276 c[8] = -3600*(-2*phiyz + fDy2plus1*fPhiTwist*v.z()) ;
277 c[7] = -7200*(phixz - 2*fDz*v.y() + (fdeltaY + fDy2minus1)*v.z()) ;
278 c[6] = 120*(-52*phiyz - 120*fDz*v.x() + 60*fdeltaX*v.z()
279 + 11*fDy2plus1*fPhiTwist*v.z()) ;
280 c[5] = 240*(16*phixz - 52*fDz*v.y() + 26*fdeltaY*v.z()
281 + 11*fDy2minus1*v.z()) ;
282 c[4] = 12*(127*phiyz + 640*fDz*v.x() - 320*fdeltaX*v.z()
283 + 4*fDy2plus1*fPhiTwist*v.z()) ;
284 c[3] = -404*phixz + 3048*fDz*v.y() - 1524*fdeltaY*v.z()
285 + 96*fDy2minus1*v.z() ;
286 c[2] = -72*phiyz + 404*(-2*fDz*v.x() + fdeltaX*v.z()) ;
287 c[1] = 12*(phixz - 12*fDz*v.y() + 6*fdeltaY*v.z()) ;
288 c[0] = 24*fDz*v.x() - 12*fdeltaX*v.z() ;
289
290
291#ifdef G4TWISTDEBUG
292 G4cout << "coef = " << c[0] << " "
293 << c[1] << " "
294 << c[2] << " "
295 << c[3] << " "
296 << c[4] << " "
297 << c[5] << " "
298 << c[6] << " "
299 << c[7] << " "
300 << c[8] << G4endl ;
301#endif
302
303 G4JTPolynomialSolver trapEq ;
304 G4int num = trapEq.FindRoots(c,8,srd,si);
305
306 for (G4int i = 0 ; i<num ; ++i ) // loop over all mathematical solutions
307 {
308 if ( si[i]==0.0 ) // only real solutions
309 {
310#ifdef G4TWISTDEBUG
311 G4cout << "Solution " << i << " : " << srd[i] << G4endl ;
312#endif
313 phi = std::fmod(srd[i] , pihalf) ;
314 u = (1/std::cos(phi)*(2*phixz + 4*fDz*phi*v.x()
315 - 2*fdeltaX*phi*v.z() + (fDy2plus1*fPhiTwist
316 + 2*fDy2minus1*phi)*v.z()* std::sin(phi)))/(2.*fPhiTwist*v.z()) ;
317
318 xbuftmp.phi = phi ;
319 xbuftmp.u = u ;
320 xbuftmp.areacode = sOutside ;
321 xbuftmp.distance = kInfinity ;
322 xbuftmp.isvalid = false ;
323
324 xbuf.push_back(xbuftmp) ; // store it to xbuf
325
326#ifdef G4TWISTDEBUG
327 G4cout << "solution " << i << " = " << phi << " , " << u << G4endl ;
328#endif
329
330 } // end if real solution
331 } // end loop i
332 } // end general case
333
334 nxx = (G4int)xbuf.size() ; // save the number of solutions
335
336 G4ThreeVector xxonsurface ; // point on surface
337 G4ThreeVector surfacenormal ; // normal vector
338 G4double deltaX ; // distance between intersection point and point on surface
339 G4double theta ; // angle between track and surfacenormal
340 G4double factor ; // a scaling factor
341 G4int maxint = 30 ; // number of iterations
342
343
344 for (auto & k : xbuf)
345 {
346#ifdef G4TWISTDEBUG
347 G4cout << "Solution " << k << " : "
348 << "reconstructed phiR = " << xbuf[k].phi
349 << ", uR = " << xbuf[k].u << G4endl ;
350#endif
351
352 phi = k.phi ; // get the stored values for phi and u
353 u = k.u ;
354
355 IsConverged = false ; // no convergence at the beginning
356
357 for ( G4int i = 1 ; i<maxint ; ++i )
358 {
359 xxonsurface = SurfacePoint(phi,u) ;
360 surfacenormal = NormAng(phi,u) ;
361 tmpdist = DistanceToPlaneWithV(p, v, xxonsurface, surfacenormal, tmpxx);
362 deltaX = ( tmpxx - xxonsurface ).mag() ;
363 theta = std::fabs(std::acos(v*surfacenormal) - pihalf) ;
364 if ( theta < 0.001 )
365 {
366 factor = 50 ;
367 IsParallel = true ;
368 }
369 else
370 {
371 factor = 1 ;
372 }
373
374#ifdef G4TWISTDEBUG
375 G4cout << "Step i = " << i << ", distance = "
376 << tmpdist << ", " << deltaX << G4endl ;
377 G4cout << "X = " << tmpxx << G4endl ;
378#endif
379
380 GetPhiUAtX(tmpxx, phi, u); // new point xx is accepted and phi/u replaced
381
382#ifdef G4TWISTDEBUG
383 G4cout << "approximated phi = " << phi << ", u = " << u << G4endl ;
384#endif
385
386 if ( deltaX <= factor*ctol ) { IsConverged = true ; break ; }
387 } // end iterative loop (i)
388
389 if ( std::fabs(tmpdist)<ctol ) { tmpdist = 0 ; }
390
391#ifdef G4TWISTDEBUG
392 G4cout << "refined solution " << phi << " , " << u << G4endl ;
393 G4cout << "distance = " << tmpdist << G4endl ;
394 G4cout << "local X = " << tmpxx << G4endl ;
395#endif
396
397 tmpisvalid = false ; // init
398
399 if ( IsConverged )
400 {
401 if (validate == kValidateWithTol)
402 {
403 tmpareacode = GetAreaCode(tmpxx);
404 if (!IsOutside(tmpareacode))
405 {
406 if (tmpdist >= 0) { tmpisvalid = true; }
407 }
408 }
409 else if (validate == kValidateWithoutTol)
410 {
411 tmpareacode = GetAreaCode(tmpxx, false);
412 if (IsInside(tmpareacode))
413 {
414 if (tmpdist >= 0) { tmpisvalid = true; }
415 }
416 }
417 else // kDontValidate
418 {
419 G4Exception("G4TwistTrapParallelSide::DistanceToSurface()",
420 "GeomSolids0001", FatalException,
421 "Feature NOT implemented !");
422 }
423 }
424 else
425 {
426 tmpdist = kInfinity; // no convergence after 10 steps
427 tmpisvalid = false ; // solution is not vaild
428 }
429
430 // store the found values
431 k.xx = tmpxx ;
432 k.distance = tmpdist ;
433 k.areacode = tmpareacode ;
434 k.isvalid = tmpisvalid ;
435 } // end loop over physical solutions (variable k)
436
437 std::sort(xbuf.begin() , xbuf.end(), DistanceSort ) ; // sorting
438
439#ifdef G4TWISTDEBUG
440 G4cout << G4endl << "list xbuf after sorting : " << G4endl ;
441 G4cout << G4endl << G4endl ;
442#endif
443
444 // erase identical intersection (within kCarTolerance)
445 xbuf.erase(std::unique(xbuf.begin(),xbuf.end(),EqualIntersection),xbuf.end());
446
447
448 // add guesses
449
450 auto nxxtmp = (G4int)xbuf.size() ;
451
452 if ( nxxtmp<2 || IsParallel )
453 {
454 // positive end
455#ifdef G4TWISTDEBUG
456 G4cout << "add guess at +z/2 .. " << G4endl ;
457#endif
458
459 phi = fPhiTwist/2 ;
460 u = 0 ;
461
462 xbuftmp.phi = phi ;
463 xbuftmp.u = u ;
464 xbuftmp.areacode = sOutside ;
465 xbuftmp.distance = kInfinity ;
466 xbuftmp.isvalid = false ;
467
468 xbuf.push_back(xbuftmp) ; // store it to xbuf
469
470#ifdef G4TWISTDEBUG
471 G4cout << "add guess at -z/2 .. " << G4endl ;
472#endif
473
474 phi = -fPhiTwist/2 ;
475 u = 0 ;
476
477 xbuftmp.phi = phi ;
478 xbuftmp.u = u ;
479 xbuftmp.areacode = sOutside ;
480 xbuftmp.distance = kInfinity ;
481 xbuftmp.isvalid = false ;
482
483 xbuf.push_back(xbuftmp) ; // store it to xbuf
484
485 for ( std::size_t k = nxxtmp ; k<xbuf.size() ; ++k )
486 {
487#ifdef G4TWISTDEBUG
488 G4cout << "Solution " << k << " : "
489 << "reconstructed phiR = " << xbuf[k].phi
490 << ", uR = " << xbuf[k].u << G4endl ;
491#endif
492
493 phi = xbuf[k].phi ; // get the stored values for phi and u
494 u = xbuf[k].u ;
495
496 IsConverged = false ; // no convergence at the beginning
497
498 for ( G4int i = 1 ; i<maxint ; ++i )
499 {
500 xxonsurface = SurfacePoint(phi,u) ;
501 surfacenormal = NormAng(phi,u) ;
502 tmpdist = DistanceToPlaneWithV(p, v, xxonsurface, surfacenormal, tmpxx);
503 deltaX = ( tmpxx - xxonsurface ).mag() ;
504 theta = std::fabs(std::acos(v*surfacenormal) - pihalf) ;
505 if ( theta < 0.001 )
506 {
507 factor = 50 ;
508 }
509 else
510 {
511 factor = 1 ;
512 }
513
514#ifdef G4TWISTDEBUG
515 G4cout << "Step i = " << i << ", distance = "
516 << tmpdist << ", " << deltaX << G4endl ;
517 G4cout << "X = " << tmpxx << G4endl ;
518#endif
519
520 GetPhiUAtX(tmpxx, phi, u) ; // new point xx accepted and phi/u replaced
521
522#ifdef G4TWISTDEBUG
523 G4cout << "approximated phi = " << phi << ", u = " << u << G4endl ;
524#endif
525
526 if ( deltaX <= factor*ctol ) { IsConverged = true ; break ; }
527 } // end iterative loop (i)
528
529 if ( std::fabs(tmpdist)<ctol ) { tmpdist = 0 ; }
530
531#ifdef G4TWISTDEBUG
532 G4cout << "refined solution " << phi << " , " << u << G4endl ;
533 G4cout << "distance = " << tmpdist << G4endl ;
534 G4cout << "local X = " << tmpxx << G4endl ;
535#endif
536
537 tmpisvalid = false ; // init
538
539 if ( IsConverged )
540 {
541 if (validate == kValidateWithTol)
542 {
543 tmpareacode = GetAreaCode(tmpxx);
544 if (!IsOutside(tmpareacode))
545 {
546 if (tmpdist >= 0) { tmpisvalid = true; }
547 }
548 }
549 else if (validate == kValidateWithoutTol)
550 {
551 tmpareacode = GetAreaCode(tmpxx, false);
552 if (IsInside(tmpareacode))
553 {
554 if (tmpdist >= 0) { tmpisvalid = true; }
555 }
556 }
557 else // kDontValidate
558 {
559 G4Exception("G4TwistedBoxSide::DistanceToSurface()",
560 "GeomSolids0001", FatalException,
561 "Feature NOT implemented !");
562 }
563
564 }
565 else
566 {
567 tmpdist = kInfinity; // no convergence after 10 steps
568 tmpisvalid = false ; // solution is not vaild
569 }
570
571 // store the found values
572 xbuf[k].xx = tmpxx ;
573 xbuf[k].distance = tmpdist ;
574 xbuf[k].areacode = tmpareacode ;
575 xbuf[k].isvalid = tmpisvalid ;
576
577 } // end loop over physical solutions
578 } // end less than 2 solutions
579
580 // sort again
581 std::sort(xbuf.begin() , xbuf.end(), DistanceSort ) ; // sorting
582
583 // erase identical intersection (within kCarTolerance)
584 xbuf.erase(std::unique(xbuf.begin(),xbuf.end(),EqualIntersection),xbuf.end());
585
586#ifdef G4TWISTDEBUG
587 G4cout << G4endl << "list xbuf after sorting : " << G4endl ;
588 G4cout << G4endl << G4endl ;
589#endif
590
591 nxx = (G4int)xbuf.size() ; // determine number of solutions again.
592
593 for ( G4int i = 0 ; i<(G4int)xbuf.size() ; ++i )
594 {
595 distance[i] = xbuf[i].distance;
596 gxx[i] = ComputeGlobalPoint(xbuf[i].xx);
597 areacode[i] = xbuf[i].areacode ;
598 isvalid[i] = xbuf[i].isvalid ;
599
600 fCurStatWithV.SetCurrentStatus(i, gxx[i], distance[i], areacode[i],
601 isvalid[i], nxx, validate, &gp, &gv);
602#ifdef G4TWISTDEBUG
603 G4cout << "element Nr. " << i
604 << ", local Intersection = " << xbuf[i].xx
605 << ", distance = " << xbuf[i].distance
606 << ", u = " << xbuf[i].u
607 << ", phi = " << xbuf[i].phi
608 << ", isvalid = " << xbuf[i].isvalid
609 << G4endl ;
610#endif
611 } // end for( i ) loop
612
613#ifdef G4TWISTDEBUG
614 G4cout << "G4TwistTrapParallelSide finished " << G4endl ;
615 G4cout << nxx << " possible physical solutions found" << G4endl ;
616 for ( G4int k= 0 ; k< nxx ; ++k )
617 {
618 G4cout << "global intersection Point found: " << gxx[k] << G4endl ;
619 G4cout << "distance = " << distance[k] << G4endl ;
620 G4cout << "isvalid = " << isvalid[k] << G4endl ;
621 }
622#endif
623
624 return nxx ;
625}
626
627//=====================================================================
628//* DistanceToSurface -------------------------------------------------
629
631 G4ThreeVector gxx[],
632 G4double distance[],
633 G4int areacode[])
634{
635 const G4double ctol = 0.5 * kCarTolerance;
636
637 fCurStat.ResetfDone(kDontValidate, &gp);
638
639 if (fCurStat.IsDone())
640 {
641 for (G4int i=0; i<fCurStat.GetNXX(); ++i)
642 {
643 gxx[i] = fCurStat.GetXX(i);
644 distance[i] = fCurStat.GetDistance(i);
645 areacode[i] = fCurStat.GetAreacode(i);
646 }
647 return fCurStat.GetNXX();
648 }
649
650 // initialize
651 for (G4int i=0; i<G4VSURFACENXX; ++i)
652 {
653 distance[i] = kInfinity;
654 areacode[i] = sOutside;
655 gxx[i].set(kInfinity, kInfinity, kInfinity);
656 }
657
659 G4ThreeVector xx; // intersection point
660 G4ThreeVector xxonsurface ; // interpolated intersection point
661
662 // the surfacenormal at that surface point
663 G4double phiR = 0 ; //
664 G4double uR = 0 ;
665
666 G4ThreeVector surfacenormal ;
667 G4double deltaX ;
668
669 G4int maxint = 20 ;
670
671 for ( G4int i = 1 ; i<maxint ; ++i )
672 {
673 xxonsurface = SurfacePoint(phiR,uR) ;
674 surfacenormal = NormAng(phiR,uR) ;
675 distance[0] = DistanceToPlane(p, xxonsurface, surfacenormal, xx); // new XX
676 deltaX = ( xx - xxonsurface ).mag() ;
677
678#ifdef G4TWISTDEBUG
679 G4cout << "i = " << i << ", distance = "
680 << distance[0] << ", " << deltaX << G4endl ;
681 G4cout << "X = " << xx << G4endl ;
682#endif
683
684 // the new point xx is accepted and phi/psi replaced
685 GetPhiUAtX(xx, phiR, uR) ;
686
687 if ( deltaX <= ctol ) { break ; }
688 }
689
690 // check validity of solution ( valid phi,psi )
691
692 G4double halfphi = 0.5*fPhiTwist ;
693 G4double uMax = GetBoundaryMax(phiR) ;
694 G4double uMin = GetBoundaryMin(phiR) ;
695
696 if ( phiR > halfphi ) { phiR = halfphi ; }
697 if ( phiR < -halfphi ) { phiR = -halfphi ; }
698 if ( uR > uMax ) { uR = uMax ; }
699 if ( uR < uMin ) { uR = uMin ; }
700
701 xxonsurface = SurfacePoint(phiR,uR) ;
702 distance[0] = ( p - xx ).mag() ;
703 if ( distance[0] <= ctol ) { distance[0] = 0 ; }
704
705 // end of validity
706
707#ifdef G4TWISTDEBUG
708 G4cout << "refined solution " << phiR << " , " << uR << " , " << G4endl ;
709 G4cout << "distance = " << distance[0] << G4endl ;
710 G4cout << "X = " << xx << G4endl ;
711#endif
712
713 G4bool isvalid = true;
714 gxx[0] = ComputeGlobalPoint(xx);
715
716#ifdef G4TWISTDEBUG
717 G4cout << "intersection Point found: " << gxx[0] << G4endl ;
718 G4cout << "distance = " << distance[0] << G4endl ;
719#endif
720
721 fCurStat.SetCurrentStatus(0, gxx[0], distance[0], areacode[0],
722 isvalid, 1, kDontValidate, &gp);
723 return 1;
724}
725
726//=====================================================================
727//* GetAreaCode -------------------------------------------------------
728
729G4int G4TwistTrapParallelSide::GetAreaCode(const G4ThreeVector& xx,
730 G4bool withTol)
731{
732 // We must use the function in local coordinate system.
733 // See the description of DistanceToSurface(p,v).
734
735 const G4double ctol = 0.5 * kCarTolerance;
736
737 G4double phi ;
738 G4double yprime ;
739 GetPhiUAtX(xx, phi,yprime ) ;
740
741 G4double fXAxisMax = GetBoundaryMax(phi) ;
742 G4double fXAxisMin = GetBoundaryMin(phi) ;
743
744#ifdef G4TWISTDEBUG
745 G4cout << "GetAreaCode: phi = " << phi << G4endl ;
746 G4cout << "GetAreaCode: yprime = " << yprime << G4endl ;
747 G4cout << "Intervall is " << fXAxisMin << " to " << fXAxisMax << G4endl ;
748#endif
749
750 G4int areacode = sInside;
751
752 if (fAxis[0] == kXAxis && fAxis[1] == kZAxis)
753 {
754 G4int zaxis = 1;
755
756 if (withTol)
757 {
758 G4bool isoutside = false;
759
760 // test boundary of xaxis
761
762 if (yprime < fXAxisMin + ctol)
763 {
764 areacode |= (sAxis0 & (sAxisX | sAxisMin)) | sBoundary;
765 if (yprime <= fXAxisMin - ctol) { isoutside = true; }
766 }
767 else if (yprime > fXAxisMax - ctol)
768 {
769 areacode |= (sAxis0 & (sAxisX | sAxisMax)) | sBoundary;
770 if (yprime >= fXAxisMax + ctol) { isoutside = true; }
771 }
772
773 // test boundary of z-axis
774
775 if (xx.z() < fAxisMin[zaxis] + ctol)
776 {
777 areacode |= (sAxis1 & (sAxisZ | sAxisMin));
778
779 if ((areacode & sBoundary) != 0)
780 {
781 areacode |= sCorner; // xx is on corner.
782 }
783 else
784 {
785 areacode |= sBoundary;
786 }
787 if (xx.z() <= fAxisMin[zaxis] - ctol) { isoutside = true; }
788 }
789 else if (xx.z() > fAxisMax[zaxis] - ctol)
790 {
791 areacode |= (sAxis1 & (sAxisZ | sAxisMax));
792
793 if ((areacode & sBoundary) != 0)
794 {
795 areacode |= sCorner; // xx is on corner.
796 }
797 else
798 {
799 areacode |= sBoundary;
800 }
801 if (xx.z() >= fAxisMax[zaxis] + ctol) { isoutside = true; }
802 }
803
804 // if isoutside = true, clear inside bit.
805 // if not on boundary, add axis information.
806
807 if (isoutside)
808 {
809 G4int tmpareacode = areacode & (~sInside);
810 areacode = tmpareacode;
811 }
812 else if ((areacode & sBoundary) != sBoundary)
813 {
814 areacode |= (sAxis0 & sAxisX) | (sAxis1 & sAxisZ);
815 }
816
817 }
818 else
819 {
820 // boundary of y-axis
821
822 if (yprime < fXAxisMin )
823 {
824 areacode |= (sAxis0 & (sAxisX | sAxisMin)) | sBoundary;
825 }
826 else if (yprime > fXAxisMax)
827 {
828 areacode |= (sAxis0 & (sAxisX | sAxisMax)) | sBoundary;
829 }
830
831 // boundary of z-axis
832
833 if (xx.z() < fAxisMin[zaxis])
834 {
835 areacode |= (sAxis1 & (sAxisZ | sAxisMin));
836 if ((areacode & sBoundary) != 0)
837 {
838 areacode |= sCorner; // xx is on corner.
839 }
840 else
841 {
842 areacode |= sBoundary;
843 }
844 }
845 else if (xx.z() > fAxisMax[zaxis])
846 {
847 areacode |= (sAxis1 & (sAxisZ | sAxisMax)) ;
848 if ((areacode & sBoundary) != 0)
849 {
850 areacode |= sCorner; // xx is on corner.
851 }
852 else
853 {
854 areacode |= sBoundary;
855 }
856 }
857
858 if ((areacode & sBoundary) != sBoundary)
859 {
860 areacode |= (sAxis0 & sAxisX) | (sAxis1 & sAxisZ);
861 }
862 }
863 return areacode;
864 }
865
866 G4Exception("G4TwistTrapParallelSide::GetAreaCode()",
867 "GeomSolids0001", FatalException,
868 "Feature NOT implemented !");
869
870 return areacode;
871}
872
873//=====================================================================
874//* SetCorners() ------------------------------------------------------
875
876void G4TwistTrapParallelSide::SetCorners()
877{
878
879 // Set Corner points in local coodinate.
880
881 if (fAxis[0] == kXAxis && fAxis[1] == kZAxis)
882 {
883 G4double x, y, z;
884
885 // corner of Axis0min and Axis1min
886
887 x = -fdeltaX/2. + (-fDx2 + fDy1*fTAlph)*std::cos(fPhiTwist/2.)
888 + fDy1*std::sin(fPhiTwist/2.) ;
889 y = -fdeltaY/2. + fDy1*std::cos(fPhiTwist/2.)
890 + (fDx2 - fDy1*fTAlph)*std::sin(fPhiTwist/2.) ;
891 z = -fDz ;
892
893 SetCorner(sC0Min1Min, x, y, z);
894
895 // corner of Axis0max and Axis1min
896
897 x = -fdeltaX/2. + (fDx2 + fDy1*fTAlph)*std::cos(fPhiTwist/2.)
898 + fDy1*std::sin(fPhiTwist/2.) ;
899 y = -fdeltaY/2. + fDy1*std::cos(fPhiTwist/2.)
900 - (fDx2 + fDy1*fTAlph)*std::sin(fPhiTwist/2.) ;
901 z = -fDz;
902
903 SetCorner(sC0Max1Min, x, y, z);
904
905 // corner of Axis0max and Axis1max
906 x = fdeltaX/2. + (fDx4 + fDy2*fTAlph)*std::cos(fPhiTwist/2.)
907 - fDy2*std::sin(fPhiTwist/2.) ;
908 y = fdeltaY/2. + fDy2*std::cos(fPhiTwist/2.)
909 + (fDx4 + fDy2*fTAlph)*std::sin(fPhiTwist/2.) ;
910 z = fDz ;
911
912 SetCorner(sC0Max1Max, x, y, z);
913
914 // corner of Axis0min and Axis1max
915 x = fdeltaX/2. + (-fDx4 + fDy2*fTAlph)*std::cos(fPhiTwist/2.)
916 - fDy2*std::sin(fPhiTwist/2.) ;
917 y = fdeltaY/2. + fDy2*std::cos(fPhiTwist/2.)
918 + (-fDx4 + fDy2*fTAlph)*std::sin(fPhiTwist/2.) ;
919 z = fDz ;
920
921 SetCorner(sC0Min1Max, x, y, z);
922 }
923 else
924 {
925 G4Exception("G4TwistTrapParallelSide::SetCorners()",
926 "GeomSolids0001", FatalException,
927 "Method NOT implemented !");
928 }
929}
930
931//=====================================================================
932//* SetBoundaries() ---------------------------------------------------
933
934void G4TwistTrapParallelSide::SetBoundaries()
935{
936 // Set direction-unit vector of boundary-lines in local coodinate.
937 //
938
939 G4ThreeVector direction;
940
941 if (fAxis[0] == kXAxis && fAxis[1] == kZAxis)
942 {
943 // sAxis0 & sAxisMin
945 direction = direction.unit();
946 SetBoundary(sAxis0 & (sAxisX | sAxisMin), direction,
948
949 // sAxis0 & sAxisMax
951 direction = direction.unit();
952 SetBoundary(sAxis0 & (sAxisX | sAxisMax), direction,
954
955 // sAxis1 & sAxisMin
957 direction = direction.unit();
958 SetBoundary(sAxis1 & (sAxisZ | sAxisMin), direction,
960
961 // sAxis1 & sAxisMax
963 direction = direction.unit();
964 SetBoundary(sAxis1 & (sAxisZ | sAxisMax), direction,
966 }
967 else
968 {
969 G4Exception("G4TwistTrapParallelSide::SetCorners()",
970 "GeomSolids0001", FatalException,
971 "Feature NOT implemented !");
972 }
973}
974
975//=====================================================================
976//* GetPhiUAtX() ------------------------------------------------------
977
978void
979G4TwistTrapParallelSide::GetPhiUAtX( const G4ThreeVector& p,
980 G4double& phi, G4double& u )
981{
982 // find closest point XX on surface for a given point p
983 // X0 is a point on the surface, d is the direction
984 // ( both for a fixed z = pz)
985
986 // phi is given by the z coordinate of p
987
988 phi = p.z()/(2*fDz)*fPhiTwist ;
989
990 u = ((-(fdeltaX*phi) + fPhiTwist*p.x())* std::cos(phi)
991 + (-(fdeltaY*phi) + fPhiTwist*p.y())*std::sin(phi))/fPhiTwist ;
992}
993
994//=====================================================================
995//* ProjectPoint() ----------------------------------------------------
996
997G4ThreeVector G4TwistTrapParallelSide::ProjectPoint(const G4ThreeVector& p,
998 G4bool isglobal)
999{
1000 // Get Rho at p.z() on Hyperbolic Surface.
1001 G4ThreeVector tmpp;
1002 if (isglobal)
1003 {
1004 tmpp = fRot.inverse()*p - fTrans;
1005 }
1006 else
1007 {
1008 tmpp = p;
1009 }
1010
1011 G4double phi ;
1012 G4double u ;
1013
1014 GetPhiUAtX( tmpp, phi, u ) ; // calculate (phi, u) for p close to surface
1015
1016 G4ThreeVector xx = SurfacePoint(phi,u) ; // transform back to Cartesian coords
1017
1018 if (isglobal)
1019 {
1020 return (fRot * xx + fTrans);
1021 }
1022 return xx;
1023}
1024
1025//=====================================================================
1026//* GetFacets() -------------------------------------------------------
1027
1028void G4TwistTrapParallelSide::GetFacets( G4int k, G4int n, G4double xyz[][3],
1029 G4int faces[][4], G4int iside )
1030{
1031 G4double phi ;
1032 G4double z, u ; // the two parameters for the surface equation
1033 G4ThreeVector p ; // a point on the surface, given by (z,u)
1034
1035 G4int nnode ;
1036 G4int nface ;
1037
1038 G4double umin, umax ;
1039
1040 // calculate the (n-1)*(k-1) vertices
1041
1042 for ( G4int i = 0 ; i<n ; ++i )
1043 {
1044 z = -fDz+i*(2.*fDz)/(n-1) ;
1045 phi = z*fPhiTwist/(2*fDz) ;
1046 umin = GetBoundaryMin(phi) ;
1047 umax = GetBoundaryMax(phi) ;
1048
1049 for ( G4int j = 0 ; j<k ; ++j )
1050 {
1051 nnode = GetNode(i,j,k,n,iside) ;
1052 u = umax - j*(umax-umin)/(k-1) ;
1053 p = SurfacePoint(phi,u,true) ; // surface point in global coords
1054
1055 xyz[nnode][0] = p.x() ;
1056 xyz[nnode][1] = p.y() ;
1057 xyz[nnode][2] = p.z() ;
1058
1059 if ( i<n-1 && j<k-1 ) // conterclock wise filling
1060 {
1061 nface = GetFace(i,j,k,n,iside) ;
1062 faces[nface][0] = GetEdgeVisibility(i,j,k,n,0,-1)
1063 * (GetNode(i ,j ,k,n,iside)+1) ; // fortran numbering
1064 faces[nface][1] = GetEdgeVisibility(i,j,k,n,1,-1)
1065 * (GetNode(i ,j+1,k,n,iside)+1) ;
1066 faces[nface][2] = GetEdgeVisibility(i,j,k,n,2,-1)
1067 * (GetNode(i+1,j+1,k,n,iside)+1) ;
1068 faces[nface][3] = GetEdgeVisibility(i,j,k,n,3,-1)
1069 * (GetNode(i+1,j ,k,n,iside)+1) ;
1070 }
1071 }
1072 }
1073}
const G4double kCarTolerance
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4VSURFACENXX
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double z() const
Hep3Vector unit() const
double x() const
double y() const
void set(double x, double y, double z)
G4int FindRoots(G4double *op, G4int degree, G4double *zeror, G4double *zeroi)
G4int DistanceToSurface(const G4ThreeVector &gp, const G4ThreeVector &gv, G4ThreeVector gxx[], G4double distance[], G4int areacode[], G4bool isvalid[], EValidate validate=kValidateWithTol) override
G4TwistTrapParallelSide(const G4String &name, G4double PhiTwist, G4double pDz, G4double pTheta, G4double pPhi, G4double pDy1, G4double pDx1, G4double pDx2, G4double pDy2, G4double pDx3, G4double pDx4, G4double pAlph, G4double AngleSide)
G4ThreeVector GetNormal(const G4ThreeVector &p, G4bool isGlobal=false) override
static const G4int sC0Min1Min
static const G4int sC0Min1Max
G4VTwistSurface(const G4String &name)
G4double DistanceToPlane(const G4ThreeVector &p, const G4ThreeVector &x0, const G4ThreeVector &n0, G4ThreeVector &xx)
G4int GetNode(G4int i, G4int j, G4int m, G4int n, G4int iside)
static const G4int sOutside
G4ThreeVector ComputeGlobalDirection(const G4ThreeVector &lp) const
static const G4int sAxisMax
static const G4int sAxis0
G4int GetFace(G4int i, G4int j, G4int m, G4int n, G4int iside)
G4RotationMatrix fRot
G4int GetEdgeVisibility(G4int i, G4int j, G4int m, G4int n, G4int number, G4int orientation)
G4ThreeVector ComputeLocalDirection(const G4ThreeVector &gp) const
static const G4int sAxisMin
static const G4int sC0Max1Max
static const G4int sAxis1
G4bool IsInside(G4int areacode, G4bool testbitmode=false) const
G4ThreeVector fTrans
virtual void SetBoundary(const G4int &axiscode, const G4ThreeVector &direction, const G4ThreeVector &x0, const G4int &boundarytype)
G4ThreeVector ComputeLocalPoint(const G4ThreeVector &gp) const
void SetCorner(G4int areacode, G4double x, G4double y, G4double z)
G4ThreeVector GetCorner(G4int areacode) const
static const G4int sBoundary
static const G4int sAxisZ
G4bool IsOutside(G4int areacode) const
static const G4int sCorner
static const G4int sC0Max1Min
static const G4int sInside
CurrentStatus fCurStatWithV
static const G4int sAxisX
G4double DistanceToPlaneWithV(const G4ThreeVector &p, const G4ThreeVector &v, const G4ThreeVector &x0, const G4ThreeVector &n0, G4ThreeVector &xx)
G4ThreeVector ComputeGlobalPoint(const G4ThreeVector &lp) const
G4SurfCurNormal fCurrentNormal
CurrentStatus fCurStat
@ kXAxis
Definition geomdefs.hh:55
@ kZAxis
Definition geomdefs.hh:57