Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParameterisationPolyhedra.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// G4ParameterisationPolyhedra[Rho/Phi/Z] implementation
27//
28// 14.10.03 - P.Arce, Initial version
29// 08.04.04 - I.Hrivnacova, Implemented reflection
30// --------------------------------------------------------------------
31
33
34#include <iomanip>
36#include "G4ThreeVector.hh"
38#include "G4RotationMatrix.hh"
39#include "G4VPhysicalVolume.hh"
40#include "G4LogicalVolume.hh"
41#include "G4ReflectedSolid.hh"
42#include "G4Polyhedra.hh"
43
44//--------------------------------------------------------------------------
47 G4double offset, G4VSolid* msolid,
48 DivisionType divType )
49 : G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid )
50{
51 std::ostringstream message;
52 /* #ifdef G4MULTITHREADED
53 message << "Divisions for G4Polyhedra currently NOT supported in MT-mode."
54 << G4endl
55 << "Sorry! Solid: " << msolid->GetName();
56 G4Exception("G4VParameterisationPolyhedra::G4VParameterisationPolyhedra()",
57 "GeomDiv0001", FatalException, message);
58 #endif */
59
60 G4Polyhedra* msol = (G4Polyhedra*)(msolid);
61 if ((msolid->GetEntityType() != "G4ReflectedSolid") && (msol->IsGeneric()))
62 {
63 message << "Generic construct for G4Polyhedra NOT supported." << G4endl
64 << "Sorry! Solid: " << msol->GetName();
65 G4Exception("G4VParameterisationPolyhedra::G4VParameterisationPolyhedra()",
66 "GeomDiv0001", FatalException, message);
67 }
68 if (msolid->GetEntityType() == "G4ReflectedSolid")
69 {
70 // Get constituent solid
71 //
72 G4VSolid* mConstituentSolid
73 = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid();
74 msol = (G4Polyhedra*)(mConstituentSolid);
75
76 // Get parameters
77 //
78 G4int nofSides = msol->GetOriginalParameters()->numSide;
79 G4int nofZplanes = msol->GetOriginalParameters()->Num_z_planes;
80 G4double* zValues = msol->GetOriginalParameters()->Z_values;
81 G4double* rminValues = msol->GetOriginalParameters()->Rmin;
82 G4double* rmaxValues = msol->GetOriginalParameters()->Rmax;
83
84 // Invert z values, convert radius parameters
85 //
86 G4double* rminValues2 = new G4double[nofZplanes];
87 G4double* rmaxValues2 = new G4double[nofZplanes];
88 G4double* zValuesRefl = new G4double[nofZplanes];
89 for (G4int i=0; i<nofZplanes; ++i)
90 {
91 rminValues2[i] = rminValues[i] * ConvertRadiusFactor(*msol);
92 rmaxValues2[i] = rmaxValues[i] * ConvertRadiusFactor(*msol);
93 zValuesRefl[i] = - zValues[i];
94 }
95
96 G4Polyhedra* newSolid
97 = new G4Polyhedra(msol->GetName(),
98 msol->GetStartPhi(),
99 msol->GetEndPhi() - msol->GetStartPhi(),
100 nofSides,
101 nofZplanes, zValuesRefl, rminValues2, rmaxValues2);
102
103 delete [] rminValues2;
104 delete [] rmaxValues2;
105 delete [] zValuesRefl;
106
107 msol = newSolid;
108 fmotherSolid = newSolid;
109 fReflectedSolid = true;
110 fDeleteSolid = true;
111 }
112}
113
114//------------------------------------------------------------------------
116{
117}
118
119//--------------------------------------------------------------------------
121G4VParameterisationPolyhedra::
122ConvertRadiusFactor(const G4Polyhedra& phedra) const
123{
124 G4double phiTotal = phedra.GetEndPhi() - phedra.GetStartPhi();
125 G4int nofSides = phedra.GetOriginalParameters()->numSide;
126
127 if ( (phiTotal <=0) || (phiTotal >
128 2*pi+G4GeometryTolerance::GetInstance()->GetAngularTolerance()) )
129 { phiTotal = 2*pi; }
130
131 return std::cos(0.5*phiTotal/nofSides);
132}
133
134//--------------------------------------------------------------------------
137 G4double width, G4double offset,
138 G4VSolid* msolid, DivisionType divType )
139 : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType )
140{
141
143 SetType( "DivisionPolyhedraRho" );
144
146 G4PolyhedraHistorical* original_pars = msol->GetOriginalParameters();
147
148 if( divType == DivWIDTH )
149 {
150 fnDiv = CalculateNDiv( original_pars->Rmax[0]
151 - original_pars->Rmin[0], width, offset );
152 }
153 else if( divType == DivNDIV )
154 {
155 fwidth = CalculateWidth( original_pars->Rmax[0]
156 - original_pars->Rmin[0], nDiv, offset );
157 }
158
159#ifdef G4DIVDEBUG
160 if( verbose >= 1 )
161 {
162 G4cout << " G4ParameterisationPolyhedraRho - # divisions " << fnDiv
163 << " = " << nDiv << G4endl
164 << " Offset " << foffset << " = " << offset << G4endl
165 << " Width " << fwidth << " = " << width << G4endl;
166 }
167#endif
168}
169
170//------------------------------------------------------------------------
172{
173}
174
175//---------------------------------------------------------------------
177{
179
181
183 {
184 std::ostringstream message;
185 message << "In solid " << msol->GetName() << G4endl
186 << "Division along R will be done with a width "
187 << "different for each solid section." << G4endl
188 << "WIDTH will not be used !";
189 G4Exception("G4ParameterisationPolyhedraRho::CheckParametersValidity()",
190 "GeomDiv1001", JustWarning, message);
191 }
192 if( foffset != 0. )
193 {
194 std::ostringstream message;
195 message << "In solid " << msol->GetName() << G4endl
196 << "Division along R will be done with a width "
197 << "different for each solid section." << G4endl
198 << "OFFSET will not be used !";
199 G4Exception("G4ParameterisationPolyhedraRho::CheckParametersValidity()",
200 "GeomDiv1001", JustWarning, message);
201 }
202}
203
204//------------------------------------------------------------------------
206{
208 G4PolyhedraHistorical* original_pars = msol->GetOriginalParameters();
209 return original_pars->Rmax[0] - original_pars->Rmin[0];
210}
211
212//--------------------------------------------------------------------------
213void
215ComputeTransformation( const G4int, G4VPhysicalVolume* physVol ) const
216{
217 //----- translation
218 G4ThreeVector origin(0.,0.,0.);
219
220 //----- set translation
221 physVol->SetTranslation( origin );
222
223 //----- calculate rotation matrix: unit
224
225#ifdef G4DIVDEBUG
226 if( verbose >= 2 )
227 {
228 G4cout << " G4ParameterisationPolyhedraRho " << G4endl
229 << " foffset: " << foffset/CLHEP::deg
230 << " - fwidth: " << fwidth/CLHEP::deg << G4endl;
231 }
232#endif
233
234 ChangeRotMatrix( physVol );
235
236#ifdef G4DIVDEBUG
237 if( verbose >= 2 )
238 {
239 G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraRho "
240 << G4endl
241 << " Position: " << origin
242 << " - Width: " << fwidth
243 << " - Axis: " << faxis << G4endl;
244 }
245#endif
246}
247
248//--------------------------------------------------------------------------
249void
251ComputeDimensions( G4Polyhedra& phedra, const G4int copyNo,
252 const G4VPhysicalVolume* ) const
253{
255
256 G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters();
257 G4PolyhedraHistorical origparam( *origparamMother );
258 G4int nZplanes = origparamMother->Num_z_planes;
259
260 G4double width = 0.;
261 for( G4int ii = 0; ii < nZplanes; ++ii )
262 {
263 width = CalculateWidth( origparamMother->Rmax[ii]
264 - origparamMother->Rmin[ii], fnDiv, foffset );
265 origparam.Rmin[ii] = origparamMother->Rmin[ii]+foffset+width*copyNo;
266 origparam.Rmax[ii] = origparamMother->Rmin[ii]+foffset+width*(copyNo+1);
267 }
268
269 phedra.SetOriginalParameters(&origparam); // copy values & transfer pointers
270 phedra.Reset(); // reset to new solid parameters
271
272#ifdef G4DIVDEBUG
273 if( verbose >= -2 )
274 {
275 G4cout << "G4ParameterisationPolyhedraRho::ComputeDimensions()" << G4endl
276 << "-- Parametrised phedra copy-number: " << copyNo << G4endl;
277 phedra.DumpInfo();
278 }
279#endif
280}
281
282//--------------------------------------------------------------------------
285 G4double width, G4double offset,
286 G4VSolid* msolid, DivisionType divType )
287 : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType )
288{
290 SetType( "DivisionPolyhedraPhi" );
291
293 G4double deltaPhi = msol->GetEndPhi() - msol->GetStartPhi();
294
295 if( divType == DivWIDTH )
296 {
297 fnDiv = msol->GetNumSide();
298 }
299
300 fwidth = CalculateWidth( deltaPhi, fnDiv, 0.0 );
301
302#ifdef G4DIVDEBUG
303 if( verbose >= 1 )
304 {
305 G4cout << " G4ParameterisationPolyhedraPhi - # divisions " << fnDiv
306 << " = " << nDiv << G4endl
307 << " Offset " << foffset << " = " << offset << G4endl
308 << " Width " << fwidth << " = " << width << G4endl;
309 }
310#endif
311}
312
313//------------------------------------------------------------------------
315{
316}
317
318//------------------------------------------------------------------------
320{
322 return msol->GetEndPhi() - msol->GetStartPhi();
323}
324
325//---------------------------------------------------------------------
327{
329
331
333 {
334 std::ostringstream message;
335 message << "In solid " << msol->GetName() << G4endl
336 << " Division along PHI will be done splitting "
337 << "in the defined numSide." << G4endl
338 << "WIDTH will not be used !";
339 G4Exception("G4ParameterisationPolyhedraPhi::CheckParametersValidity()",
340 "GeomDiv1001", JustWarning, message);
341 }
342 if( foffset != 0. )
343 {
344 std::ostringstream message;
345 message << "In solid " << msol->GetName() << G4endl
346 << "Division along PHI will be done splitting "
347 << "in the defined numSide." << G4endl
348 << "OFFSET will not be used !";
349 G4Exception("G4ParameterisationPolyhedraPhi::CheckParametersValidity()",
350 "GeomDiv1001", JustWarning, message);
351 }
352
353 G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters();
354
355 if( origparamMother->numSide != fnDiv && fDivisionType != DivWIDTH)
356 {
357 std::ostringstream message;
358 message << "Configuration not supported." << G4endl
359 << "Division along PHI will be done splitting in the defined"
360 << G4endl
361 << "numSide, i.e, the number of division would be :"
362 << origparamMother->numSide << " instead of " << fnDiv << " !";
363 G4Exception("G4ParameterisationPolyhedraPhi::CheckParametersValidity()",
364 "GeomDiv0001", FatalException, message);
365 }
366}
367
368//--------------------------------------------------------------------------
369void
371ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol ) const
372{
373 //----- translation
374 G4ThreeVector origin(0.,0.,0.);
375 //----- set translation
376 physVol->SetTranslation( origin );
377
378 //----- calculate rotation matrix (so that all volumes point to the centre)
379 G4double posi = copyNo*fwidth;
380
381#ifdef G4DIVDEBUG
382 if( verbose >= 2 )
383 {
384 G4cout << " G4ParameterisationPolyhedraPhi - position: " << posi/CLHEP::deg
385 << G4endl
386 << " copyNo: " << copyNo
387 << " - fwidth: " << fwidth/CLHEP::deg << G4endl;
388 }
389#endif
390
391 ChangeRotMatrix( physVol, -posi );
392
393#ifdef G4DIVDEBUG
394 if( verbose >= 2 )
395 {
396 G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraPhi "
397 << copyNo << G4endl
398 << " Position: " << origin << " - Width: " << fwidth
399 << " - Axis: " << faxis << G4endl;
400 }
401#endif
402}
403
404//--------------------------------------------------------------------------
405void
407ComputeDimensions( G4Polyhedra& phedra, const G4int,
408 const G4VPhysicalVolume* ) const
409{
411
412 G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters();
413 G4PolyhedraHistorical origparam( *origparamMother );
414
415 origparam.numSide = 1;
416 origparam.Start_angle = origparamMother->Start_angle;
417 origparam.Opening_angle = fwidth;
418
419 phedra.SetOriginalParameters(&origparam); // copy values & transfer pointers
420 phedra.Reset(); // reset to new solid parameters
421
422#ifdef G4DIVDEBUG
423 if( verbose >= 2 )
424 {
425 G4cout << "G4ParameterisationPolyhedraPhi::ComputeDimensions():" << G4endl;
426 phedra.DumpInfo();
427 }
428#endif
429}
430
431//--------------------------------------------------------------------------
434 G4double width, G4double offset,
435 G4VSolid* msolid, DivisionType divType )
436 : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType ),
437 fOrigParamMother(((G4Polyhedra*)fmotherSolid)->GetOriginalParameters())
438{
440 SetType( "DivisionPolyhedraZ" );
441
442 if( divType == DivWIDTH )
443 {
444 fnDiv =
445 CalculateNDiv(fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
446 - fOrigParamMother->Z_values[0] , width, offset);
447 }
448 else if( divType == DivNDIV )
449 {
450 fwidth =
451 CalculateNDiv(fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
452 - fOrigParamMother->Z_values[0] , nDiv, offset);
453 }
454
455#ifdef G4DIVDEBUG
456 if( verbose >= 1 )
457 {
458 G4cout << " G4ParameterisationPolyhedraZ - # divisions " << fnDiv << " = "
459 << nDiv << G4endl
460 << " Offset " << foffset << " = " << offset << G4endl
461 << " Width " << fwidth << " = " << width << G4endl;
462 }
463#endif
464}
465
466//---------------------------------------------------------------------
468{
469}
470
471//------------------------------------------------------------------------
472G4double G4ParameterisationPolyhedraZ::GetR(G4double z,
473 G4double z1, G4double r1,
474 G4double z2, G4double r2) const
475{
476 // Linear parameterisation:
477 // r = az + b
478 // a = (r1 - r2)/(z1-z2)
479 // b = r1 - a*z1
480
481 return (r1-r2)/(z1-z2)*z + ( r1 - (r1-r2)/(z1-z2)*z1 ) ;
482}
483
484//------------------------------------------------------------------------
485G4double G4ParameterisationPolyhedraZ::GetRmin(G4double z, G4int nseg) const
486{
487// Get Rmin in the given z position for the given polyhedra segment
488
489 return GetR(z,
490 fOrigParamMother->Z_values[nseg],
491 fOrigParamMother->Rmin[nseg],
492 fOrigParamMother->Z_values[nseg+1],
493 fOrigParamMother->Rmin[nseg+1]);
494}
495
496//------------------------------------------------------------------------
497G4double G4ParameterisationPolyhedraZ::GetRmax(G4double z, G4int nseg) const
498{
499// Get Rmax in the given z position for the given polyhedra segment
500
501 return GetR(z,
502 fOrigParamMother->Z_values[nseg],
503 fOrigParamMother->Rmax[nseg],
504 fOrigParamMother->Z_values[nseg+1],
505 fOrigParamMother->Rmax[nseg+1]);
506}
507
508//------------------------------------------------------------------------
510{
511 return std::abs(fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
512 -fOrigParamMother->Z_values[0]);
513}
514
515//---------------------------------------------------------------------
517{
519
520 // Division will be following the mother polyhedra segments
521 //
522 if( fDivisionType == DivNDIV )
523 {
524 if( fOrigParamMother->Num_z_planes-1 != fnDiv )
525 {
526 std::ostringstream message;
527 message << "Configuration not supported." << G4endl
528 << "Division along Z will be done splitting in the defined"
529 << G4endl
530 << "Z planes, i.e, the number of division would be :"
531 << fOrigParamMother->Num_z_planes-1 << " instead of "
532 << fnDiv << " !";
533 G4Exception("G4ParameterisationPolyhedraZ::CheckParametersValidity()",
534 "GeomDiv0001", FatalException, message);
535 }
536 }
537
538 // Division will be done within one polyhedra segment
539 // with applying given width and offset
540 //
542 {
543 // Check if divided region does not span over more
544 // than one z segment
545
546 G4int isegstart = -1; // number of the segment containing start position
547 G4int isegend = -1; // number of the segment containing end position
548
549 if ( !fReflectedSolid )
550 {
551 // The start/end position of the divided region
552 //
553 G4double zstart = fOrigParamMother->Z_values[0] + foffset;
554 G4double zend = fOrigParamMother->Z_values[0]
555 + foffset + fnDiv*fwidth;
556
557 G4int counter = 0;
558 while ( isegend < 0 && counter < fOrigParamMother->Num_z_planes-1 )
559 {
560 // first segment
561 if ( zstart >= fOrigParamMother->Z_values[counter] &&
562 zstart < fOrigParamMother->Z_values[counter+1] )
563 {
564 isegstart = counter;
565 }
566 // last segment
567 if ( zend > fOrigParamMother->Z_values[counter] &&
568 zend <= fOrigParamMother->Z_values[counter+1] )
569 {
570 isegend = counter;
571 }
572 ++counter;
573 } // Loop checking, 06.08.2015, G.Cosmo
574 }
575 else
576 {
577 // The start/end position of the divided region
578 //
579 G4double zstart = fOrigParamMother->Z_values[0] - foffset;
580 G4double zend = fOrigParamMother->Z_values[0]
581 - (foffset + fnDiv* fwidth);
582
583 G4int counter = 0;
584 while ( isegend < 0 && counter < fOrigParamMother->Num_z_planes-1 )
585 {
586 // first segment
587 if ( zstart <= fOrigParamMother->Z_values[counter] &&
588 zstart > fOrigParamMother->Z_values[counter+1] )
589 {
590 isegstart = counter;
591 }
592 // last segment
593 if ( zend < fOrigParamMother->Z_values[counter] &&
594 zend >= fOrigParamMother->Z_values[counter+1] )
595 {
596 isegend = counter;
597 }
598 ++counter;
599 } // Loop checking, 06.08.2015, G.Cosmo
600 }
601
602 if ( isegstart != isegend )
603 {
604 std::ostringstream message;
605 message << "Configuration not supported." << G4endl
606 << "Division with user defined width." << G4endl
607 << "Solid " << fmotherSolid->GetName() << G4endl
608 << "Divided region is not between two Z planes.";
609 G4Exception("G4ParameterisationPolyhedraZ::CheckParametersValidity()",
610 "GeomDiv0001", FatalException, message);
611 }
612
613 fNSegment = isegstart;
614 }
615}
616
617//---------------------------------------------------------------------
618void
620ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol) const
621{
622 G4double posi;
623 if ( fDivisionType == DivNDIV )
624 {
625 // The position of the centre of copyNo-th mother polycone segment
626
627 posi = ( fOrigParamMother->Z_values[copyNo]
628 + fOrigParamMother->Z_values[copyNo+1])/2;
629 physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
630 }
631
633 {
634 // The position of the centre of copyNo-th division
635
636 posi = fOrigParamMother->Z_values[0];
637
638 if ( !fReflectedSolid )
639 posi += foffset + (2*copyNo + 1) * fwidth/2.;
640 else
641 posi -= foffset + (2*copyNo + 1) * fwidth/2.;
642
643 physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
644 }
645
646 //----- calculate rotation matrix: unit
647
648#ifdef G4DIVDEBUG
649 if( verbose >= 2 )
650 {
651 G4cout << " G4ParameterisationPolyhedraZ - position: " << posi << G4endl
652 << " copyNo: " << copyNo << " - foffset: " << foffset/CLHEP::deg
653 << " - fwidth: " << fwidth/CLHEP::deg << G4endl;
654 }
655#endif
656
657 ChangeRotMatrix( physVol );
658
659#ifdef G4DIVDEBUG
660 if( verbose >= 2 )
661 {
662 G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraZ "
663 << copyNo << G4endl
664 << " Position: (0,0,0) - Width: " << fwidth
665 << " - Axis: " << faxis << G4endl;
666 }
667#endif
668}
669
670//---------------------------------------------------------------------
671void
673ComputeDimensions( G4Polyhedra& phedra, const G4int copyNo,
674 const G4VPhysicalVolume* ) const
675{
676 // Define division solid
677 //
678 G4PolyhedraHistorical origparam;
679 G4int nz = 2;
680 origparam.Num_z_planes = nz;
681 origparam.numSide = fOrigParamMother->numSide;
682 origparam.Start_angle = fOrigParamMother->Start_angle;
683 origparam.Opening_angle = fOrigParamMother->Opening_angle;
684
685 // Define division solid z sections
686 //
687 origparam.Z_values = new G4double[nz];
688 origparam.Rmin = new G4double[nz];
689 origparam.Rmax = new G4double[nz];
690 origparam.Z_values[0] = - fwidth/2.;
691 origparam.Z_values[1] = fwidth/2.;
692
693 if ( fDivisionType == DivNDIV )
694 {
695 // The position of the centre of copyNo-th mother polycone segment
696 //
697 G4double posi = ( fOrigParamMother->Z_values[copyNo]
698 + fOrigParamMother->Z_values[copyNo+1])/2;
699
700 origparam.Z_values[0] = fOrigParamMother->Z_values[copyNo] - posi;
701 origparam.Z_values[1] = fOrigParamMother->Z_values[copyNo+1] - posi;
702 origparam.Rmin[0] = fOrigParamMother->Rmin[copyNo];
703 origparam.Rmin[1] = fOrigParamMother->Rmin[copyNo+1];
704 origparam.Rmax[0] = fOrigParamMother->Rmax[copyNo];
705 origparam.Rmax[1] = fOrigParamMother->Rmax[copyNo+1];
706 }
707
709 {
710 if ( !fReflectedSolid )
711 {
712 origparam.Z_values[0] = -fwidth/2.;
713 origparam.Z_values[1] = fwidth/2.;
714
715 // The position of the centre of copyNo-th division
716 //
717 G4double posi = fOrigParamMother->Z_values[0]
718 + foffset + (2*copyNo + 1) * fwidth/2.;
719
720 // The first and last z sides z values
721 G4double zstart = posi - fwidth/2.;
722 G4double zend = posi + fwidth/2.;
723 origparam.Rmin[0] = GetRmin(zstart, fNSegment);
724 origparam.Rmax[0] = GetRmax(zstart, fNSegment);
725 origparam.Rmin[1] = GetRmin(zend, fNSegment);
726 origparam.Rmax[1] = GetRmax(zend, fNSegment);
727 }
728 else
729 {
730 origparam.Z_values[0] = fwidth/2.;
731 origparam.Z_values[1] = -fwidth/2.;
732
733 // The position of the centre of copyNo-th division
734 //
735 G4double posi = fOrigParamMother->Z_values[0]
736 - ( foffset + (2*copyNo + 1) * fwidth/2.);
737
738 // The first and last z sides z values
739 //
740 G4double zstart = posi + fwidth/2.;
741 G4double zend = posi - fwidth/2.;
742 origparam.Rmin[0] = GetRmin(zstart, fNSegment);
743 origparam.Rmax[0] = GetRmax(zstart, fNSegment);
744 origparam.Rmin[1] = GetRmin(zend, fNSegment);
745 origparam.Rmax[1] = GetRmax(zend, fNSegment);
746 }
747
748 // It can happen due to rounding errors
749 //
750 if ( origparam.Rmin[0] < 0.0 ) origparam.Rmin[0] = 0.0;
751 if ( origparam.Rmin[nz-1] < 0.0 ) origparam.Rmin[1] = 0.0;
752 }
753
754 phedra.SetOriginalParameters(&origparam); // copy values & transfer pointers
755 phedra.Reset(); // reset to new solid parameters
756
757#ifdef G4DIVDEBUG
758 if( verbose >= 2 )
759 {
760 G4cout << "G4ParameterisationPolyhedraZ::ComputeDimensions()" << G4endl
761 << "-- Parametrised phedra copy-number: " << copyNo << G4endl;
762 phedra.DumpInfo();
763 }
764#endif
765}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4GeometryTolerance * GetInstance()
G4ParameterisationPolyhedraPhi(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
void ComputeDimensions(G4Polyhedra &phedra, const G4int copyNo, const G4VPhysicalVolume *physVol) const
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
void ComputeDimensions(G4Polyhedra &phedra, const G4int copyNo, const G4VPhysicalVolume *physVol) const
G4ParameterisationPolyhedraRho(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
G4ParameterisationPolyhedraZ(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
void ComputeDimensions(G4Polyhedra &phedra, const G4int copyNo, const G4VPhysicalVolume *physVol) const
G4double GetEndPhi() const
void SetOriginalParameters(G4PolyhedraHistorical *pars)
G4int GetNumSide() const
G4PolyhedraHistorical * GetOriginalParameters() const
G4bool Reset()
Definition: G4Polyhedra.cc:463
G4double GetStartPhi() const
G4bool IsGeneric() const
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.0) const
G4VParameterisationPolyhedra(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
void SetTranslation(const G4ThreeVector &v)
G4String GetName() const
void DumpInfo() const
virtual G4GeometryType GetEntityType() const =0
EAxis
Definition: geomdefs.hh:54