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