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

#include <G4ParameterisationTrd.hh>

+ Inheritance diagram for G4ParameterisationTrdX:

Public Member Functions

 G4ParameterisationTrdX (EAxis axis, G4int nCopies, G4double width, G4double offset, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationTrdX ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Trd &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const
 
void ComputeDimensions (G4Trap &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const
 
- Public Member Functions inherited from G4VParameterisationTrd
 G4VParameterisationTrd (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
virtual ~G4VParameterisationTrd ()
 
- 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 ()=default
 
virtual ~G4VPVParameterisation ()=default
 
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 G4VParameterisationTrd
G4bool bDivInTrap = false
 
- 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 75 of file G4ParameterisationTrd.hh.

Constructor & Destructor Documentation

◆ G4ParameterisationTrdX()

G4ParameterisationTrdX::G4ParameterisationTrdX ( EAxis  axis,
G4int  nCopies,
G4double  width,
G4double  offset,
G4VSolid motherSolid,
DivisionType  divType 
)

Definition at line 78 of file G4ParameterisationTrd.cc.

82 : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
83{
85 SetType( "DivisionTrdX" );
86
87 G4Trd* msol = (G4Trd*)(fmotherSolid);
88 if( divType == DivWIDTH )
89 {
91 width, offset );
92 }
93 else if( divType == DivNDIV )
94 {
96 nDiv, offset );
97 }
98
99#ifdef G4DIVDEBUG
100 if( verbose >= 1 )
101 {
102 G4cout << " G4ParameterisationTrdX - ## divisions " << fnDiv << " = "
103 << nDiv << G4endl
104 << " Offset " << foffset << " = " << offset << G4endl
105 << " Width " << fwidth << " = " << width << G4endl;
106 }
107#endif
108
109 G4double mpDx1 = msol->GetXHalfLength1();
110 G4double mpDx2 = msol->GetXHalfLength2();
111 if( std::fabs(mpDx1 - mpDx2) > kCarTolerance )
112 {
113 bDivInTrap = true;
114 }
115}
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Definition: G4Trd.hh:63
G4double GetXHalfLength2() const
G4double GetXHalfLength1() const
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const

◆ ~G4ParameterisationTrdX()

G4ParameterisationTrdX::~G4ParameterisationTrdX ( )

Definition at line 118 of file G4ParameterisationTrd.cc.

119{
120}

Member Function Documentation

◆ ComputeDimensions() [1/2]

void G4ParameterisationTrdX::ComputeDimensions ( G4Trap trd,
const G4int  copyNo,
const G4VPhysicalVolume pv 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 192 of file G4ParameterisationTrd.cc.

194{
195 G4Trd* msol = (G4Trd*)(fmotherSolid);
196 G4double pDy1 = msol->GetYHalfLength1();
197 G4double pDy2 = msol->GetYHalfLength2();
198 G4double pDz = msol->GetZHalfLength();
199 //fwidth is at Z=0
200 G4double pDx1 = msol->GetXHalfLength1();
201 G4double pDx2 = msol->GetXHalfLength2();
202 // G4double pDxAVE = (pDx1+pDx2)/2.;
203 G4double xChangeRatio = (pDx2-pDx1) / (pDx2+pDx1);
204 G4double fWidChange = xChangeRatio*fwidth;
205 G4double fWid1 = fwidth - fWidChange;
206 G4double fWid2 = fwidth + fWidChange;
207 G4double fOffsetChange = xChangeRatio*foffset/2.;
208 G4double fOffset1 = foffset - fOffsetChange;
209 G4double fOffset2 = foffset + fOffsetChange;
210 G4double cxy1 = -pDx1+fOffset1 + (copyNo+0.5)*fWid1;// centre of the side at y=-pDy1
211 G4double cxy2 = -pDx2+fOffset2 + (copyNo+0.5)*fWid2;// centre of the side at y=+pDy1
212 G4double alp = std::atan( (cxy2-cxy1)/(pDz*2.) );
213
214 pDx1 = fwidth/2. - fWidChange/2.;
215 pDx2 = fwidth/2. + fWidChange/2.;
216
217
218 trap.SetAllParameters ( pDz,
219 alp,
220 0.,
221 pDy1,
222 pDx1,
223 pDx1,
224 0.,
225 pDy2,
226 pDx2 - fhgap,
227 pDx2 - fhgap * pDx2/pDx1,
228 0.);
229
230#ifdef G4DIVDEBUG
231 if( verbose >= 2 )
232 {
233 G4cout << " G4ParameterisationTrdX::ComputeDimensions():"
234 << copyNo << G4endl;
235 trap.DumpInfo();
236 }
237#endif
238}
G4double GetYHalfLength2() const
G4double GetYHalfLength1() const
G4double GetZHalfLength() const

◆ ComputeDimensions() [2/2]

void G4ParameterisationTrdX::ComputeDimensions ( G4Trd trd,
const G4int  copyNo,
const G4VPhysicalVolume pv 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 169 of file G4ParameterisationTrd.cc.

171{
172 G4Trd* msol = (G4Trd*)(fmotherSolid);
173 G4double pDy1 = msol->GetYHalfLength1();
174 G4double pDy2 = msol->GetYHalfLength2();
175 G4double pDz = msol->GetZHalfLength();
176 G4double pDx = fwidth/2. - fhgap;
177
178 trd.SetAllParameters ( pDx, pDx, pDy1, pDy2, pDz );
179
180#ifdef G4DIVDEBUG
181 if( verbose >= 2 )
182 {
183 G4cout << " G4ParameterisationTrdX::ComputeDimensions():"
184 << copyNo << G4endl;
185 trd.DumpInfo();
186 }
187#endif
188}
void SetAllParameters(G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition: G4Trd.cc:128
void DumpInfo() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 131 of file G4ParameterisationTrd.cc.

134{
135 G4Trd* msol = (G4Trd*)(fmotherSolid );
136 G4double mdx = ( msol->GetXHalfLength1() + msol->GetXHalfLength2() ) / 2.;
137 //----- translation
138 G4ThreeVector origin(0.,0.,0.);
139 G4double posi;
140 posi = -mdx + foffset + (copyNo+0.5)*fwidth;
141 if( faxis == kXAxis )
142 {
143 origin.setX( posi );
144 }
145 else
146 {
147 std::ostringstream message;
148 message << "Only axes along X are allowed ! Axis: " << faxis;
149 G4Exception("G4ParameterisationTrdX::ComputeTransformation()",
150 "GeomDiv0002", FatalException, message);
151 }
152
153#ifdef G4DIVDEBUG
154 if( verbose >= 2 )
155 {
156 G4cout << std::setprecision(8)
157 << " G4ParameterisationTrdX::ComputeTransformation() "
158 << copyNo << G4endl
159 << " Position: " << origin << " - Axis: " << faxis << G4endl;
160 }
161#endif
162
163 //----- set translation
164 physVol->SetTranslation( origin );
165}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
void SetTranslation(const G4ThreeVector &v)
@ kXAxis
Definition: geomdefs.hh:55

◆ GetMaxParameter()

G4double G4ParameterisationTrdX::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 123 of file G4ParameterisationTrd.cc.

124{
125 G4Trd* msol = (G4Trd*)(fmotherSolid);
126 return (msol->GetXHalfLength1()+msol->GetXHalfLength2());
127}

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