Geant4 9.6.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=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 133 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 289 of file G4ParameterisationTrd.cc.

293 : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
294{
296 SetType( "DivisionTrdY" );
297
298 G4Trd* msol = (G4Trd*)(fmotherSolid);
299 if( divType == DivWIDTH )
300 {
302 width, offset );
303 }
304 else if( divType == DivNDIV )
305 {
307 nDiv, offset );
308 }
309
310#ifdef G4DIVDEBUG
311 if( verbose >= 1 )
312 {
313 G4cout << " G4ParameterisationTrdY no divisions " << fnDiv
314 << " = " << nDiv << G4endl
315 << " Offset " << foffset << " = " << offset << G4endl
316 << " width " << fwidth << " = " << width << G4endl;
317 }
318#endif
319}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT 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 322 of file G4ParameterisationTrd.cc.

323{
324}

Member Function Documentation

◆ CheckParametersValidity()

void G4ParameterisationTrdY::CheckParametersValidity ( )
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 395 of file G4ParameterisationTrd.cc.

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

Referenced by G4ParameterisationTrdY().

◆ ComputeDimensions()

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

Reimplemented from G4VPVParameterisation.

Definition at line 371 of file G4ParameterisationTrd.cc.

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

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

◆ GetMaxParameter()

G4double G4ParameterisationTrdY::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 327 of file G4ParameterisationTrd.cc.

328{
329 G4Trd* msol = (G4Trd*)(fmotherSolid);
330 return 2*msol->GetYHalfLength1();
331}

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