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

#include <G4ParameterisationTrd.hh>

+ Inheritance diagram for G4ParameterisationTrdY:

Public Member Functions

 G4ParameterisationTrdY (EAxis axis, G4int nCopies, G4double width, G4double offset, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationTrdY () 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 &, const G4int, const G4VPhysicalVolume *) 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 122 of file G4ParameterisationTrd.hh.

Constructor & Destructor Documentation

◆ G4ParameterisationTrdY()

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

Definition at line 239 of file G4ParameterisationTrd.cc.

243 : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
244{
246 SetType( "DivisionTrdY" );
247
248 auto msol = (G4Trd*)(fmotherSolid);
249 if( divType == DivWIDTH )
250 {
251 fnDiv = CalculateNDiv( 2*msol->GetYHalfLength1(),
252 width, offset );
253 }
254 else if( divType == DivNDIV )
255 {
256 fwidth = CalculateWidth( 2*msol->GetYHalfLength1(),
257 nDiv, offset );
258 }
259
260#ifdef G4DIVDEBUG
261 if( verbose >= 1 )
262 {
263 G4cout << " G4ParameterisationTrdY no divisions " << fnDiv
264 << " = " << nDiv << G4endl
265 << " Offset " << foffset << " = " << offset << G4endl
266 << " width " << fwidth << " = " << width << G4endl;
267 }
268#endif
269}
#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)

◆ ~G4ParameterisationTrdY()

G4ParameterisationTrdY::~G4ParameterisationTrdY ( )
overridedefault

Member Function Documentation

◆ ComputeDimensions() [1/2]

void G4ParameterisationTrdY::ComputeDimensions ( G4Trap & trap,
const G4int copyNo,
const G4VPhysicalVolume *  ) const
overridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 345 of file G4ParameterisationTrd.cc.

347{
348 auto msol = (G4Trd*)(fmotherSolid);
349 G4double pDx1 = msol->GetXHalfLength1();
350 G4double pDx2 = msol->GetXHalfLength2();
351 G4double pDz = msol->GetZHalfLength();
352 //fwidth is at Z=0
353 G4double pDy1 = msol->GetYHalfLength1();
354 G4double pDy2 = msol->GetYHalfLength2();
355 // G4double pDxAVE = (pDy1+pDy2)/2.;
356 G4double yChangeRatio = (pDy2-pDy1) / (pDy2+pDy1);
357 G4double fWidChange = yChangeRatio*fwidth;
358 G4double fWid1 = fwidth - fWidChange;
359 G4double fWid2 = fwidth + fWidChange;
360 G4double fOffsetChange = yChangeRatio*foffset/2.;
361 G4double fOffset1 = foffset - fOffsetChange;
362 G4double fOffset2 = foffset + fOffsetChange;
363 G4double cyx1 = -pDy1+fOffset1 + (copyNo+0.5)*fWid1;// centre of the side at y=-pDy1
364 G4double cyx2 = -pDy2+fOffset2 + (copyNo+0.5)*fWid2;// centre of the side at y=+pDy1
365 G4double alp = std::atan( (cyx2-cyx1)/(pDz*2.) );
366
367 pDy1 = fwidth/2. - fWidChange/2.;
368 pDy2 = fwidth/2. + fWidChange/2.;
369
370
371 trap.SetAllParameters ( pDz,
372 alp,
373 90*CLHEP::deg,
374 pDy1,
375 pDx1,
376 pDx1,
377 0.,
378 pDy2,
379 pDx2 - fhgap,
380 pDx2 - fhgap * pDx2/pDx1,
381 0.);
382
383#ifdef G4DIVDEBUG
384 if( verbose >= 2 )
385 {
386 G4cout << " G4ParameterisationTrdY::ComputeDimensions():"
387 << copyNo << G4endl;
388 trap.DumpInfo();
389 }
390#endif
391}
double G4double
Definition G4Types.hh:83
void SetAllParameters(G4double pDz, G4double pTheta, G4double pPhi, G4double pDy1, G4double pDx1, G4double pDx2, G4double pAlp1, G4double pDy2, G4double pDx3, G4double pDx4, G4double pAlp2)
Definition G4Trap.cc:280
void DumpInfo() const

◆ ComputeDimensions() [2/2]

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

Reimplemented from G4VPVParameterisation.

Definition at line 320 of file G4ParameterisationTrd.cc.

322{
323 //---- The division along Y of a Trd will result a Trd, only
324 //--- if Y at -Z and +Z are equal, else use the G4Trap version
325 auto msol = (G4Trd*)(fmotherSolid);
326
327 G4double pDx1 = msol->GetXHalfLength1();
328 G4double pDx2 = msol->GetXHalfLength2();
329 G4double pDz = msol->GetZHalfLength();
330 G4double pDy = fwidth/2. - fhgap;
331
332 trd.SetAllParameters ( pDx1, pDx2, pDy, pDy, pDz );
333
334#ifdef G4DIVDEBUG
335 if( verbose >= 2 )
336 {
337 G4cout << " G4ParameterisationTrdY::ComputeDimensions():" << G4endl;
338 trd.DumpInfo();
339 }
340#endif
341}
void SetAllParameters(G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition G4Trd.cc:126

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 283 of file G4ParameterisationTrd.cc.

285{
286 auto msol = (G4Trd*)(fmotherSolid );
287 G4double mdy = ( msol->GetYHalfLength1() + msol->GetYHalfLength2() ) / 2.;
288
289 //----- translation
290 G4ThreeVector origin(0.,0.,0.);
291 G4double posi = -mdy + foffset + (copyNo+0.5)*fwidth;
292
293 if( faxis == kYAxis )
294 {
295 origin.setY( posi );
296 }
297 else
298 {
299 std::ostringstream message;
300 message << "Only axes along Y are allowed ! Axis: " << faxis;
301 G4Exception("G4ParameterisationTrdY::ComputeTransformation()",
302 "GeomDiv0002", FatalException, message);
303 }
304
305#ifdef G4DIVDEBUG
306 if( verbose >= 2 )
307 {
308 G4cout << std::setprecision(8)
309 << " G4ParameterisationTrdY::ComputeTransformation " << copyNo
310 << " pos " << origin << " rot mat " << " axis " << faxis << G4endl;
311 }
312#endif
313
314 //----- set translation
315 physVol->SetTranslation( origin );
316}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
void SetTranslation(const G4ThreeVector &v)
@ kYAxis
Definition geomdefs.hh:56

◆ GetMaxParameter()

G4double G4ParameterisationTrdY::GetMaxParameter ( ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 275 of file G4ParameterisationTrd.cc.

276{
277 auto msol = (G4Trd*)(fmotherSolid);
278 return (msol->GetYHalfLength1()+msol->GetYHalfLength2());
279}

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