Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4ParameterisationTrdX Class Reference

#include <G4ParameterisationTrd.hh>

Inheritance diagram for G4ParameterisationTrdX:
G4VParameterisationTrd G4VDivisionParameterisation G4VPVParameterisation

Public Member Functions

 G4ParameterisationTrdX (EAxis axis, G4int nCopies, G4double width, G4double offset, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationTrdX ()
 
void CheckParametersValidity ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Trd &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const
 
void ComputeDimensions (G4Trap &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const
 
G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *)
 
- 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 ()
 
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 G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=0)
 
virtual G4bool IsNested () const
 
virtual G4VVolumeMaterialScannerGetMaterialScanner ()
 

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
 
void CheckOffset (G4double maxPar)
 
void CheckNDivAndWidth (G4double maxPar)
 
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 G4ThreadLocal G4int verbose = 5
 

Detailed Description

Definition at line 82 of file G4ParameterisationTrd.hh.

Constructor & Destructor Documentation

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

Definition at line 83 of file G4ParameterisationTrd.cc.

References G4VParameterisationTrd::bDivInTrap, G4VDivisionParameterisation::CalculateNDiv(), G4VDivisionParameterisation::CalculateWidth(), CheckParametersValidity(), DivNDIV, DivWIDTH, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fnDiv, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Trd::GetXHalfLength1(), G4Trd::GetXHalfLength2(), G4VDivisionParameterisation::kCarTolerance, G4VDivisionParameterisation::SetType(), G4VDivisionParameterisation::verbose, and width.

86  : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
87 {
89  SetType( "DivisionTrdX" );
90 
91  G4Trd* msol = (G4Trd*)(fmotherSolid);
92  if( divType == DivWIDTH )
93  {
95  width, offset );
96  }
97  else if( divType == DivNDIV )
98  {
100  nDiv, offset );
101  }
102 
103 #ifdef G4DIVDEBUG
104  if( verbose >= 1 )
105  {
106  G4cout << " G4ParameterisationTrdX - ## divisions " << fnDiv << " = "
107  << nDiv << G4endl
108  << " Offset " << foffset << " = " << offset << G4endl
109  << " Width " << fwidth << " = " << width << G4endl;
110  }
111 #endif
112 
113  G4double mpDx1 = msol->GetXHalfLength1();
114  G4double mpDx2 = msol->GetXHalfLength2();
115  if( std::fabs(mpDx1 - mpDx2) > kCarTolerance )
116  {
117  bDivInTrap = true;
118  }
119 }
void SetType(const G4String &type)
#define width
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
Definition: G4Trd.hh:71
G4double GetXHalfLength2() const
G4GLOB_DLL std::ostream G4cout
G4VParameterisationTrd(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
#define G4endl
Definition: G4ios.hh:61
G4double GetXHalfLength1() const
double G4double
Definition: G4Types.hh:76
G4ParameterisationTrdX::~G4ParameterisationTrdX ( )

Definition at line 122 of file G4ParameterisationTrd.cc.

123 {
124 }

Member Function Documentation

void G4ParameterisationTrdX::CheckParametersValidity ( )
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 259 of file G4ParameterisationTrd.cc.

References G4VDivisionParameterisation::CheckParametersValidity().

Referenced by G4ParameterisationTrdX().

260 {
262 /*
263  G4Trd* msol = (G4Trd*)(fmotherSolid);
264 
265  G4double mpDx1 = msol->GetXHalfLength1();
266  G4double mpDx2 = msol->GetXHalfLength2();
267  bDivInTrap = false;
268 
269  if( std::fabs(mpDx1 - mpDx2) > kCarTolerance )
270  {
271  std::ostringstream message;
272  message << "Invalid solid specification. NOT supported." << G4endl
273  << "Making a division of a TRD along axis X," << G4endl
274  << "while the X half lengths are not equal," << G4endl
275  << "is not (yet) supported. It will result" << G4endl
276  << "in non-equal division solids.";
277  G4Exception("G4ParameterisationTrdX::CheckParametersValidity()",
278  "GeomDiv0001", FatalException, message);
279  }
280 */
281 }
void G4ParameterisationTrdX::ComputeDimensions ( G4Trd trd,
const G4int  copyNo,
const G4VPhysicalVolume pv 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 183 of file G4ParameterisationTrd.cc.

References G4VSolid::DumpInfo(), G4VDivisionParameterisation::fhgap, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Trd::GetYHalfLength1(), G4Trd::GetYHalfLength2(), G4Trd::GetZHalfLength(), G4Trd::SetAllParameters(), and G4VDivisionParameterisation::verbose.

184 {
185  G4Trd* msol = (G4Trd*)(fmotherSolid);
186 
187  G4double pDy1 = msol->GetYHalfLength1();
188  G4double pDy2 = msol->GetYHalfLength2();
189  G4double pDz = msol->GetZHalfLength();
190  G4double pDx = fwidth/2. - fhgap;
191 
192  trd.SetAllParameters ( pDx, pDx, pDy1, pDy2, pDz );
193 
194 #ifdef G4DIVDEBUG
195  if( verbose >= 2 )
196  {
197  G4cout << " G4ParameterisationTrdX::ComputeDimensions():"
198  << copyNo << G4endl;
199  trd.DumpInfo();
200  }
201 #endif
202 }
void SetAllParameters(G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition: G4Trd.cc:169
G4double GetYHalfLength1() const
Definition: G4Trd.hh:71
G4double GetZHalfLength() const
void DumpInfo() const
G4GLOB_DLL std::ostream G4cout
G4double GetYHalfLength2() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void G4ParameterisationTrdX::ComputeDimensions ( G4Trap trd,
const G4int  copyNo,
const G4VPhysicalVolume pv 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 220 of file G4ParameterisationTrd.cc.

References G4VSolid::DumpInfo(), G4VDivisionParameterisation::fhgap, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fnDiv, G4VDivisionParameterisation::foffset, G4cout, G4endl, G4Trd::GetXHalfLength1(), G4Trd::GetXHalfLength2(), G4Trd::GetYHalfLength1(), G4Trd::GetYHalfLength2(), G4Trd::GetZHalfLength(), G4Trap::SetAllParameters(), and G4VDivisionParameterisation::verbose.

222 {
223  G4Trd* msol = (G4Trd*)(fmotherSolid);
224  G4double pDy1 = msol->GetYHalfLength1();
225  G4double pDy2 = msol->GetYHalfLength2();
226  G4double pDz = msol->GetZHalfLength();
227  G4double pDx1 = msol->GetXHalfLength1()/fnDiv;
228  G4double pDx2 = msol->GetXHalfLength2()/fnDiv;
229 
230  G4double cxy1 = -msol->GetXHalfLength1() + foffset
231  + (copyNo+0.5)*pDx1*2;// centre of the side at y=-pDy1
232  G4double cxy2 = -msol->GetXHalfLength2() + foffset
233  + (copyNo+0.5)*pDx2*2;// centre of the side at y=+pDy1
234  G4double alp = std::atan( (cxy2-cxy1)/pDz );
235 
236  trap.SetAllParameters ( pDz,
237  0.,
238  0.,
239  pDy1,
240  pDx1,
241  pDx2,
242  alp,
243  pDy2,
244  pDx1 - fhgap,
245  pDx2 - fhgap * pDx2/pDx1,
246  alp);
247 
248 #ifdef G4DIVDEBUG
249  if( verbose >= 2 )
250  {
251  G4cout << " G4ParameterisationTrdX::ComputeDimensions():"
252  << copyNo << G4endl;
253  trap.DumpInfo();
254  }
255 #endif
256 }
G4double GetYHalfLength1() const
Definition: G4Trd.hh:71
G4double GetZHalfLength() const
G4double GetXHalfLength2() const
G4GLOB_DLL std::ostream G4cout
G4double GetYHalfLength2() const
#define G4endl
Definition: G4ios.hh:61
G4double GetXHalfLength1() const
double G4double
Definition: G4Types.hh:76
G4VSolid * G4ParameterisationTrdX::ComputeSolid ( const G4int  i,
G4VPhysicalVolume pv 
)
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 205 of file G4ParameterisationTrd.cc.

References G4VParameterisationTrd::bDivInTrap, G4VDivisionParameterisation::ComputeSolid(), and G4VDivisionParameterisation::fmotherSolid.

206 {
207  if( bDivInTrap )
208  {
210  }
211  else
212  {
213  return fmotherSolid;
214  }
215 }
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
void G4ParameterisationTrdX::ComputeTransformation ( const G4int  copyNo,
G4VPhysicalVolume physVol 
) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 136 of file G4ParameterisationTrd.cc.

References G4VParameterisationTrd::bDivInTrap, FatalException, G4VDivisionParameterisation::faxis, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fnDiv, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Exception(), G4Trd::GetXHalfLength1(), G4Trd::GetXHalfLength2(), kXAxis, G4VPhysicalVolume::SetTranslation(), CLHEP::Hep3Vector::setX(), and G4VDivisionParameterisation::verbose.

138 {
139  G4Trd* msol = (G4Trd*)(fmotherSolid );
140  G4double mdx = ( msol->GetXHalfLength1() + msol->GetXHalfLength2() ) / 2.;
141 
142  //----- translation
143  G4ThreeVector origin(0.,0.,0.);
144  G4double posi;
145  if( !bDivInTrap )
146  {
147  posi = -mdx + foffset + (copyNo+0.5)*fwidth;
148  }
149  else
150  {
151  G4double aveHL = (msol->GetXHalfLength1()+msol->GetXHalfLength2())/2.;
152  posi = - aveHL + foffset + (copyNo+0.5)*aveHL/fnDiv*2;
153  }
154  if( faxis == kXAxis )
155  {
156  origin.setX( posi );
157  }
158  else
159  {
160  std::ostringstream message;
161  message << "Only axes along X are allowed ! Axis: " << faxis;
162  G4Exception("G4ParameterisationTrdX::ComputeTransformation()",
163  "GeomDiv0002", FatalException, message);
164  }
165 
166 #ifdef G4DIVDEBUG
167  if( verbose >= 2 )
168  {
169  G4cout << std::setprecision(8)
170  << " G4ParameterisationTrdX::ComputeTransformation() "
171  << copyNo << G4endl
172  << " Position: " << origin << " - Axis: " << faxis << G4endl;
173  }
174 #endif
175 
176  //----- set translation
177  physVol->SetTranslation( origin );
178 }
Definition: G4Trd.hh:71
G4double GetXHalfLength2() const
G4GLOB_DLL std::ostream G4cout
void SetTranslation(const G4ThreeVector &v)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
G4double GetXHalfLength1() const
double G4double
Definition: G4Types.hh:76
G4double G4ParameterisationTrdX::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 127 of file G4ParameterisationTrd.cc.

References G4VDivisionParameterisation::fmotherSolid, G4Trd::GetXHalfLength1(), and G4Trd::GetXHalfLength2().

128 {
129  G4Trd* msol = (G4Trd*)(fmotherSolid);
130  return (msol->GetXHalfLength1()+msol->GetXHalfLength2());
131 }
Definition: G4Trd.hh:71
G4double GetXHalfLength2() const
G4double GetXHalfLength1() const

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