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

#include <G4ParameterisationPolyhedra.hh>

+ Inheritance diagram for G4ParameterisationPolyhedraZ:

Public Member Functions

 G4ParameterisationPolyhedraZ (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationPolyhedraZ ()
 
void CheckParametersValidity ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Polyhedra &phedra, const G4int copyNo, const G4VPhysicalVolume *physVol) const
 
- Public Member Functions inherited from G4VParameterisationPolyhedra
 G4VParameterisationPolyhedra (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
virtual ~G4VParameterisationPolyhedra ()
 
- Public Member Functions inherited from G4VDivisionParameterisation
 G4VDivisionParameterisation (EAxis axis, G4int nDiv, G4double width, G4double offset, DivisionType divType, G4VSolid *motherSolid=0)
 
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=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 (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.) 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
 
G4double fwidth
 
G4double foffset
 
DivisionType fDivisionType
 
G4VSolidfmotherSolid
 
G4bool fReflectedSolid
 
G4bool fDeleteSolid
 
G4int theVoluFirstCopyNo
 
G4double kCarTolerance
 
G4double fhgap
 
- Static Protected Attributes inherited from G4VDivisionParameterisation
static G4int verbose = 5
 

Detailed Description

Definition at line 183 of file G4ParameterisationPolyhedra.hh.

Constructor & Destructor Documentation

◆ G4ParameterisationPolyhedraZ()

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

Definition at line 425 of file G4ParameterisationPolyhedra.cc.

429 : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType ),
430 fNSegment(0),
431 fOrigParamMother(((G4Polyhedra*)fmotherSolid)->GetOriginalParameters())
432{
434 SetType( "DivisionPolyhedraZ" );
435
436 if( divType == DivWIDTH )
437 {
438 fnDiv =
439 CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
440 - fOrigParamMother->Z_values[0] , width, offset );
441 }
442 else if( divType == DivNDIV )
443 {
444 fwidth =
445 CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
446 - fOrigParamMother->Z_values[0] , nDiv, offset );
447 }
448
449#ifdef G4DIVDEBUG
450 if( verbose >= 1 )
451 {
452 G4cout << " G4ParameterisationPolyhedraZ - # divisions " << fnDiv << " = "
453 << nDiv << G4endl
454 << " Offset " << foffset << " = " << offset << G4endl
455 << " Width " << fwidth << " = " << width << G4endl;
456 }
457#endif
458}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void SetType(const G4String &type)
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const

◆ ~G4ParameterisationPolyhedraZ()

G4ParameterisationPolyhedraZ::~G4ParameterisationPolyhedraZ ( )

Definition at line 461 of file G4ParameterisationPolyhedra.cc.

462{
463}

Member Function Documentation

◆ CheckParametersValidity()

void G4ParameterisationPolyhedraZ::CheckParametersValidity ( )
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 510 of file G4ParameterisationPolyhedra.cc.

511{
513
514 // Division will be following the mother polyhedra segments
515 if( fDivisionType == DivNDIV ) {
516 if( fOrigParamMother->Num_z_planes-1 != fnDiv ) {
517 std::ostringstream message;
518 message << "Configuration not supported." << G4endl
519 << "Division along Z will be done splitting in the defined"
520 << G4endl
521 << "Z planes, i.e, the number of division would be :"
522 << fOrigParamMother->Num_z_planes-1 << " instead of "
523 << fnDiv << " !";
524 G4Exception("G4ParameterisationPolyhedraZ::CheckParametersValidity()",
525 "GeomDiv0001", FatalException, message);
526 }
527 }
528
529 // Division will be done within one polyhedra segment
530 // with applying given width and offset
532 // Check if divided region does not span over more
533 // than one z segment
534
535 G4int isegstart = -1; // number of the segment containing start position
536 G4int isegend = -1; // number of the segment containing end position
537
538 if ( ! fReflectedSolid ) {
539 // The start/end position of the divided region
540 G4double zstart
541 = fOrigParamMother->Z_values[0] + foffset;
542 G4double zend
543 = fOrigParamMother->Z_values[0] + foffset + fnDiv* fwidth;
544
545 G4int counter = 0;
546 while ( isegend < 0 && counter < fOrigParamMother->Num_z_planes-1 ) {
547 // first segment
548 if ( zstart >= fOrigParamMother->Z_values[counter] &&
549 zstart < fOrigParamMother->Z_values[counter+1] ) {
550 isegstart = counter;
551 }
552 // last segment
553 if ( zend > fOrigParamMother->Z_values[counter] &&
554 zend <= fOrigParamMother->Z_values[counter+1] ) {
555 isegend = counter;
556 }
557 ++counter;
558 }
559 }
560 else {
561 // The start/end position of the divided region
562 G4double zstart
563 = fOrigParamMother->Z_values[0] - foffset;
564 G4double zend
565 = fOrigParamMother->Z_values[0] - ( foffset + fnDiv* fwidth);
566
567 G4int counter = 0;
568 while ( isegend < 0 && counter < fOrigParamMother->Num_z_planes-1 ) {
569 // first segment
570 if ( zstart <= fOrigParamMother->Z_values[counter] &&
571 zstart > fOrigParamMother->Z_values[counter+1] ) {
572 isegstart = counter;
573 }
574 // last segment
575 if ( zend < fOrigParamMother->Z_values[counter] &&
576 zend >= fOrigParamMother->Z_values[counter+1] ) {
577 isegend = counter;
578 }
579 ++counter;
580 }
581 }
582
583 if ( isegstart != isegend ) {
584 std::ostringstream message;
585 message << "Configuration not supported." << G4endl
586 << "Division with user defined width." << G4endl
587 << "Solid " << fmotherSolid->GetName() << G4endl
588 << "Divided region is not between two Z planes.";
589 G4Exception("G4ParameterisationPolyhedraZ::CheckParametersValidity()",
590 "GeomDiv0001", FatalException, message);
591 }
592
593 fNSegment = isegstart;
594 }
595}
@ FatalException
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4String GetName() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Referenced by G4ParameterisationPolyhedraZ().

◆ ComputeDimensions()

void G4ParameterisationPolyhedraZ::ComputeDimensions ( G4Polyhedra phedra,
const G4int  copyNo,
const G4VPhysicalVolume physVol 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 648 of file G4ParameterisationPolyhedra.cc.

651{
652 // Define division solid
653 G4PolyhedraHistorical origparam;
654 G4int nz = 2;
655 origparam.Num_z_planes = nz;
656 origparam.numSide = fOrigParamMother->numSide;
657 origparam.Start_angle = fOrigParamMother->Start_angle;
658 origparam.Opening_angle = fOrigParamMother->Opening_angle;
659
660 // Define division solid z sections
661 origparam.Z_values = new G4double[nz];
662 origparam.Rmin = new G4double[nz];
663 origparam.Rmax = new G4double[nz];
664 origparam.Z_values[0] = - fwidth/2.;
665 origparam.Z_values[1] = fwidth/2.;
666
667 if ( fDivisionType == DivNDIV ) {
668 // The position of the centre of copyNo-th mother polycone segment
669 G4double posi = ( fOrigParamMother->Z_values[copyNo]
670 + fOrigParamMother->Z_values[copyNo+1])/2;
671
672 origparam.Z_values[0] = fOrigParamMother->Z_values[copyNo] - posi;
673 origparam.Z_values[1] = fOrigParamMother->Z_values[copyNo+1] - posi;
674 origparam.Rmin[0] = fOrigParamMother->Rmin[copyNo];
675 origparam.Rmin[1] = fOrigParamMother->Rmin[copyNo+1];
676 origparam.Rmax[0] = fOrigParamMother->Rmax[copyNo];
677 origparam.Rmax[1] = fOrigParamMother->Rmax[copyNo+1];
678 }
679
681 if ( ! fReflectedSolid ) {
682 origparam.Z_values[0] = - fwidth/2.;
683 origparam.Z_values[1] = fwidth/2.;
684
685 // The position of the centre of copyNo-th division
686 G4double posi
687 = fOrigParamMother->Z_values[0] + foffset + (2*copyNo + 1) * fwidth/2.;
688
689 // The first and last z sides z values
690 G4double zstart = posi - fwidth/2.;
691 G4double zend = posi + fwidth/2.;
692 origparam.Rmin[0] = GetRmin(zstart, fNSegment);
693 origparam.Rmax[0] = GetRmax(zstart, fNSegment);
694 origparam.Rmin[1] = GetRmin(zend, fNSegment);
695 origparam.Rmax[1] = GetRmax(zend, fNSegment);
696 }
697 else {
698 origparam.Z_values[0] = fwidth/2.;
699 origparam.Z_values[1] = - fwidth/2.;
700
701 // The position of the centre of copyNo-th division
702 G4double posi
703 = fOrigParamMother->Z_values[0] - ( foffset + (2*copyNo + 1) * fwidth/2.);
704
705 // The first and last z sides z values
706 G4double zstart = posi + fwidth/2.;
707 G4double zend = posi - fwidth/2.;
708 origparam.Rmin[0] = GetRmin(zstart, fNSegment);
709 origparam.Rmax[0] = GetRmax(zstart, fNSegment);
710 origparam.Rmin[1] = GetRmin(zend, fNSegment);
711 origparam.Rmax[1] = GetRmax(zend, fNSegment);
712 }
713
714 // It can happen due to rounding errors
715 if ( origparam.Rmin[0] < 0.0 ) origparam.Rmin[0] = 0.0;
716 if ( origparam.Rmin[nz-1] < 0.0 ) origparam.Rmin[1] = 0.0;
717 }
718
719 phedra.SetOriginalParameters(&origparam); // copy values & transfer pointers
720 phedra.Reset(); // reset to new solid parameters
721
722#ifdef G4DIVDEBUG
723 if( verbose >= 2 )
724 {
725 G4cout << "G4ParameterisationPolyhedraZ::ComputeDimensions()" << G4endl
726 << "-- Parametrised phedra copy-number: " << copyNo << G4endl;
727 phedra.DumpInfo();
728 }
729#endif
730}
void SetOriginalParameters(G4PolyhedraHistorical *pars)
G4bool Reset()
Definition: G4Polyhedra.cc:465
void DumpInfo() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 599 of file G4ParameterisationPolyhedra.cc.

601{
602 if ( fDivisionType == DivNDIV ) {
603 // The position of the centre of copyNo-th mother polycone segment
604 G4double posi = ( fOrigParamMother->Z_values[copyNo]
605 + fOrigParamMother->Z_values[copyNo+1])/2;
606 physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
607 }
608
610 // The position of the centre of copyNo-th division
611
612 G4double posi = fOrigParamMother->Z_values[0];
613
614 if ( ! fReflectedSolid )
615 posi += foffset + (2*copyNo + 1) * fwidth/2.;
616 else
617 posi -= foffset + (2*copyNo + 1) * fwidth/2.;
618
619 physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
620 }
621
622 //----- calculate rotation matrix: unit
623
624#ifdef G4DIVDEBUG
625 if( verbose >= 2 )
626 {
627 G4cout << " G4ParameterisationPolyhedraZ - position: " << posi << G4endl
628 << " copyNo: " << copyNo << " - foffset: " << foffset/deg
629 << " - fwidth: " << fwidth/deg << G4endl;
630 }
631#endif
632
633 ChangeRotMatrix( physVol );
634
635#ifdef G4DIVDEBUG
636 if( verbose >= 2 )
637 {
638 G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraZ "
639 << copyNo << G4endl
640 << " Position: " << origin << " - Width: " << fwidth
641 << " - Axis: " << faxis << G4endl;
642 }
643#endif
644}
CLHEP::Hep3Vector G4ThreeVector
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaxParameter()

G4double G4ParameterisationPolyhedraZ::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 503 of file G4ParameterisationPolyhedra.cc.

504{
505 return std::abs (fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
506 -fOrigParamMother->Z_values[0]);
507}

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