Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParameterisationPolyconeRho Class Reference

#include <G4ParameterisationPolycone.hh>

+ Inheritance diagram for G4ParameterisationPolyconeRho:

Public Member Functions

 G4ParameterisationPolyconeRho (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationPolyconeRho () override
 
void CheckParametersValidity () override
 
G4double GetMaxParameter () const override
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const override
 
void ComputeDimensions (G4Polycone &pcone, const G4int copyNo, const G4VPhysicalVolume *physVol) const override
 
- Public Member Functions inherited from G4VParameterisationPolycone
 G4VParameterisationPolycone (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
 ~G4VParameterisationPolycone () override
 
- Public Member Functions inherited from G4VDivisionParameterisation
 G4VDivisionParameterisation (EAxis axis, G4int nDiv, G4double width, G4double offset, DivisionType divType, G4VSolid *motherSolid=nullptr)
 
 ~G4VDivisionParameterisation () override
 
G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *) override
 
const G4StringGetType () const
 
EAxis GetAxis () const
 
G4int GetNoDiv () const
 
G4double GetWidth () const
 
G4double GetOffset () const
 
G4VSolidGetMotherSolid () const
 
void SetType (const G4String &type)
 
G4int VolumeFirstCopyNo () const
 
void SetHalfGap (G4double hg)
 
G4double GetHalfGap () const
 
- Public Member Functions inherited from G4VPVParameterisation
 G4VPVParameterisation ()=default
 
virtual ~G4VPVParameterisation ()=default
 
virtual G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr)
 
virtual G4bool IsNested () const
 
virtual G4VVolumeMaterialScannerGetMaterialScanner ()
 

Additional Inherited Members

- Protected Member Functions inherited from G4VDivisionParameterisation
void ChangeRotMatrix (G4VPhysicalVolume *physVol, G4double rotZ=0.0) const
 
G4int CalculateNDiv (G4double motherDim, G4double width, G4double offset) const
 
G4double CalculateWidth (G4double motherDim, G4int nDiv, G4double offset) const
 
void CheckOffset (G4double maxPar)
 
void CheckNDivAndWidth (G4double maxPar)
 
G4double OffsetZ () const
 
- Protected Attributes inherited from G4VDivisionParameterisation
G4String ftype
 
EAxis faxis
 
G4int fnDiv = 0
 
G4double fwidth = 0.0
 
G4double foffset = 0.0
 
DivisionType fDivisionType
 
G4VSolidfmotherSolid = nullptr
 
G4bool fReflectedSolid = false
 
G4bool fDeleteSolid = false
 
G4int theVoluFirstCopyNo = 1
 
G4double kCarTolerance
 
G4double fhgap = 0.0
 
- Static Protected Attributes inherited from G4VDivisionParameterisation
static G4ThreadLocal G4RotationMatrixfRot = nullptr
 
static const G4int verbose = 5
 

Detailed Description

Definition at line 73 of file G4ParameterisationPolycone.hh.

Constructor & Destructor Documentation

◆ G4ParameterisationPolyconeRho()

G4ParameterisationPolyconeRho::G4ParameterisationPolyconeRho ( EAxis axis,
G4int nCopies,
G4double offset,
G4double step,
G4VSolid * motherSolid,
DivisionType divType )

Definition at line 96 of file G4ParameterisationPolycone.cc.

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}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
G4VParameterisationPolycone(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)

◆ ~G4ParameterisationPolyconeRho()

G4ParameterisationPolyconeRho::~G4ParameterisationPolyconeRho ( )
overridedefault

Member Function Documentation

◆ CheckParametersValidity()

void G4ParameterisationPolyconeRho::CheckParametersValidity ( )
overridevirtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 134 of file G4ParameterisationPolycone.cc.

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}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)

Referenced by G4ParameterisationPolyconeRho().

◆ ComputeDimensions()

void G4ParameterisationPolyconeRho::ComputeDimensions ( G4Polycone & pcone,
const G4int copyNo,
const G4VPhysicalVolume * physVol ) const
overridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 208 of file G4ParameterisationPolycone.cc.

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}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
void SetOriginalParameters(G4PolyconeHistorical *pars)
G4bool Reset()
void DumpInfo() const

◆ ComputeTransformation()

void G4ParameterisationPolyconeRho::ComputeTransformation ( const G4int copyNo,
G4VPhysicalVolume * physVol ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 173 of file G4ParameterisationPolycone.cc.

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}
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.0) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaxParameter()

G4double G4ParameterisationPolyconeRho::GetMaxParameter ( ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 163 of file G4ParameterisationPolycone.cc.

164{
165 auto msol = (G4Polycone*)(fmotherSolid);
166 G4PolyconeHistorical* original_pars = msol->GetOriginalParameters();
167 return original_pars->Rmax[0] - original_pars->Rmin[0];
168}

The documentation for this class was generated from the following files: