Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParameterisationPolyconeZ Class Reference

#include <G4ParameterisationPolycone.hh>

+ Inheritance diagram for G4ParameterisationPolyconeZ:

Public Member Functions

 G4ParameterisationPolyconeZ (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationPolyconeZ ()
 
void CheckParametersValidity ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Polycone &pcone, const G4int copyNo, const G4VPhysicalVolume *physVol) const
 
- Public Member Functions inherited from G4VParameterisationPolycone
 G4VParameterisationPolycone (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
virtual ~G4VParameterisationPolycone ()
 
- Public Member Functions inherited from G4VDivisionParameterisation
 G4VDivisionParameterisation (EAxis axis, G4int nDiv, G4double width, G4double offset, DivisionType divType, G4VSolid *motherSolid=nullptr)
 
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=nullptr)
 
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 Member Functions inherited from G4VDivisionParameterisation
void ChangeRotMatrix (G4VPhysicalVolume *physVol, G4double rotZ=0.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 G4VDivisionParameterisation
G4String ftype
 
EAxis faxis
 
G4int fnDiv = 0
 
G4double fwidth = 0.0
 
G4double foffset = 0.0
 
DivisionType fDivisionType
 
G4VSolidfmotherSolid = nullptr
 
G4bool fReflectedSolid = false
 
G4bool fDeleteSolid = false
 
G4int theVoluFirstCopyNo = 1
 
G4double kCarTolerance
 
G4double fhgap = 0.0
 
- Static Protected Attributes inherited from G4VDivisionParameterisation
static G4ThreadLocal G4RotationMatrixfRot = nullptr
 
static const G4int verbose = 5
 

Detailed Description

Definition at line 173 of file G4ParameterisationPolycone.hh.

Constructor & Destructor Documentation

◆ G4ParameterisationPolyconeZ()

G4ParameterisationPolyconeZ::G4ParameterisationPolyconeZ ( EAxis  axis,
G4int  nCopies,
G4double  offset,
G4double  step,
G4VSolid motherSolid,
DivisionType  divType 
)

Definition at line 351 of file G4ParameterisationPolycone.cc.

355 : G4VParameterisationPolycone( axis, nDiv, width, offset, msolid, divType ),
356 fOrigParamMother(((G4Polycone*)fmotherSolid)->GetOriginalParameters())
357{
358
360 SetType( "DivisionPolyconeZ" );
361
362 if( divType == DivWIDTH )
363 {
364 fnDiv =
365 CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
366 - fOrigParamMother->Z_values[0] , width, offset );
367 }
368 else if( divType == DivNDIV )
369 {
370 fwidth =
371 CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
372 - fOrigParamMother->Z_values[0] , nDiv, offset );
373 }
374
375#ifdef G4DIVDEBUG
376 if( verbose >= 1 )
377 {
378 G4cout << " G4ParameterisationPolyconeZ - # divisions " << fnDiv << " = "
379 << nDiv << G4endl
380 << " Offset " << foffset << " = " << offset << G4endl
381 << " Width " << fwidth << " = " << width << G4endl;
382 }
383#endif
384}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void SetType(const G4String &type)
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const

◆ ~G4ParameterisationPolyconeZ()

G4ParameterisationPolyconeZ::~G4ParameterisationPolyconeZ ( )

Definition at line 387 of file G4ParameterisationPolycone.cc.

388{
389}

Member Function Documentation

◆ CheckParametersValidity()

void G4ParameterisationPolyconeZ::CheckParametersValidity ( )
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 436 of file G4ParameterisationPolycone.cc.

437{
439
440 // Division will be following the mother polycone segments
441 //
442 if( fDivisionType == DivNDIV )
443 {
444 if( fnDiv > fOrigParamMother->Num_z_planes-1 )
445 {
446 std::ostringstream error;
447 error << "Configuration not supported." << G4endl
448 << "Division along Z will be done by splitting in the defined"
449 << G4endl
450 << "Z planes, i.e, the number of division would be: "
451 << fOrigParamMother->Num_z_planes-1
452 << ", instead of: " << fnDiv << " !";
453 G4Exception("G4ParameterisationPolyconeZ::CheckParametersValidity()",
454 "GeomDiv0001", FatalException, error);
455 }
456 }
457
458 // Division will be done within one polycone segment
459 // with applying given width and offset
460 //
462 {
463 // Check if divided region does not span over more
464 // than one z segment
465
466 G4int isegstart = -1; // number of the segment containing start position
467 G4int isegend = -1; // number of the segment containing end position
468
469 if ( !fReflectedSolid )
470 {
471 // The start/end position of the divided region
472 //
473 G4double zstart
474 = fOrigParamMother->Z_values[0] + foffset;
475 G4double zend
476 = fOrigParamMother->Z_values[0] + foffset + fnDiv* fwidth;
477
478 G4int counter = 0;
479 while ( isegend < 0 && counter < fOrigParamMother->Num_z_planes-1 )
480 {
481 // first segment
482 if ( zstart >= fOrigParamMother->Z_values[counter] &&
483 zstart < fOrigParamMother->Z_values[counter+1] )
484 {
485 isegstart = counter;
486 }
487 // last segment
488 if ( zend > fOrigParamMother->Z_values[counter] &&
489 zend <= fOrigParamMother->Z_values[counter+1] )
490 {
491 isegend = counter;
492 }
493 ++counter;
494 } // Loop checking, 06.08.2015, G.Cosmo
495 }
496 else
497 {
498 // The start/end position of the divided region
499 //
500 G4double zstart
501 = fOrigParamMother->Z_values[0] - foffset;
502 G4double zend
503 = fOrigParamMother->Z_values[0] - ( foffset + fnDiv* fwidth);
504
505 G4int counter = 0;
506 while ( isegend < 0 && counter < fOrigParamMother->Num_z_planes-1 )
507 {
508 // first segment
509 if ( zstart <= fOrigParamMother->Z_values[counter] &&
510 zstart > fOrigParamMother->Z_values[counter+1] )
511 {
512 isegstart = counter;
513 }
514 // last segment
515 if ( zend < fOrigParamMother->Z_values[counter] &&
516 zend >= fOrigParamMother->Z_values[counter+1] )
517 {
518 isegend = counter;
519 }
520 ++counter;
521 } // Loop checking, 06.08.2015, G.Cosmo
522 }
523
524
525 if ( isegstart != isegend )
526 {
527 std::ostringstream message;
528 message << "Condiguration not supported." << G4endl
529 << "Division with user defined width." << G4endl
530 << "Solid " << fmotherSolid->GetName() << G4endl
531 << "Divided region is not between two z planes.";
532 G4Exception("G4ParameterisationPolyconeZ::CheckParametersValidity()",
533 "GeomDiv0001", FatalException, message);
534 }
535
536 fNSegment = isegstart;
537 }
538}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4String GetName() const

Referenced by G4ParameterisationPolyconeZ().

◆ ComputeDimensions()

void G4ParameterisationPolyconeZ::ComputeDimensions ( G4Polycone pcone,
const G4int  copyNo,
const G4VPhysicalVolume physVol 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 594 of file G4ParameterisationPolycone.cc.

597{
598
599 // Define division solid
600 //
601 G4PolyconeHistorical origparam;
602 G4int nz = 2;
603 origparam.Num_z_planes = nz;
604 origparam.Start_angle = fOrigParamMother->Start_angle;
605 origparam.Opening_angle = fOrigParamMother->Opening_angle;
606
607 // Define division solid z sections
608 //
609 origparam.Z_values = new G4double[nz];
610 origparam.Rmin = new G4double[nz];
611 origparam.Rmax = new G4double[nz];
612
613 if ( fDivisionType == DivNDIV )
614 {
615 // The position of the centre of copyNo-th mother polycone segment
616 G4double posi = (fOrigParamMother->Z_values[copyNo]
617 + fOrigParamMother->Z_values[copyNo+1])/2;
618
619 origparam.Z_values[0] = fOrigParamMother->Z_values[copyNo] - posi;
620 origparam.Z_values[1] = fOrigParamMother->Z_values[copyNo+1] - posi;
621 origparam.Rmin[0] = fOrigParamMother->Rmin[copyNo];
622 origparam.Rmin[1] = fOrigParamMother->Rmin[copyNo+1];
623 origparam.Rmax[0] = fOrigParamMother->Rmax[copyNo];
624 origparam.Rmax[1] = fOrigParamMother->Rmax[copyNo+1];
625 }
626
628 {
629 if ( !fReflectedSolid )
630 {
631 origparam.Z_values[0] = - fwidth/2.;
632 origparam.Z_values[1] = fwidth/2.;
633
634 // The position of the centre of copyNo-th division
635 //
636 G4double posi = fOrigParamMother->Z_values[0]
637 + foffset + (2*copyNo + 1) * fwidth/2.;
638
639 // The first and last z sides z values
640 //
641 G4double zstart = posi - fwidth/2.;
642 G4double zend = posi + fwidth/2.;
643 origparam.Rmin[0] = GetRmin(zstart, fNSegment);
644 origparam.Rmax[0] = GetRmax(zstart, fNSegment);
645 origparam.Rmin[1] = GetRmin(zend, fNSegment);
646 origparam.Rmax[1] = GetRmax(zend, fNSegment);
647 }
648 else
649 {
650 origparam.Z_values[0] = fwidth/2.;
651 origparam.Z_values[1] = - fwidth/2.;
652
653 // The position of the centre of copyNo-th division
654 //
655 G4double posi = fOrigParamMother->Z_values[0]
656 - ( foffset + (2*copyNo + 1) * fwidth/2.);
657
658 // The first and last z sides z values
659 //
660 G4double zstart = posi + fwidth/2.;
661 G4double zend = posi - fwidth/2.;
662 origparam.Rmin[0] = GetRmin(zstart, fNSegment);
663 origparam.Rmax[0] = GetRmax(zstart, fNSegment);
664 origparam.Rmin[1] = GetRmin(zend, fNSegment);
665 origparam.Rmax[1] = GetRmax(zend, fNSegment);
666 }
667
668 // It can happen due to rounding errors
669 //
670 if ( origparam.Rmin[0] < 0.0 ) origparam.Rmin[0] = 0.0;
671 if ( origparam.Rmin[nz-1] < 0.0 ) origparam.Rmin[1] = 0.0;
672 }
673
674 pcone.SetOriginalParameters(&origparam); // copy values & transfer pointers
675 pcone.Reset(); // reset to new solid parameters
676
677#ifdef G4DIVDEBUG
678 if( verbose >= 2 )
679 {
680 G4cout << "G4ParameterisationPolyconeZ::ComputeDimensions()" << G4endl
681 << "-- Parametrised pcone copy-number: " << copyNo << G4endl;
682 pcone.DumpInfo();
683 }
684#endif
685}
void SetOriginalParameters(G4PolyconeHistorical *pars)
G4bool Reset()
Definition: G4Polycone.cc:443
void DumpInfo() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 542 of file G4ParameterisationPolycone.cc.

544{
545 if ( fDivisionType == DivNDIV )
546 {
547 // The position of the centre of copyNo-th mother polycone segment
548 //
549 G4double posi = ( fOrigParamMother->Z_values[copyNo]
550 + fOrigParamMother->Z_values[copyNo+1])/2;
551 physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
552 }
553
555 {
556 // The position of the centre of copyNo-th division
557 //
558 G4double posi = fOrigParamMother->Z_values[0];
559
560 if ( !fReflectedSolid )
561 posi += foffset + (2*copyNo + 1) * fwidth/2.;
562 else
563 posi -= foffset + (2*copyNo + 1) * fwidth/2.;
564
565 physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
566 }
567
568 //----- calculate rotation matrix: unit
569
570#ifdef G4DIVDEBUG
571 if( verbose >= 2 )
572 {
573 G4cout << " G4ParameterisationPolyconeZ - position: " << posi << G4endl
574 << " copyNo: " << copyNo << " - foffset: " << foffset/deg
575 << " - fwidth: " << fwidth/deg << G4endl;
576 }
577#endif
578
579 ChangeRotMatrix( physVol );
580
581#ifdef G4DIVDEBUG
582 if( verbose >= 2 )
583 {
584 G4cout << std::setprecision(8) << " G4ParameterisationPolyconeZ "
585 << copyNo << G4endl
586 << " Position: " << origin << " - Width: " << fwidth
587 << " - Axis: " << faxis << G4endl;
588 }
589#endif
590}
CLHEP::Hep3Vector G4ThreeVector
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.0) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaxParameter()

G4double G4ParameterisationPolyconeZ::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 429 of file G4ParameterisationPolycone.cc.

430{
431 return std::abs (fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
432 -fOrigParamMother->Z_values[0]);
433}

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