Geant4 10.7.0
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 ()
 
void CheckParametersValidity ()
 
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
 
G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *)
 
- 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 ()
 
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 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

◆ CheckParametersValidity()

void G4ParameterisationTrdX::CheckParametersValidity ( )
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 255 of file G4ParameterisationTrd.cc.

256{
258/*
259 G4Trd* msol = (G4Trd*)(fmotherSolid);
260
261 G4double mpDx1 = msol->GetXHalfLength1();
262 G4double mpDx2 = msol->GetXHalfLength2();
263 bDivInTrap = false;
264
265 if( std::fabs(mpDx1 - mpDx2) > kCarTolerance )
266 {
267 std::ostringstream message;
268 message << "Invalid solid specification. NOT supported." << G4endl
269 << "Making a division of a TRD along axis X," << G4endl
270 << "while the X half lengths are not equal," << G4endl
271 << "is not (yet) supported. It will result" << G4endl
272 << "in non-equal division solids.";
273 G4Exception("G4ParameterisationTrdX::CheckParametersValidity()",
274 "GeomDiv0001", FatalException, message);
275 }
276*/
277}

Referenced by G4ParameterisationTrdX().

◆ ComputeDimensions() [1/2]

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

Reimplemented from G4VPVParameterisation.

Definition at line 216 of file G4ParameterisationTrd.cc.

218{
219 G4Trd* msol = (G4Trd*)(fmotherSolid);
220 G4double pDy1 = msol->GetYHalfLength1();
221 G4double pDy2 = msol->GetYHalfLength2();
222 G4double pDz = msol->GetZHalfLength();
223 G4double pDx1 = msol->GetXHalfLength1()/fnDiv;
224 G4double pDx2 = msol->GetXHalfLength2()/fnDiv;
225
226 G4double cxy1 = -msol->GetXHalfLength1() + foffset
227 + (copyNo+0.5)*pDx1*2;// centre of the side at y=-pDy1
228 G4double cxy2 = -msol->GetXHalfLength2() + foffset
229 + (copyNo+0.5)*pDx2*2;// centre of the side at y=+pDy1
230 G4double alp = std::atan( (cxy2-cxy1)/pDz );
231
232 trap.SetAllParameters ( pDz,
233 0.,
234 0.,
235 pDy1,
236 pDx1,
237 pDx2,
238 alp,
239 pDy2,
240 pDx1 - fhgap,
241 pDx2 - fhgap * pDx2/pDx1,
242 alp);
243
244#ifdef G4DIVDEBUG
245 if( verbose >= 2 )
246 {
247 G4cout << " G4ParameterisationTrdX::ComputeDimensions():"
248 << copyNo << G4endl;
249 trap.DumpInfo();
250 }
251#endif
252}
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 178 of file G4ParameterisationTrd.cc.

180{
181 G4Trd* msol = (G4Trd*)(fmotherSolid);
182
183 G4double pDy1 = msol->GetYHalfLength1();
184 G4double pDy2 = msol->GetYHalfLength2();
185 G4double pDz = msol->GetZHalfLength();
186 G4double pDx = fwidth/2. - fhgap;
187
188 trd.SetAllParameters ( pDx, pDx, pDy1, pDy2, pDz );
189
190#ifdef G4DIVDEBUG
191 if( verbose >= 2 )
192 {
193 G4cout << " G4ParameterisationTrdX::ComputeDimensions():"
194 << copyNo << G4endl;
195 trd.DumpInfo();
196 }
197#endif
198}
void SetAllParameters(G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition: G4Trd.cc:128
void DumpInfo() const

◆ ComputeSolid()

G4VSolid * G4ParameterisationTrdX::ComputeSolid ( const G4int  i,
G4VPhysicalVolume pv 
)
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 200 of file G4ParameterisationTrd.cc.

202{
203 if( bDivInTrap )
204 {
206 }
207 else
208 {
209 return fmotherSolid;
210 }
211}
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)

◆ 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
138 //----- translation
139 G4ThreeVector origin(0.,0.,0.);
140 G4double posi;
141 if( !bDivInTrap )
142 {
143 posi = -mdx + foffset + (copyNo+0.5)*fwidth;
144 }
145 else
146 {
147 G4double aveHL = (msol->GetXHalfLength1()+msol->GetXHalfLength2())/2.;
148 posi = - aveHL + foffset + (copyNo+0.5)*aveHL/fnDiv*2;
149 }
150 if( faxis == kXAxis )
151 {
152 origin.setX( posi );
153 }
154 else
155 {
156 std::ostringstream message;
157 message << "Only axes along X are allowed ! Axis: " << faxis;
158 G4Exception("G4ParameterisationTrdX::ComputeTransformation()",
159 "GeomDiv0002", FatalException, message);
160 }
161
162#ifdef G4DIVDEBUG
163 if( verbose >= 2 )
164 {
165 G4cout << std::setprecision(8)
166 << " G4ParameterisationTrdX::ComputeTransformation() "
167 << copyNo << G4endl
168 << " Position: " << origin << " - Axis: " << faxis << G4endl;
169 }
170#endif
171
172 //----- set translation
173 physVol->SetTranslation( origin );
174}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
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: