Geant4 10.7.0
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 ()
 
void CheckParametersValidity ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Trd &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 ()
 
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 129 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 286 of file G4ParameterisationTrd.cc.

290 : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
291{
293 SetType( "DivisionTrdY" );
294
295 G4Trd* msol = (G4Trd*)(fmotherSolid);
296 if( divType == DivWIDTH )
297 {
299 width, offset );
300 }
301 else if( divType == DivNDIV )
302 {
304 nDiv, offset );
305 }
306
307#ifdef G4DIVDEBUG
308 if( verbose >= 1 )
309 {
310 G4cout << " G4ParameterisationTrdY no divisions " << fnDiv
311 << " = " << nDiv << G4endl
312 << " Offset " << foffset << " = " << offset << G4endl
313 << " width " << fwidth << " = " << width << G4endl;
314 }
315#endif
316}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Definition: G4Trd.hh:63
G4double GetYHalfLength1() const
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const

◆ ~G4ParameterisationTrdY()

G4ParameterisationTrdY::~G4ParameterisationTrdY ( )

Definition at line 319 of file G4ParameterisationTrd.cc.

320{
321}

Member Function Documentation

◆ CheckParametersValidity()

void G4ParameterisationTrdY::CheckParametersValidity ( )
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 393 of file G4ParameterisationTrd.cc.

394{
396
397 G4Trd* msol = (G4Trd*)(fmotherSolid);
398
399 G4double mpDy1 = msol->GetYHalfLength1();
400 G4double mpDy2 = msol->GetYHalfLength2();
401
402 if( std::fabs(mpDy1 - mpDy2) > kCarTolerance )
403 {
404 std::ostringstream message;
405 message << "Invalid solid specification. NOT supported." << G4endl
406 << "Making a division of a TRD along axis Y while" << G4endl
407 << "the Y half lengths are not equal is not (yet)" << G4endl
408 << "supported. It will result in non-equal" << G4endl
409 << "division solids.";
410 G4Exception("G4ParameterisationTrdY::CheckParametersValidity()",
411 "GeomDiv0001", FatalException, message);
412 }
413}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
G4double GetYHalfLength2() const

Referenced by G4ParameterisationTrdY().

◆ ComputeDimensions()

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

Reimplemented from G4VPVParameterisation.

Definition at line 369 of file G4ParameterisationTrd.cc.

371{
372 //---- The division along Y of a Trd will result a Trd, only
373 //--- if Y at -Z and +Z are equal, else use the G4Trap version
374 G4Trd* msol = (G4Trd*)(fmotherSolid);
375
376 G4double pDx1 = msol->GetXHalfLength1();
377 G4double pDx2 = msol->GetXHalfLength2();
378 G4double pDz = msol->GetZHalfLength();
379 G4double pDy = fwidth/2. - fhgap;
380
381 trd.SetAllParameters ( pDx1, pDx2, pDy, pDy, pDz );
382
383#ifdef G4DIVDEBUG
384 if( verbose >= 2 )
385 {
386 G4cout << " G4ParameterisationTrdY::ComputeDimensions():" << G4endl;
387 trd.DumpInfo();
388 }
389#endif
390}
void SetAllParameters(G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition: G4Trd.cc:128
G4double GetXHalfLength2() const
G4double GetXHalfLength1() const
G4double GetZHalfLength() const
void DumpInfo() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 332 of file G4ParameterisationTrd.cc.

334{
335 G4Trd* msol = (G4Trd*)(fmotherSolid );
336 G4double mdy = msol->GetYHalfLength1();
337
338 //----- translation
339 G4ThreeVector origin(0.,0.,0.);
340 G4double posi = -mdy + foffset + (copyNo+0.5)*fwidth;
341
342 if( faxis == kYAxis )
343 {
344 origin.setY( posi );
345 }
346 else
347 {
348 std::ostringstream message;
349 message << "Only axes along Y are allowed ! Axis: " << faxis;
350 G4Exception("G4ParameterisationTrdY::ComputeTransformation()",
351 "GeomDiv0002", FatalException, message);
352 }
353
354#ifdef G4DIVDEBUG
355 if( verbose >= 2 )
356 {
357 G4cout << std::setprecision(8)
358 << " G4ParameterisationTrdY::ComputeTransformation " << copyNo
359 << " pos " << origin << " rot mat " << " axis " << faxis << G4endl;
360 }
361#endif
362
363 //----- set translation
364 physVol->SetTranslation( origin );
365}
void SetTranslation(const G4ThreeVector &v)
@ kYAxis
Definition: geomdefs.hh:56

◆ GetMaxParameter()

G4double G4ParameterisationTrdY::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 324 of file G4ParameterisationTrd.cc.

325{
326 G4Trd* msol = (G4Trd*)(fmotherSolid);
327 return 2*msol->GetYHalfLength1();
328}

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