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