Geant4-11
Public Member Functions | Protected Attributes | Private Attributes
G4tgbPlaceParamLinear Class Reference

#include <G4tgbPlaceParamLinear.hh>

Inheritance diagram for G4tgbPlaceParamLinear:
G4tgbPlaceParameterisation G4VPVParameterisation

Public Member Functions

void CheckNExtraData (G4tgrPlaceParameterisation *tgrParam, G4int nWcheck, WLSIZEtype st, const G4String &methodName)
 
virtual void ComputeDimensions (G4Box &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Cons &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Ellipsoid &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Hype &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Orb &, 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 (G4Sphere &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Torus &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Trap &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Trd &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Tubs &, const G4int, const G4VPhysicalVolume *) const
 
virtual G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr)
 
virtual G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *)
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
 G4tgbPlaceParamLinear (G4tgrPlaceParameterisation *)
 
EAxis GetAxis () const
 
virtual G4VVolumeMaterialScannerGetMaterialScanner ()
 
G4int GetNCopies () const
 
virtual G4bool IsNested () const
 
 ~G4tgbPlaceParamLinear ()
 

Protected Attributes

EAxis theAxis = kUndefined
 
G4int theNCopies = 0
 
G4RotationMatrixtheRotationMatrix = nullptr
 
G4ThreeVector theTranslation
 

Private Attributes

G4ThreeVector theDirection
 
G4double theOffset = 0.0
 
G4double theStep = 0.0
 

Detailed Description

Definition at line 45 of file G4tgbPlaceParamLinear.hh.

Constructor & Destructor Documentation

◆ G4tgbPlaceParamLinear()

G4tgbPlaceParamLinear::G4tgbPlaceParamLinear ( G4tgrPlaceParameterisation tgrParam)

Definition at line 43 of file G4tgbPlaceParamLinear.cc.

46{
47 //---- Get translation and rotation
48 if(tgrParam->GetParamType() == "LINEAR")
49 {
50 CheckNExtraData(tgrParam, 6, WLSIZE_EQ, "G4tgbPlaceParamLinear:");
52 G4ThreeVector(tgrParam->GetExtraData()[3], tgrParam->GetExtraData()[4],
53 tgrParam->GetExtraData()[5]);
55 }
56 else
57 {
58 CheckNExtraData(tgrParam, 3, WLSIZE_EQ, "G4tgbPlaceParamLinear:");
59 if(tgrParam->GetParamType() == "LINEAR_X")
60 {
61 theDirection = G4ThreeVector(1., 0., 0.);
63 }
64 else if(tgrParam->GetParamType() == "LINEAR_Y")
65 {
66 theDirection = G4ThreeVector(0., 1., 0.);
68 }
69 else if(tgrParam->GetParamType() == "LINEAR_Z")
70 {
71 theDirection = G4ThreeVector(0., 0., 1.);
73 }
74 }
75
76 if(theDirection.mag() == 0.)
77 {
78 G4Exception("G4tgbPlaceParamLinear::G4tgbPlaceParamLinear()",
79 "InvalidSetup", FatalException, "Direction is zero !");
80 }
81 else
82 {
84 }
85
86 theNCopies = G4int(tgrParam->GetExtraData()[0]);
87 theStep = tgrParam->GetExtraData()[1];
88 theOffset = tgrParam->GetExtraData()[2];
89
91
92#ifdef G4VERBOSE
94 {
95 G4cout << " G4tgbPlaceParamLinear::G4tgbPlaceParamLinear(): "
96 << " param type " << tgrParam->GetParamType() << G4endl
97 << " N copies " << theNCopies << G4endl << " step " << theStep
98 << G4endl << " offset " << theOffset << G4endl << " translation "
99 << theTranslation << G4endl << " direction " << theDirection
100 << G4endl << " axis " << theAxis << G4endl;
101 }
102#endif
103}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
CLHEP::Hep3Vector G4ThreeVector
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
@ WLSIZE_EQ
Definition: G4tgrUtils.hh:47
double mag() const
G4tgbPlaceParameterisation(G4tgrPlaceParameterisation *tgrParam)
void CheckNExtraData(G4tgrPlaceParameterisation *tgrParam, G4int nWcheck, WLSIZEtype st, const G4String &methodName)
static G4int GetVerboseLevel()
std::vector< G4double > GetExtraData() const
const G4String & GetParamType() const
@ kYAxis
Definition: geomdefs.hh:56
@ kXAxis
Definition: geomdefs.hh:55
@ kZAxis
Definition: geomdefs.hh:57

References G4tgbPlaceParameterisation::CheckNExtraData(), FatalException, G4cout, G4endl, G4Exception(), G4tgrPlaceParameterisation::GetExtraData(), G4tgrPlaceParameterisation::GetParamType(), G4tgrMessenger::GetVerboseLevel(), kXAxis, kYAxis, kZAxis, CLHEP::Hep3Vector::mag(), G4tgbPlaceParameterisation::theAxis, theDirection, G4tgbPlaceParameterisation::theNCopies, theOffset, theStep, G4tgbPlaceParameterisation::theTranslation, and WLSIZE_EQ.

◆ ~G4tgbPlaceParamLinear()

G4tgbPlaceParamLinear::~G4tgbPlaceParamLinear ( )

Definition at line 38 of file G4tgbPlaceParamLinear.cc.

39{
40}

Member Function Documentation

◆ CheckNExtraData()

void G4tgbPlaceParameterisation::CheckNExtraData ( G4tgrPlaceParameterisation tgrParam,
G4int  nWcheck,
WLSIZEtype  st,
const G4String methodName 
)
inherited

Definition at line 61 of file G4tgbPlaceParameterisation.cc.

64{
65 std::vector<G4double> extraData = tgrParam->GetExtraData();
66 G4int ndata = extraData.size();
67
68 G4String outStr = methodName + " " + tgrParam->GetType() + " ";
69 G4bool isOK = G4tgrUtils::CheckListSize(ndata, nWcheck, st, outStr);
70
71 if(!isOK)
72 {
73 G4String chartmp = G4UIcommand::ConvertToString(nWcheck);
74 outStr += chartmp + G4String(" words");
75 G4cerr << outStr;
76 G4cerr << " NUMBER OF WORDS " << ndata << G4endl;
77 G4Exception("G4tgbPlaceParameterisation::CheckNExtraData", "InvalidData",
78 FatalException, "Invalid data size.");
79 }
80}
bool G4bool
Definition: G4Types.hh:86
G4GLOB_DLL std::ostream G4cerr
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:445
const G4String & GetType() const
Definition: G4tgrPlace.hh:55
static G4bool CheckListSize(unsigned int nWreal, unsigned int nWcheck, WLSIZEtype st, G4String &outstr)
Definition: G4tgrUtils.cc:512

References G4tgrUtils::CheckListSize(), G4UIcommand::ConvertToString(), FatalException, G4cerr, G4endl, G4Exception(), G4tgrPlaceParameterisation::GetExtraData(), and G4tgrPlace::GetType().

Referenced by G4tgbPlaceParamCircle::G4tgbPlaceParamCircle(), G4tgbPlaceParamLinear(), and G4tgbPlaceParamSquare::G4tgbPlaceParamSquare().

◆ ComputeDimensions() [1/13]

virtual void G4VPVParameterisation::ComputeDimensions ( G4Box ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtualinherited

◆ ComputeDimensions() [2/13]

virtual void G4VPVParameterisation::ComputeDimensions ( G4Cons ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtualinherited

◆ ComputeDimensions() [3/13]

virtual void G4VPVParameterisation::ComputeDimensions ( G4Ellipsoid ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtualinherited

◆ ComputeDimensions() [4/13]

virtual void G4VPVParameterisation::ComputeDimensions ( G4Hype ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtualinherited

◆ ComputeDimensions() [5/13]

virtual void G4VPVParameterisation::ComputeDimensions ( G4Orb ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtualinherited

◆ ComputeDimensions() [6/13]

virtual void G4VPVParameterisation::ComputeDimensions ( G4Para ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtualinherited

◆ ComputeDimensions() [7/13]

virtual void G4VPVParameterisation::ComputeDimensions ( G4Polycone ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtualinherited

◆ ComputeDimensions() [8/13]

virtual void G4VPVParameterisation::ComputeDimensions ( G4Polyhedra ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtualinherited

◆ ComputeDimensions() [9/13]

virtual void G4VPVParameterisation::ComputeDimensions ( G4Sphere ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtualinherited

◆ ComputeDimensions() [10/13]

virtual void G4VPVParameterisation::ComputeDimensions ( G4Torus ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtualinherited

◆ ComputeDimensions() [11/13]

virtual void G4VPVParameterisation::ComputeDimensions ( G4Trap ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtualinherited

◆ ComputeDimensions() [12/13]

virtual void G4VPVParameterisation::ComputeDimensions ( G4Trd ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtualinherited

◆ ComputeDimensions() [13/13]

virtual void G4VPVParameterisation::ComputeDimensions ( G4Tubs ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtualinherited

◆ ComputeMaterial()

G4Material * G4VPVParameterisation::ComputeMaterial ( const G4int  repNo,
G4VPhysicalVolume currentVol,
const G4VTouchable parentTouch = nullptr 
)
virtualinherited

◆ ComputeSolid()

G4VSolid * G4VPVParameterisation::ComputeSolid ( const  G4int,
G4VPhysicalVolume pPhysicalVol 
)
virtualinherited

◆ ComputeTransformation()

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

Reimplemented from G4tgbPlaceParameterisation.

Definition at line 106 of file G4tgbPlaceParamLinear.cc.

108{
109 G4ThreeVector origin = theTranslation + copyNo * theStep * theDirection;
110
111#ifdef G4VERBOSE
113 {
114 G4cout << " G4tgbPlaceParamLinear::ComputeTransformation() -"
115 << physVol->GetName() << G4endl << " copyNo " << copyNo << " pos "
116 << origin << G4endl;
117 }
118#endif
119 //----- Set traslation and rotation
120 physVol->SetTranslation(origin);
121 physVol->SetCopyNo(copyNo);
123}
virtual void SetCopyNo(G4int CopyNo)=0
const G4String & GetName() const
void SetTranslation(const G4ThreeVector &v)
void SetRotation(G4RotationMatrix *)

References G4cout, G4endl, G4VPhysicalVolume::GetName(), G4tgrMessenger::GetVerboseLevel(), G4VPhysicalVolume::SetCopyNo(), G4VPhysicalVolume::SetRotation(), G4VPhysicalVolume::SetTranslation(), theDirection, G4tgbPlaceParameterisation::theRotationMatrix, theStep, and G4tgbPlaceParameterisation::theTranslation.

◆ GetAxis()

EAxis G4tgbPlaceParameterisation::GetAxis ( ) const
inlineinherited

Definition at line 60 of file G4tgbPlaceParameterisation.hh.

60{ return theAxis; }

References G4tgbPlaceParameterisation::theAxis.

Referenced by G4tgbVolume::ConstructG4PhysVol().

◆ GetMaterialScanner()

G4VVolumeMaterialScanner * G4VPVParameterisation::GetMaterialScanner ( )
virtualinherited

Reimplemented in G4VNestedParameterisation.

Definition at line 62 of file G4VPVParameterisation.cc.

63{
64 return nullptr;
65}

Referenced by G4Region::ScanVolumeTree().

◆ GetNCopies()

G4int G4tgbPlaceParameterisation::GetNCopies ( ) const
inlineinherited

Definition at line 59 of file G4tgbPlaceParameterisation.hh.

59{ return theNCopies; }

References G4tgbPlaceParameterisation::theNCopies.

Referenced by G4tgbVolume::ConstructG4PhysVol().

◆ IsNested()

G4bool G4VPVParameterisation::IsNested ( ) const
virtualinherited

Field Documentation

◆ theAxis

EAxis G4tgbPlaceParameterisation::theAxis = kUndefined
protectedinherited

◆ theDirection

G4ThreeVector G4tgbPlaceParamLinear::theDirection
private

Definition at line 56 of file G4tgbPlaceParamLinear.hh.

Referenced by ComputeTransformation(), and G4tgbPlaceParamLinear().

◆ theNCopies

G4int G4tgbPlaceParameterisation::theNCopies = 0
protectedinherited

◆ theOffset

G4double G4tgbPlaceParamLinear::theOffset = 0.0
private

Definition at line 58 of file G4tgbPlaceParamLinear.hh.

Referenced by G4tgbPlaceParamLinear().

◆ theRotationMatrix

G4RotationMatrix* G4tgbPlaceParameterisation::theRotationMatrix = nullptr
protectedinherited

◆ theStep

G4double G4tgbPlaceParamLinear::theStep = 0.0
private

Definition at line 57 of file G4tgbPlaceParamLinear.hh.

Referenced by ComputeTransformation(), and G4tgbPlaceParamLinear().

◆ theTranslation

G4ThreeVector G4tgbPlaceParameterisation::theTranslation
protectedinherited

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