Geant4 9.6.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=0)
 
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=0)
 
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 (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.) 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
 
- Protected Attributes inherited from G4VDivisionParameterisation
G4String ftype
 
EAxis faxis
 
G4int fnDiv
 
G4double fwidth
 
G4double foffset
 
DivisionType fDivisionType
 
G4VSolidfmotherSolid
 
G4bool fReflectedSolid
 
G4bool fDeleteSolid
 
G4int theVoluFirstCopyNo
 
G4double kCarTolerance
 
G4double fhgap
 
- Static Protected Attributes inherited from G4VDivisionParameterisation
static G4int verbose = 5
 

Detailed Description

Definition at line 81 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 82 of file G4ParameterisationTrd.cc.

86 : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
87{
89 SetType( "DivisionTrdX" );
90
91 G4Trd* msol = (G4Trd*)(fmotherSolid);
92 if( divType == DivWIDTH )
93 {
95 width, offset );
96 }
97 else if( divType == DivNDIV )
98 {
100 nDiv, offset );
101 }
102
103#ifdef G4DIVDEBUG
104 if( verbose >= 1 )
105 {
106 G4cout << " G4ParameterisationTrdX - ## divisions " << fnDiv << " = "
107 << nDiv << G4endl
108 << " Offset " << foffset << " = " << offset << G4endl
109 << " Width " << fwidth << " = " << width << G4endl;
110 }
111#endif
112
113 G4double mpDx1 = msol->GetXHalfLength1();
114 G4double mpDx2 = msol->GetXHalfLength2();
115 if( std::fabs(mpDx1 - mpDx2) > kCarTolerance )
116 {
117 bDivInTrap = true;
118 }
119}
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT 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 122 of file G4ParameterisationTrd.cc.

123{
124}

Member Function Documentation

◆ CheckParametersValidity()

void G4ParameterisationTrdX::CheckParametersValidity ( )
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 258 of file G4ParameterisationTrd.cc.

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

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 219 of file G4ParameterisationTrd.cc.

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

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

◆ ComputeSolid()

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

Reimplemented from G4VDivisionParameterisation.

Definition at line 203 of file G4ParameterisationTrd.cc.

205{
206 if( bDivInTrap )
207 {
209 }
210 else
211 {
212 return fmotherSolid;
213 }
214}
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 135 of file G4ParameterisationTrd.cc.

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

◆ GetMaxParameter()

G4double G4ParameterisationTrdX::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 127 of file G4ParameterisationTrd.cc.

128{
129 G4Trd* msol = (G4Trd*)(fmotherSolid);
130 return (msol->GetXHalfLength1()+msol->GetXHalfLength2());
131}

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