Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParameterisationPolycone.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// G4ParameterisationPolycone[Rho/Phi/Z] implementation
27//
28// 26.05.03 - P.Arce, Initial version
29// 08.04.04 - I.Hrivnacova, Implemented reflection
30//---------------------------------------------------------------------
31
33
34#include <iomanip>
35#include "G4ThreeVector.hh"
36#include "G4RotationMatrix.hh"
37#include "G4VPhysicalVolume.hh"
38#include "G4LogicalVolume.hh"
39#include "G4ReflectedSolid.hh"
40
41//-----------------------------------------------------------------------
44 G4double offset, G4VSolid* msolid,
45 DivisionType divType )
46 : G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid )
47{
48 /*#ifdef G4MULTITHREADED
49 std::ostringstream message;
50 message << "Divisions for G4Polycone currently NOT supported in MT-mode."
51 << G4endl
52 << "Sorry! Solid: " << msolid->GetName();
53 G4Exception("G4VParameterisationPolycone::G4VParameterisationPolycone()",
54 "GeomDiv0001", FatalException, message);
55 #endif */
56 auto msol = (G4Polycone*)(msolid);
57 if (msolid->GetEntityType() == "G4ReflectedSolid")
58 {
59 // Get constituent solid
60 //
61 G4VSolid* mConstituentSolid
62 = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid();
63 msol = (G4Polycone*)(mConstituentSolid);
64
65 // Get parameters
66 //
67 G4int nofZplanes = msol->GetOriginalParameters()->Num_z_planes;
68 G4double* zValues = msol->GetOriginalParameters()->Z_values;
69 G4double* rminValues = msol->GetOriginalParameters()->Rmin;
70 G4double* rmaxValues = msol->GetOriginalParameters()->Rmax;
71
72 // Invert z values
73 //
74 auto zValuesRefl = new G4double[nofZplanes];
75 for (G4int i=0; i<nofZplanes; ++i) { zValuesRefl[i] = - zValues[i]; }
76
77 auto newSolid
78 = new G4Polycone(msol->GetName(),
79 msol->GetStartPhi(),
80 msol->GetEndPhi() - msol->GetStartPhi(),
81 nofZplanes, zValuesRefl, rminValues, rmaxValues);
82
83 delete [] zValuesRefl;
84
85 msol = newSolid;
86 fmotherSolid = newSolid;
87 fReflectedSolid = true;
88 fDeleteSolid = true;
89 }
90}
91
92//---------------------------------------------------------------------
94
95//---------------------------------------------------------------------
98 G4double width, G4double offset,
99 G4VSolid* msolid, DivisionType divType )
100 : G4VParameterisationPolycone( axis, nDiv, width, offset, msolid, divType )
101{
103 SetType( "DivisionPolyconeRho" );
104
105 auto msol = (G4Polycone*)(fmotherSolid);
106 G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
107
108 if( divType == DivWIDTH )
109 {
110 fnDiv = CalculateNDiv( origparamMother->Rmax[0]
111 - origparamMother->Rmin[0], width, offset );
112 }
113 else if( divType == DivNDIV )
114 {
115 fwidth = CalculateWidth( origparamMother->Rmax[0]
116 - origparamMother->Rmin[0], nDiv, offset );
117 }
118
119#ifdef G4DIVDEBUG
120 if( verbose >= -1 )
121 {
122 G4cout << " G4ParameterisationPolyconeRho - # divisions " << fnDiv
123 << " = " << nDiv << G4endl
124 << " Offset " << foffset << " = " << offset << G4endl
125 << " Width " << fwidth << " = " << width << G4endl;
126 }
127#endif
128}
129
130//---------------------------------------------------------------------
132
133//---------------------------------------------------------------------
135{
137
138 auto msol = (G4Polycone*)(fmotherSolid);
139
141 {
142 std::ostringstream message;
143 message << "In solid " << msol->GetName() << G4endl
144 << "Division along R will be done with a width "
145 << "different for each solid section." << G4endl
146 << "WIDTH will not be used !";
147 G4Exception("G4VParameterisationPolycone::CheckParametersValidity()",
148 "GeomDiv1001", JustWarning, message);
149 }
150 if( foffset != 0. )
151 {
152 std::ostringstream message;
153 message << "In solid " << msol->GetName() << G4endl
154 << "Division along R will be done with a width "
155 << "different for each solid section." << G4endl
156 << "OFFSET will not be used !";
157 G4Exception("G4VParameterisationPolycone::CheckParametersValidity()",
158 "GeomDiv1001", JustWarning, message);
159 }
160}
161
162//------------------------------------------------------------------------
164{
165 auto msol = (G4Polycone*)(fmotherSolid);
166 G4PolyconeHistorical* original_pars = msol->GetOriginalParameters();
167 return original_pars->Rmax[0] - original_pars->Rmin[0];
168}
169
170
171//---------------------------------------------------------------------
172void
174ComputeTransformation( const G4int, G4VPhysicalVolume* physVol ) const
175{
176 //----- translation
177 G4ThreeVector origin(0.,0.,0.);
178 //----- set translation
179 physVol->SetTranslation( origin );
180
181 //----- calculate rotation matrix: unit
182
183#ifdef G4DIVDEBUG
184 if( verbose >= 2 )
185 {
186 G4cout << " G4ParameterisationPolyconeRho " << G4endl
187 << " foffset: " << foffset
188 << " - fwidth: " << fwidth << G4endl;
189 }
190#endif
191
192 ChangeRotMatrix( physVol );
193
194#ifdef G4DIVDEBUG
195 if( verbose >= 2 )
196 {
197 G4cout << std::setprecision(8) << " G4ParameterisationPolyconeRho "
198 << G4endl
199 << " Position: (0,0,0)"
200 << " - Width: " << fwidth/CLHEP::deg
201 << " - Axis: " << faxis << G4endl;
202 }
203#endif
204}
205
206//---------------------------------------------------------------------
207void
209ComputeDimensions( G4Polycone& pcone, const G4int copyNo,
210 const G4VPhysicalVolume* ) const
211{
212 auto msol = (G4Polycone*)(fmotherSolid);
213
214 G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
215 G4PolyconeHistorical origparam( *origparamMother );
216 G4int nZplanes = origparamMother->Num_z_planes;
217
218 G4double width = 0.;
219 for( G4int ii = 0; ii < nZplanes; ++ii )
220 {
221 width = CalculateWidth( origparamMother->Rmax[ii]
222 - origparamMother->Rmin[ii], fnDiv, foffset );
223 origparam.Rmin[ii] = origparamMother->Rmin[ii]+foffset+width*copyNo;
224 origparam.Rmax[ii] = origparamMother->Rmin[ii]+foffset+width*(copyNo+1);
225 }
226
227 pcone.SetOriginalParameters(&origparam); // copy values & transfer pointers
228 pcone.Reset(); // reset to new solid parameters
229
230#ifdef G4DIVDEBUG
231 if( verbose >= -2 )
232 {
233 G4cout << "G4ParameterisationPolyconeRho::ComputeDimensions()" << G4endl
234 << "-- Parametrised pcone copy-number: " << copyNo << G4endl;
235 pcone.DumpInfo();
236 }
237#endif
238}
239
240//---------------------------------------------------------------------
243 G4double width, G4double offset,
244 G4VSolid* msolid, DivisionType divType )
245 : G4VParameterisationPolycone( axis, nDiv, width, offset, msolid, divType )
246{
248 SetType( "DivisionPolyconePhi" );
249
250 auto msol = (G4Polycone*)(fmotherSolid);
251 G4double deltaPhi = msol->GetEndPhi() - msol->GetStartPhi();
252
253 if( divType == DivWIDTH )
254 {
255 fnDiv = CalculateNDiv( deltaPhi, width, offset );
256 }
257 else if( divType == DivNDIV )
258 {
259 fwidth = CalculateWidth( deltaPhi, nDiv, offset );
260 }
261
262#ifdef G4DIVDEBUG
263 if( verbose >= 1 )
264 {
265 G4cout << " G4ParameterisationPolyconePhi - # divisions " << fnDiv
266 << " = " << nDiv << G4endl
267 << " Offset " << foffset/CLHEP::deg << " = " << offset/CLHEP::deg << G4endl
268 << " Width " << fwidth/CLHEP::deg << " = " << width/CLHEP::deg << G4endl;
269 }
270#endif
271}
272
273//---------------------------------------------------------------------
275
276//------------------------------------------------------------------------
278{
279 auto msol = (G4Polycone*)(fmotherSolid);
280 return msol->GetEndPhi() - msol->GetStartPhi();
281}
282
283//---------------------------------------------------------------------
284void
286ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol ) const
287{
288 //----- translation
289 G4ThreeVector origin(0.,0.,0.);
290 //----- set translation
291 physVol->SetTranslation( origin );
292
293 //----- calculate rotation matrix (so that all volumes point to the centre)
294 G4double posi = foffset + copyNo*fwidth;
295
296#ifdef G4DIVDEBUG
297 if( verbose >= 2 )
298 {
299 G4cout << " G4ParameterisationPolyconePhi - position: " << posi/CLHEP::deg
300 << G4endl
301 << " copyNo: " << copyNo << " - foffset: " << foffset/CLHEP::deg
302 << " - fwidth: " << fwidth/CLHEP::deg << G4endl;
303 }
304#endif
305
306 ChangeRotMatrix( physVol, -posi );
307
308#ifdef G4DIVDEBUG
309 if( verbose >= 2 )
310 {
311 G4cout << std::setprecision(8) << " G4ParameterisationPolyconePhi "
312 << copyNo << G4endl
313 << " Position: (0,0,0) - Width: " << fwidth
314 << " - Axis: " << faxis << G4endl;
315 }
316#endif
317}
318
319//---------------------------------------------------------------------
320void
322ComputeDimensions( G4Polycone& pcone, const G4int,
323 const G4VPhysicalVolume* ) const
324{
325 auto msol = (G4Polycone*)(fmotherSolid);
326
327 G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
328 G4PolyconeHistorical origparam( *origparamMother );
329 origparam.Start_angle = origparamMother->Start_angle;
330 origparam.Opening_angle = fwidth;
331
332 pcone.SetOriginalParameters(&origparam); // copy values & transfer pointers
333 pcone.Reset(); // reset to new solid parameters
334
335#ifdef G4DIVDEBUG
336 if( verbose >= 2 )
337 {
338 G4cout << "G4ParameterisationPolyconePhi::ComputeDimensions():" << G4endl;
339 pcone.DumpInfo();
340 }
341#endif
342}
343
344//---------------------------------------------------------------------
347 G4double width, G4double offset,
348 G4VSolid* msolid, DivisionType divType)
349 : G4VParameterisationPolycone( axis, nDiv, width, offset, msolid, divType ),
350 fOrigParamMother(((G4Polycone*)fmotherSolid)->GetOriginalParameters())
351{
352
354 SetType( "DivisionPolyconeZ" );
355
356 if( divType == DivWIDTH )
357 {
358 fnDiv =
359 CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
360 - fOrigParamMother->Z_values[0] , width, offset );
361 }
362 else if( divType == DivNDIV )
363 {
364 fwidth =
365 CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
366 - fOrigParamMother->Z_values[0] , nDiv, offset );
367 }
368
369#ifdef G4DIVDEBUG
370 if( verbose >= 1 )
371 {
372 G4cout << " G4ParameterisationPolyconeZ - # divisions " << fnDiv << " = "
373 << nDiv << G4endl
374 << " Offset " << foffset << " = " << offset << G4endl
375 << " Width " << fwidth << " = " << width << G4endl;
376 }
377#endif
378}
379
380//---------------------------------------------------------------------
382
383//------------------------------------------------------------------------
384G4double G4ParameterisationPolyconeZ::GetR(G4double z,
385 G4double z1, G4double r1,
386 G4double z2, G4double r2) const
387{
388 // Linear parameterisation:
389 // r = az + b
390 // a = (r1 - r2)/(z1-z2)
391 // b = r1 - a*z1
392
393 return (r1-r2)/(z1-z2)*z + ( r1 - (r1-r2)/(z1-z2)*z1 ) ;
394}
395
396//------------------------------------------------------------------------
397G4double G4ParameterisationPolyconeZ::GetRmin(G4double z, G4int nseg) const
398{
399// Get Rmin in the given z position for the given polycone segment
400
401 return GetR(z,
402 fOrigParamMother->Z_values[nseg],
403 fOrigParamMother->Rmin[nseg],
404 fOrigParamMother->Z_values[nseg+1],
405 fOrigParamMother->Rmin[nseg+1]);
406}
407
408//------------------------------------------------------------------------
409G4double G4ParameterisationPolyconeZ::GetRmax(G4double z, G4int nseg) const
410{
411// Get Rmax in the given z position for the given polycone segment
412
413 return GetR(z,
414 fOrigParamMother->Z_values[nseg],
415 fOrigParamMother->Rmax[nseg],
416 fOrigParamMother->Z_values[nseg+1],
417 fOrigParamMother->Rmax[nseg+1]);
418}
419
420//------------------------------------------------------------------------
422{
423 return std::abs (fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
424 -fOrigParamMother->Z_values[0]);
425}
426
427//---------------------------------------------------------------------
429{
431
432 // Division will be following the mother polycone segments
433 //
434 if( fDivisionType == DivNDIV )
435 {
436 if( fnDiv > fOrigParamMother->Num_z_planes-1 )
437 {
438 std::ostringstream error;
439 error << "Configuration not supported." << G4endl
440 << "Division along Z will be done by splitting in the defined"
441 << G4endl
442 << "Z planes, i.e, the number of division would be: "
443 << fOrigParamMother->Num_z_planes-1
444 << ", instead of: " << fnDiv << " !";
445 G4Exception("G4ParameterisationPolyconeZ::CheckParametersValidity()",
446 "GeomDiv0001", FatalException, error);
447 }
448 }
449
450 // Division will be done within one polycone segment
451 // with applying given width and offset
452 //
454 {
455 // Check if divided region does not span over more
456 // than one z segment
457
458 G4int isegstart = -1; // number of the segment containing start position
459 G4int isegend = -1; // number of the segment containing end position
460
461 if ( !fReflectedSolid )
462 {
463 // The start/end position of the divided region
464 //
465 G4double zstart
466 = fOrigParamMother->Z_values[0] + foffset;
467 G4double zend
468 = fOrigParamMother->Z_values[0] + foffset + fnDiv* fwidth;
469
470 G4int counter = 0;
471 while ( isegend < 0 && counter < fOrigParamMother->Num_z_planes-1 )
472 {
473 // first segment
474 if ( zstart >= fOrigParamMother->Z_values[counter] &&
475 zstart < fOrigParamMother->Z_values[counter+1] )
476 {
477 isegstart = counter;
478 }
479 // last segment
480 if ( zend > fOrigParamMother->Z_values[counter] &&
481 zend <= fOrigParamMother->Z_values[counter+1] )
482 {
483 isegend = counter;
484 }
485 ++counter;
486 } // Loop checking, 06.08.2015, G.Cosmo
487 }
488 else
489 {
490 // The start/end position of the divided region
491 //
492 G4double zstart
493 = fOrigParamMother->Z_values[0] - foffset;
494 G4double zend
495 = fOrigParamMother->Z_values[0] - ( foffset + fnDiv* fwidth);
496
497 G4int counter = 0;
498 while ( isegend < 0 && counter < fOrigParamMother->Num_z_planes-1 )
499 {
500 // first segment
501 if ( zstart <= fOrigParamMother->Z_values[counter] &&
502 zstart > fOrigParamMother->Z_values[counter+1] )
503 {
504 isegstart = counter;
505 }
506 // last segment
507 if ( zend < fOrigParamMother->Z_values[counter] &&
508 zend >= fOrigParamMother->Z_values[counter+1] )
509 {
510 isegend = counter;
511 }
512 ++counter;
513 } // Loop checking, 06.08.2015, G.Cosmo
514 }
515
516
517 if ( isegstart != isegend )
518 {
519 std::ostringstream message;
520 message << "Condiguration not supported." << G4endl
521 << "Division with user defined width." << G4endl
522 << "Solid " << fmotherSolid->GetName() << G4endl
523 << "Divided region is not between two z planes.";
524 G4Exception("G4ParameterisationPolyconeZ::CheckParametersValidity()",
525 "GeomDiv0001", FatalException, message);
526 }
527
528 fNSegment = isegstart;
529 }
530}
531
532//---------------------------------------------------------------------
533void
535ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol) const
536{
537 G4double posi = 0.;
538 if ( fDivisionType == DivNDIV )
539 {
540 // The position of the centre of copyNo-th mother polycone segment
541 //
542 posi = ( fOrigParamMother->Z_values[copyNo]
543 + fOrigParamMother->Z_values[copyNo+1])/2;
544 physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
545 }
546
548 {
549 // The position of the centre of copyNo-th division
550 //
551 posi = fOrigParamMother->Z_values[0];
552
553 if ( !fReflectedSolid )
554 posi += foffset + (2*copyNo + 1) * fwidth/2.;
555 else
556 posi -= foffset + (2*copyNo + 1) * fwidth/2.;
557
558 physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
559 }
560
561 //----- calculate rotation matrix: unit
562
563#ifdef G4DIVDEBUG
564 if( verbose >= 2 )
565 {
566 G4cout << " G4ParameterisationPolyconeZ - position: " << posi << G4endl
567 << " copyNo: " << copyNo << " - foffset: " << foffset/CLHEP::deg
568 << " - fwidth: " << fwidth/CLHEP::deg << G4endl;
569 }
570#endif
571
572 ChangeRotMatrix( physVol );
573
574#ifdef G4DIVDEBUG
575 if( verbose >= 2 )
576 {
577 G4cout << std::setprecision(8) << " G4ParameterisationPolyconeZ "
578 << copyNo << G4endl
579 << " Position: (0,0,0) - Width: " << fwidth
580 << " - Axis: " << faxis << G4endl;
581 }
582#endif
583}
584
585//---------------------------------------------------------------------
586void
588ComputeDimensions( G4Polycone& pcone, const G4int copyNo,
589 const G4VPhysicalVolume* ) const
590{
591
592 // Define division solid
593 //
594 G4PolyconeHistorical origparam;
595 G4int nz = 2;
596 origparam.Num_z_planes = nz;
597 origparam.Start_angle = fOrigParamMother->Start_angle;
598 origparam.Opening_angle = fOrigParamMother->Opening_angle;
599
600 // Define division solid z sections
601 //
602 origparam.Z_values = new G4double[nz];
603 origparam.Rmin = new G4double[nz];
604 origparam.Rmax = new G4double[nz];
605
606 if ( fDivisionType == DivNDIV )
607 {
608 // The position of the centre of copyNo-th mother polycone segment
609 G4double posi = (fOrigParamMother->Z_values[copyNo]
610 + fOrigParamMother->Z_values[copyNo+1])/2;
611
612 origparam.Z_values[0] = fOrigParamMother->Z_values[copyNo] - posi;
613 origparam.Z_values[1] = fOrigParamMother->Z_values[copyNo+1] - posi;
614 origparam.Rmin[0] = fOrigParamMother->Rmin[copyNo];
615 origparam.Rmin[1] = fOrigParamMother->Rmin[copyNo+1];
616 origparam.Rmax[0] = fOrigParamMother->Rmax[copyNo];
617 origparam.Rmax[1] = fOrigParamMother->Rmax[copyNo+1];
618 }
619
621 {
622 if ( !fReflectedSolid )
623 {
624 origparam.Z_values[0] = - fwidth/2.;
625 origparam.Z_values[1] = fwidth/2.;
626
627 // The position of the centre of copyNo-th division
628 //
629 G4double posi = fOrigParamMother->Z_values[0]
630 + foffset + (2*copyNo + 1) * fwidth/2.;
631
632 // The first and last z sides z values
633 //
634 G4double zstart = posi - fwidth/2.;
635 G4double zend = posi + fwidth/2.;
636 origparam.Rmin[0] = GetRmin(zstart, fNSegment);
637 origparam.Rmax[0] = GetRmax(zstart, fNSegment);
638 origparam.Rmin[1] = GetRmin(zend, fNSegment);
639 origparam.Rmax[1] = GetRmax(zend, fNSegment);
640 }
641 else
642 {
643 origparam.Z_values[0] = fwidth/2.;
644 origparam.Z_values[1] = - fwidth/2.;
645
646 // The position of the centre of copyNo-th division
647 //
648 G4double posi = fOrigParamMother->Z_values[0]
649 - ( foffset + (2*copyNo + 1) * fwidth/2.);
650
651 // The first and last z sides z values
652 //
653 G4double zstart = posi + fwidth/2.;
654 G4double zend = posi - fwidth/2.;
655 origparam.Rmin[0] = GetRmin(zstart, fNSegment);
656 origparam.Rmax[0] = GetRmax(zstart, fNSegment);
657 origparam.Rmin[1] = GetRmin(zend, fNSegment);
658 origparam.Rmax[1] = GetRmax(zend, fNSegment);
659 }
660
661 // It can happen due to rounding errors
662 //
663 if ( origparam.Rmin[0] < 0.0 ) origparam.Rmin[0] = 0.0;
664 if ( origparam.Rmin[nz-1] < 0.0 ) origparam.Rmin[1] = 0.0;
665 }
666
667 pcone.SetOriginalParameters(&origparam); // copy values & transfer pointers
668 pcone.Reset(); // reset to new solid parameters
669
670#ifdef G4DIVDEBUG
671 if( verbose >= 2 )
672 {
673 G4cout << "G4ParameterisationPolyconeZ::ComputeDimensions()" << G4endl
674 << "-- Parametrised pcone copy-number: " << copyNo << G4endl;
675 pcone.DumpInfo();
676 }
677#endif
678}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const override
void ComputeDimensions(G4Polycone &pcone, const G4int copyNo, const G4VPhysicalVolume *physVol) const override
~G4ParameterisationPolyconePhi() override
G4ParameterisationPolyconePhi(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
void ComputeDimensions(G4Polycone &pcone, const G4int copyNo, const G4VPhysicalVolume *physVol) const override
G4ParameterisationPolyconeRho(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
~G4ParameterisationPolyconeRho() override
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const override
void ComputeDimensions(G4Polycone &pcone, const G4int copyNo, const G4VPhysicalVolume *physVol) const override
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const override
~G4ParameterisationPolyconeZ() override
G4ParameterisationPolyconeZ(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
G4double GetMaxParameter() const override
void SetOriginalParameters(G4PolyconeHistorical *pars)
G4bool Reset()
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
~G4VParameterisationPolycone() override
G4VParameterisationPolycone(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