Geant4 10.7.0
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 ()
 
void CheckParametersValidity ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Polycone &pcone, const G4int copyNo, const G4VPhysicalVolume *physVol) const
 
- Public Member Functions inherited from G4VParameterisationPolycone
 G4VParameterisationPolycone (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
virtual ~G4VParameterisationPolycone ()
 
- Public Member Functions inherited from G4VDivisionParameterisation
 G4VDivisionParameterisation (EAxis axis, G4int nDiv, G4double width, G4double offset, DivisionType divType, G4VSolid *motherSolid=nullptr)
 
virtual ~G4VDivisionParameterisation ()
 
virtual G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *)
 
virtual void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const =0
 
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 ()
 
virtual ~G4VPVParameterisation ()
 
virtual void ComputeTransformation (const G4int, G4VPhysicalVolume *) const =0
 
virtual G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *)
 
virtual G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr)
 
virtual G4bool IsNested () const
 
virtual G4VVolumeMaterialScannerGetMaterialScanner ()
 
virtual void ComputeDimensions (G4Box &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Tubs &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Trd &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Trap &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Cons &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Sphere &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Orb &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Ellipsoid &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Torus &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Para &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Polycone &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Polyhedra &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Hype &, const G4int, const G4VPhysicalVolume *) const
 

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
 
virtual void CheckParametersValidity ()
 
void CheckOffset (G4double maxPar)
 
void CheckNDivAndWidth (G4double maxPar)
 
virtual G4double GetMaxParameter () const =0
 
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 98 of file G4ParameterisationPolycone.cc.

102 : G4VParameterisationPolycone( axis, nDiv, width, offset, msolid, divType )
103{
105 SetType( "DivisionPolyconeRho" );
106
108 G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
109
110 if( divType == DivWIDTH )
111 {
112 fnDiv = CalculateNDiv( origparamMother->Rmax[0]
113 - origparamMother->Rmin[0], width, offset );
114 }
115 else if( divType == DivNDIV )
116 {
117 fwidth = CalculateWidth( origparamMother->Rmax[0]
118 - origparamMother->Rmin[0], nDiv, offset );
119 }
120
121#ifdef G4DIVDEBUG
122 if( verbose >= -1 )
123 {
124 G4cout << " G4ParameterisationPolyconeRho - # divisions " << fnDiv
125 << " = " << nDiv << G4endl
126 << " Offset " << foffset << " = " << offset << G4endl
127 << " Width " << fwidth << " = " << width << G4endl;
128 }
129#endif
130}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4PolyconeHistorical * GetOriginalParameters() const
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const

◆ ~G4ParameterisationPolyconeRho()

G4ParameterisationPolyconeRho::~G4ParameterisationPolyconeRho ( )

Definition at line 133 of file G4ParameterisationPolycone.cc.

134{
135}

Member Function Documentation

◆ CheckParametersValidity()

void G4ParameterisationPolyconeRho::CheckParametersValidity ( )
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 138 of file G4ParameterisationPolycone.cc.

139{
141
143
145 {
146 std::ostringstream message;
147 message << "In solid " << msol->GetName() << G4endl
148 << "Division along R will be done with a width "
149 << "different for each solid section." << G4endl
150 << "WIDTH will not be used !";
151 G4Exception("G4VParameterisationPolycone::CheckParametersValidity()",
152 "GeomDiv1001", JustWarning, message);
153 }
154 if( foffset != 0. )
155 {
156 std::ostringstream message;
157 message << "In solid " << msol->GetName() << G4endl
158 << "Division along R will be done with a width "
159 << "different for each solid section." << G4endl
160 << "OFFSET will not be used !";
161 G4Exception("G4VParameterisationPolycone::CheckParametersValidity()",
162 "GeomDiv1001", JustWarning, message);
163 }
164}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
G4String GetName() const

Referenced by G4ParameterisationPolyconeRho().

◆ ComputeDimensions()

void G4ParameterisationPolyconeRho::ComputeDimensions ( G4Polycone pcone,
const G4int  copyNo,
const G4VPhysicalVolume physVol 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 212 of file G4ParameterisationPolycone.cc.

215{
217
218 G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
219 G4PolyconeHistorical origparam( *origparamMother );
220 G4int nZplanes = origparamMother->Num_z_planes;
221
222 G4double width = 0.;
223 for( G4int ii = 0; ii < nZplanes; ++ii )
224 {
225 width = CalculateWidth( origparamMother->Rmax[ii]
226 - origparamMother->Rmin[ii], fnDiv, foffset );
227 origparam.Rmin[ii] = origparamMother->Rmin[ii]+foffset+width*copyNo;
228 origparam.Rmax[ii] = origparamMother->Rmin[ii]+foffset+width*(copyNo+1);
229 }
230
231 pcone.SetOriginalParameters(&origparam); // copy values & transfer pointers
232 pcone.Reset(); // reset to new solid parameters
233
234#ifdef G4DIVDEBUG
235 if( verbose >= -2 )
236 {
237 G4cout << "G4ParameterisationPolyconeRho::ComputeDimensions()" << G4endl
238 << "-- Parametrised pcone copy-number: " << copyNo << G4endl;
239 pcone.DumpInfo();
240 }
241#endif
242}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
void SetOriginalParameters(G4PolyconeHistorical *pars)
G4bool Reset()
Definition: G4Polycone.cc:443
void DumpInfo() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 177 of file G4ParameterisationPolycone.cc.

179{
180 //----- translation
181 G4ThreeVector origin(0.,0.,0.);
182 //----- set translation
183 physVol->SetTranslation( origin );
184
185 //----- calculate rotation matrix: unit
186
187#ifdef G4DIVDEBUG
188 if( verbose >= 2 )
189 {
190 G4cout << " G4ParameterisationPolyconeRho " << G4endl
191 << " foffset: " << foffset
192 << " - fwidth: " << fwidth << G4endl;
193 }
194#endif
195
196 ChangeRotMatrix( physVol );
197
198#ifdef G4DIVDEBUG
199 if( verbose >= 2 )
200 {
201 G4cout << std::setprecision(8) << " G4ParameterisationPolyconeRho "
202 << G4endl
203 << " Position: " << origin/mm
204 << " - Width: " << fwidth/deg
205 << " - Axis: " << faxis << G4endl;
206 }
207#endif
208}
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.0) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaxParameter()

G4double G4ParameterisationPolyconeRho::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 167 of file G4ParameterisationPolycone.cc.

168{
170 G4PolyconeHistorical* original_pars = msol->GetOriginalParameters();
171 return original_pars->Rmax[0] - original_pars->Rmin[0];
172}

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