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

#include <G4ParameterisationCons.hh>

+ Inheritance diagram for G4ParameterisationConsRho:

Public Member Functions

 G4ParameterisationConsRho (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationConsRho ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Cons &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const
 
- Public Member Functions inherited from G4VParameterisationCons
 G4VParameterisationCons (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
virtual ~G4VParameterisationCons ()
 
- 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 70 of file G4ParameterisationCons.hh.

Constructor & Destructor Documentation

◆ G4ParameterisationConsRho()

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

Definition at line 78 of file G4ParameterisationCons.cc.

82 : G4VParameterisationCons( axis, nDiv, width, offset, msolid, divType )
83{
85 SetType( "DivisionConsRho" );
86
87 G4Cons* msol = (G4Cons*)(fmotherSolid);
88 if( msol->GetInnerRadiusPlusZ() == 0. )
89 {
90 std::ostringstream message;
91 message << "OuterRadiusMinusZ = 0" << G4endl
92 << "Width is calculated as that of OuterRadiusMinusZ !";
93 G4Exception("G4ParameterisationConsRho::G4ParameterisationConsRho()",
94 "GeomDiv1001", JustWarning, message);
95 }
96
97 if( divType == DivWIDTH )
98 {
100 - msol->GetInnerRadiusMinusZ(), width, offset );
101 }
102 else if( divType == DivNDIV )
103 {
104 G4Cons* mconsol = (G4Cons*)(msolid);
106 - mconsol->GetInnerRadiusMinusZ(), nDiv, offset );
107 }
108
109#ifdef G4DIVDEBUG
110 if( verbose >= 1 )
111 {
112 G4cout << " G4ParameterisationConsRho - no divisions " << fnDiv << " = "
113 << nDiv << G4endl
114 << " Offset " << foffset << " = " << offset
115 << " - Width " << fwidth << " = " << width << G4endl;
116 }
117#endif
118}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Definition: G4Cons.hh:78
G4double GetInnerRadiusMinusZ() const
G4double GetInnerRadiusPlusZ() const
G4double GetOuterRadiusMinusZ() const
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const

◆ ~G4ParameterisationConsRho()

G4ParameterisationConsRho::~G4ParameterisationConsRho ( )

Definition at line 121 of file G4ParameterisationCons.cc.

122{
123}

Member Function Documentation

◆ ComputeDimensions()

void G4ParameterisationConsRho::ComputeDimensions ( G4Cons tubs,
const G4int  copyNo,
const G4VPhysicalVolume physVol 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 167 of file G4ParameterisationCons.cc.

170{
171 G4Cons* msol = (G4Cons*)(fmotherSolid);
172
173 G4double pRMin1 = msol->GetInnerRadiusMinusZ() + foffset + fwidth*copyNo;
174 G4double pRMax1 = msol->GetInnerRadiusMinusZ() + foffset + fwidth*(copyNo+1);
175
176 //width at Z Plus
177 //- G4double fwidthPlus =
178 // fwidth * ( msol->GetOuterRadiusPlusZ()/ msol->GetInnerRadiusPlusZ())
179 //- / ( msol->GetOuterRadiusMinusZ() - msol->GetInnerRadiusMinusZ());
180 G4double fwidthPlus = CalculateWidth( msol->GetOuterRadiusPlusZ()
181 - msol->GetInnerRadiusPlusZ(), fnDiv, foffset );
182 G4double pRMin2 = msol->GetInnerRadiusPlusZ()
183 + foffset + fwidthPlus * copyNo;
184 G4double pRMax2 = msol->GetInnerRadiusPlusZ()
185 + foffset + fwidthPlus * (copyNo+1);
186 G4double pDz = msol->GetZHalfLength();
187
188 G4double d_half_gap = fhgap * pRMax2 / pRMax1;
189 //- already rotated double pSR = foffset + copyNo*fwidth;
190 G4double pSPhi = msol->GetStartPhiAngle();
191 G4double pDPhi = msol->GetDeltaPhiAngle();;
192
193 cons.SetInnerRadiusMinusZ( pRMin1 + fhgap );
194 cons.SetOuterRadiusMinusZ( pRMax1 - fhgap );
195 cons.SetInnerRadiusPlusZ( pRMin2 + d_half_gap );
196 cons.SetOuterRadiusPlusZ( pRMax2 - d_half_gap );
197 cons.SetZHalfLength( pDz );
198 cons.SetStartPhiAngle( pSPhi, false );
199 cons.SetDeltaPhiAngle( pDPhi );
200
201#ifdef G4DIVDEBUG
202 if( verbose >= 2 )
203 {
204 G4cout << " G4ParameterisationConsRho::ComputeDimensions()" << G4endl
205 << " pRMin: " << pRMin1 << " - pRMax: " << pRMax1 << G4endl;
206 if( verbose >= 4 ) cons.DumpInfo();
207 }
208#endif
209}
double G4double
Definition: G4Types.hh:83
G4double GetOuterRadiusPlusZ() const
G4double GetStartPhiAngle() const
G4double GetDeltaPhiAngle() const
G4double GetZHalfLength() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 134 of file G4ParameterisationCons.cc.

136{
137 //----- translation
138 G4ThreeVector origin(0.,0.,0.);
139 //----- set translation
140 physVol->SetTranslation( origin );
141
142 //----- calculate rotation matrix: unit
143
144#ifdef G4DIVDEBUG
145 if( verbose >= 2 )
146 {
147 G4cout << " G4ParameterisationConsRho " << G4endl
148 << " Offset: " << foffset
149 << " - Width: " << fwidth << G4endl;
150 }
151#endif
152
153 ChangeRotMatrix( physVol );
154
155#ifdef G4DIVDEBUG
156 if( verbose >= 2 )
157 {
158 G4cout << std::setprecision(8) << " G4ParameterisationConsRho" << G4endl
159 << " Position: " << origin << " - Width: " << fwidth
160 << " - Axis: " << faxis << G4endl;
161 }
162#endif
163}
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.0) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaxParameter()

G4double G4ParameterisationConsRho::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 126 of file G4ParameterisationCons.cc.

127{
128 G4Cons* msol = (G4Cons*)(fmotherSolid);
129 return msol->GetOuterRadiusMinusZ() - msol->GetInnerRadiusMinusZ();
130}

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