#include <G4ParameterisationTrd.hh>
|
| G4ParameterisationTrdY (EAxis axis, G4int nCopies, G4double width, G4double offset, G4VSolid *motherSolid, DivisionType divType) |
|
| ~G4ParameterisationTrdY () override |
|
G4double | GetMaxParameter () const override |
|
void | ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const override |
|
void | ComputeDimensions (G4Trd &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const override |
|
void | ComputeDimensions (G4Trap &, const G4int, const G4VPhysicalVolume *) const override |
|
| G4VParameterisationTrd (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType) |
|
| ~G4VParameterisationTrd () override |
|
| G4VDivisionParameterisation (EAxis axis, G4int nDiv, G4double width, G4double offset, DivisionType divType, G4VSolid *motherSolid=nullptr) |
|
| ~G4VDivisionParameterisation () override |
|
G4VSolid * | ComputeSolid (const G4int, G4VPhysicalVolume *) override |
|
const G4String & | GetType () const |
|
EAxis | GetAxis () const |
|
G4int | GetNoDiv () const |
|
G4double | GetWidth () const |
|
G4double | GetOffset () const |
|
G4VSolid * | GetMotherSolid () const |
|
void | SetType (const G4String &type) |
|
G4int | VolumeFirstCopyNo () const |
|
void | SetHalfGap (G4double hg) |
|
G4double | GetHalfGap () const |
|
| G4VPVParameterisation ()=default |
|
virtual | ~G4VPVParameterisation ()=default |
|
virtual G4Material * | ComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr) |
|
virtual G4bool | IsNested () const |
|
virtual G4VVolumeMaterialScanner * | GetMaterialScanner () |
|
Definition at line 122 of file G4ParameterisationTrd.hh.
◆ G4ParameterisationTrdY()
Definition at line 239 of file G4ParameterisationTrd.cc.
244{
247
250 {
252 width, offset );
253 }
255 {
257 nDiv, offset );
258 }
259
260#ifdef G4DIVDEBUG
262 {
263 G4cout <<
" G4ParameterisationTrdY no divisions " <<
fnDiv
264 <<
" = " << nDiv <<
G4endl
267 }
268#endif
269}
G4GLOB_DLL std::ostream G4cout
virtual void CheckParametersValidity()
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
static const G4int verbose
G4VParameterisationTrd(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
◆ ~G4ParameterisationTrdY()
G4ParameterisationTrdY::~G4ParameterisationTrdY |
( |
| ) |
|
|
overridedefault |
◆ ComputeDimensions() [1/2]
Reimplemented from G4VPVParameterisation.
Definition at line 345 of file G4ParameterisationTrd.cc.
347{
349 G4double pDx1 = msol->GetXHalfLength1();
350 G4double pDx2 = msol->GetXHalfLength2();
351 G4double pDz = msol->GetZHalfLength();
352
353 G4double pDy1 = msol->GetYHalfLength1();
354 G4double pDy2 = msol->GetYHalfLength2();
355
356 G4double yChangeRatio = (pDy2-pDy1) / (pDy2+pDy1);
363 G4double cyx1 = -pDy1+fOffset1 + (copyNo+0.5)*fWid1;
364 G4double cyx2 = -pDy2+fOffset2 + (copyNo+0.5)*fWid2;
365 G4double alp = std::atan( (cyx2-cyx1)/(pDz*2.) );
366
367 pDy1 =
fwidth/2. - fWidChange/2.;
368 pDy2 =
fwidth/2. + fWidChange/2.;
369
370
372 alp,
373 90*CLHEP::deg,
374 pDy1,
375 pDx1,
376 pDx1,
377 0.,
378 pDy2,
380 pDx2 -
fhgap * pDx2/pDx1,
381 0.);
382
383#ifdef G4DIVDEBUG
385 {
386 G4cout <<
" G4ParameterisationTrdY::ComputeDimensions():"
389 }
390#endif
391}
void SetAllParameters(G4double pDz, G4double pTheta, G4double pPhi, G4double pDy1, G4double pDx1, G4double pDx2, G4double pAlp1, G4double pDy2, G4double pDx3, G4double pDx4, G4double pAlp2)
◆ ComputeDimensions() [2/2]
Reimplemented from G4VPVParameterisation.
Definition at line 320 of file G4ParameterisationTrd.cc.
322{
323
324
326
327 G4double pDx1 = msol->GetXHalfLength1();
328 G4double pDx2 = msol->GetXHalfLength2();
329 G4double pDz = msol->GetZHalfLength();
331
333
334#ifdef G4DIVDEBUG
336 {
337 G4cout <<
" G4ParameterisationTrdY::ComputeDimensions():" <<
G4endl;
339 }
340#endif
341}
void SetAllParameters(G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
◆ ComputeTransformation()
Implements G4VDivisionParameterisation.
Definition at line 283 of file G4ParameterisationTrd.cc.
285{
287 G4double mdy = ( msol->GetYHalfLength1() + msol->GetYHalfLength2() ) / 2.;
288
289
292
294 {
295 origin.setY( posi );
296 }
297 else
298 {
299 std::ostringstream message;
300 message <<
"Only axes along Y are allowed ! Axis: " <<
faxis;
301 G4Exception(
"G4ParameterisationTrdY::ComputeTransformation()",
303 }
304
305#ifdef G4DIVDEBUG
307 {
308 G4cout << std::setprecision(8)
309 << " G4ParameterisationTrdY::ComputeTransformation " << copyNo
310 <<
" pos " << origin <<
" rot mat " <<
" axis " <<
faxis <<
G4endl;
311 }
312#endif
313
314
316}
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
void SetTranslation(const G4ThreeVector &v)
◆ GetMaxParameter()
G4double G4ParameterisationTrdY::GetMaxParameter |
( |
| ) |
const |
|
overridevirtual |
The documentation for this class was generated from the following files: