Geant4 11.2.2
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 () override
 
G4double GetMaxParameter () const override
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const override
 
void ComputeDimensions (G4Trd &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const override
 
void ComputeDimensions (G4Trap &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const override
 
- Public Member Functions inherited from G4VParameterisationTrd
 G4VParameterisationTrd (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
 ~G4VParameterisationTrd () override
 
- Public Member Functions inherited from G4VDivisionParameterisation
 G4VDivisionParameterisation (EAxis axis, G4int nDiv, G4double width, G4double offset, DivisionType divType, G4VSolid *motherSolid=nullptr)
 
 ~G4VDivisionParameterisation () override
 
G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *) override
 
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 G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr)
 
virtual G4bool IsNested () const
 
virtual G4VVolumeMaterialScannerGetMaterialScanner ()
 

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

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

◆ ~G4ParameterisationTrdX()

G4ParameterisationTrdX::~G4ParameterisationTrdX ( )
overridedefault

Member Function Documentation

◆ ComputeDimensions() [1/2]

void G4ParameterisationTrdX::ComputeDimensions ( G4Trap & trd,
const G4int copyNo,
const G4VPhysicalVolume * pv ) const
overridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 189 of file G4ParameterisationTrd.cc.

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

◆ ComputeDimensions() [2/2]

void G4ParameterisationTrdX::ComputeDimensions ( G4Trd & trd,
const G4int copyNo,
const G4VPhysicalVolume * pv ) const
overridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 166 of file G4ParameterisationTrd.cc.

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

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 128 of file G4ParameterisationTrd.cc.

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

◆ GetMaxParameter()

G4double G4ParameterisationTrdX::GetMaxParameter ( ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 120 of file G4ParameterisationTrd.cc.

121{
122 auto msol = (G4Trd*)(fmotherSolid);
123 return (msol->GetXHalfLength1()+msol->GetXHalfLength2());
124}

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