Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4tgbPlaceParamCircle Class Reference

#include <G4tgbPlaceParamCircle.hh>

Inheritance diagram for G4tgbPlaceParamCircle:
Collaboration diagram for G4tgbPlaceParamCircle:

Public Member Functions

 G4tgbPlaceParamCircle (G4tgrPlaceParameterisation *)
 
 ~G4tgbPlaceParamCircle ()
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
- Public Member Functions inherited from G4tgbPlaceParameterisation
 G4tgbPlaceParameterisation (G4tgrPlaceParameterisation *tgrParam)
 
virtual ~G4tgbPlaceParameterisation ()
 
void CheckNExtraData (G4tgrPlaceParameterisation *tgrParam, G4int nWcheck, WLSIZEtype st, const G4String &methodName)
 
G4int GetNCopies () const
 
EAxis GetAxis () const
 
- Public Member Functions inherited from G4VPVParameterisation
 G4VPVParameterisation ()
 
virtual ~G4VPVParameterisation ()
 
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 (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 Attributes inherited from G4tgbPlaceParameterisation
G4int theNCopies
 
EAxis theAxis
 
G4ThreeVector theTranslation
 
G4RotationMatrixtheRotationMatrix
 

Detailed Description

Definition at line 52 of file G4tgbPlaceParamCircle.hh.

Constructor & Destructor Documentation

G4tgbPlaceParamCircle::G4tgbPlaceParamCircle ( G4tgrPlaceParameterisation tgrParam)

Definition at line 52 of file G4tgbPlaceParamCircle.cc.

53  : G4tgbPlaceParameterisation(tgrParam)
54 {
55  //---- Get translation and rotation
56  if( tgrParam->GetParamType() == "CIRCLE" )
57  {
58  CheckNExtraData( tgrParam, 7, WLSIZE_EQ, "G4tgbPlaceParamCircle:");
59  theCircleAxis = G4ThreeVector( tgrParam->GetExtraData()[4],
60  tgrParam->GetExtraData()[5],
61  tgrParam->GetExtraData()[6] );
62 
63  G4ThreeVector zaxis(0.,0.,-1.);
64  if( zaxis.cross(theCircleAxis).mag() > 1.E-6 )
65  {
66  theDirInPlane = zaxis.cross(theCircleAxis);
67  }
68  else
69  {
70  theDirInPlane = theCircleAxis.cross(G4ThreeVector(0.,-1.,0.));
71  }
72  theAxis = kZAxis;
73  }
74  else
75  {
76  CheckNExtraData( tgrParam, 4, WLSIZE_EQ, "G4tgbPlaceParamCircle:");
77  if( tgrParam->GetParamType() == "CIRCLE_XY" ) {
78  theCircleAxis = G4ThreeVector(0.,0.,1.);
79  theDirInPlane = G4ThreeVector(1.,0.,0.);
80  theAxis = kZAxis;
81  } else if( tgrParam->GetParamType() == "CIRCLE_XZ" ) {
82  theCircleAxis = G4ThreeVector(0.,1.,0.);
83  theDirInPlane = G4ThreeVector(1.,0.,0.);
84  theAxis = kYAxis;
85  } else if( tgrParam->GetParamType() == "CIRCLE_YZ" ) {
86  theCircleAxis = G4ThreeVector(1.,0.,0.);
87  theDirInPlane = G4ThreeVector(0.,1.,0.);
88  theAxis = kXAxis;
89  }
90  }
91 
92  if( theCircleAxis.mag() == 0. )
93  {
94  G4Exception("G4tgbPlaceParamCircle::G4tgbPlaceParamCircle()",
95  "InvalidSetup", FatalException, "Circle axis is zero !");
96  }
97  theCircleAxis /= theCircleAxis.mag();
98 
99  theAxis = kZAxis;
100 
101  theNCopies = G4int(tgrParam->GetExtraData()[0]);
102  theStep = tgrParam->GetExtraData()[1];
103  theOffset = tgrParam->GetExtraData()[2];
104  theRadius = tgrParam->GetExtraData()[3];
105 
106 #ifdef G4VERBOSE
108  {
109  G4cout << " G4tgbPlaceParamCircle::G4tgbPlaceParamCircle():" << G4endl
110  << " param type " << tgrParam->GetParamType() << G4endl
111  << " no copies - " << theNCopies << G4endl
112  << " step - " << theStep << G4endl
113  << " offset - " << theOffset << G4endl
114  << " radius - " << theRadius << G4endl
115  << " circle axis - " << theCircleAxis << G4endl
116  << " dir in plane - " << theDirInPlane << G4endl;
117  }
118 #endif
119 }
void CheckNExtraData(G4tgrPlaceParameterisation *tgrParam, G4int nWcheck, WLSIZEtype st, const G4String &methodName)
CLHEP::Hep3Vector G4ThreeVector
int G4int
Definition: G4Types.hh:78
std::vector< G4double > GetExtraData() const
G4tgbPlaceParameterisation(G4tgrPlaceParameterisation *tgrParam)
G4GLOB_DLL std::ostream G4cout
static G4int GetVerboseLevel()
const G4String & GetParamType() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
Hep3Vector cross(const Hep3Vector &) const
double mag() const

Here is the call graph for this function:

G4tgbPlaceParamCircle::~G4tgbPlaceParamCircle ( )

Definition at line 45 of file G4tgbPlaceParamCircle.cc.

46 {
47 }

Member Function Documentation

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

Reimplemented from G4tgbPlaceParameterisation.

Definition at line 124 of file G4tgbPlaceParamCircle.cc.

125 {
126  G4double posi = theOffset + copyNo*theStep;
127  G4ThreeVector origin = theDirInPlane * theRadius;
128  origin.rotate( posi, theCircleAxis );
129 
130  //----- Calculate rotation matrix (so that all volumes point to the centre)
131  G4RotationMatrix rm;
132  rm.rotate( -posi, theCircleAxis );
133 
134  //----- Set translation and rotation
135  physVol->SetTranslation(origin);
136  G4RotationMatrix* pvRm = physVol->GetRotation();
137  if( pvRm == 0 )
138  {
139  pvRm = new G4RotationMatrix;
140  }
141  *pvRm = *theRotationMatrix * rm;
142  physVol->SetRotation(pvRm);
143  physVol->SetCopyNo( copyNo );
144 
145 #ifdef G4VERBOSE
147  {
148  G4cout << " G4tgbPlaceParamCircle::ComputeTransformation():"
149  << physVol->GetName() << G4endl
150  << " no copies - " << theNCopies << G4endl
151  << " centre - " << origin << G4endl
152  << " rotation-matrix - " << *pvRm << G4endl;
153  }
154 #endif
155 }
CLHEP::HepRotation G4RotationMatrix
const G4RotationMatrix * GetRotation() const
void SetRotation(G4RotationMatrix *)
HepRotation & rotate(double delta, const Hep3Vector &axis)
Definition: Rotation.cc:47
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
static G4int GetVerboseLevel()
virtual void SetCopyNo(G4int CopyNo)=0
void SetTranslation(const G4ThreeVector &v)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:28

Here is the call graph for this function:


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