Geant4 11.2.2
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 () override
 
void CheckParametersValidity () override
 
G4double GetMaxParameter () const override
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const override
 
void ComputeDimensions (G4Polyhedra &phedra, const G4int copyNo, const G4VPhysicalVolume *physVol) const override
 
- Public Member Functions inherited from G4VParameterisationPolyhedra
 G4VParameterisationPolyhedra (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
 ~G4VParameterisationPolyhedra () override
 
- Public Member Functions inherited from G4VDivisionParameterisation
 G4VDivisionParameterisation (EAxis axis, G4int nDiv, G4double width, G4double offset, DivisionType divType, G4VSolid *motherSolid=nullptr)
 
 ~G4VDivisionParameterisation () override
 
G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *) override
 
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 ()=default
 
virtual ~G4VPVParameterisation ()=default
 
virtual G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr)
 
virtual G4bool IsNested () const
 
virtual G4VVolumeMaterialScannerGetMaterialScanner ()
 

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
 
void CheckOffset (G4double maxPar)
 
void CheckNDivAndWidth (G4double maxPar)
 
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 181 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 426 of file G4ParameterisationPolyhedra.cc.

430 : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType ),
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:67
G4GLOB_DLL std::ostream G4cout
void SetType(const G4String &type)
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
G4VParameterisationPolyhedra(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)

◆ ~G4ParameterisationPolyhedraZ()

G4ParameterisationPolyhedraZ::~G4ParameterisationPolyhedraZ ( )
overridedefault

Member Function Documentation

◆ CheckParametersValidity()

void G4ParameterisationPolyhedraZ::CheckParametersValidity ( )
overridevirtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 508 of file G4ParameterisationPolyhedra.cc.

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

Referenced by G4ParameterisationPolyhedraZ().

◆ ComputeDimensions()

void G4ParameterisationPolyhedraZ::ComputeDimensions ( G4Polyhedra & phedra,
const G4int copyNo,
const G4VPhysicalVolume * physVol ) const
overridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 664 of file G4ParameterisationPolyhedra.cc.

667{
668 // Define division solid
669 //
670 G4PolyhedraHistorical origparam;
671 G4int nz = 2;
672 origparam.Num_z_planes = nz;
673 origparam.numSide = fOrigParamMother->numSide;
674 origparam.Start_angle = fOrigParamMother->Start_angle;
675 origparam.Opening_angle = fOrigParamMother->Opening_angle;
676
677 // Define division solid z sections
678 //
679 origparam.Z_values = new G4double[nz];
680 origparam.Rmin = new G4double[nz];
681 origparam.Rmax = new G4double[nz];
682 origparam.Z_values[0] = - fwidth/2.;
683 origparam.Z_values[1] = fwidth/2.;
684
685 if ( fDivisionType == DivNDIV )
686 {
687 // The position of the centre of copyNo-th mother polycone segment
688 //
689 G4double posi = ( fOrigParamMother->Z_values[copyNo]
690 + fOrigParamMother->Z_values[copyNo+1])/2;
691
692 origparam.Z_values[0] = fOrigParamMother->Z_values[copyNo] - posi;
693 origparam.Z_values[1] = fOrigParamMother->Z_values[copyNo+1] - posi;
694 origparam.Rmin[0] = fOrigParamMother->Rmin[copyNo];
695 origparam.Rmin[1] = fOrigParamMother->Rmin[copyNo+1];
696 origparam.Rmax[0] = fOrigParamMother->Rmax[copyNo];
697 origparam.Rmax[1] = fOrigParamMother->Rmax[copyNo+1];
698 }
699
701 {
702 if ( !fReflectedSolid )
703 {
704 origparam.Z_values[0] = -fwidth/2.;
705 origparam.Z_values[1] = fwidth/2.;
706
707 // The position of the centre of copyNo-th division
708 //
709 G4double posi = fOrigParamMother->Z_values[0]
710 + foffset + (2*copyNo + 1) * fwidth/2.;
711
712 // The first and last z sides z values
713 G4double zstart = posi - fwidth/2.;
714 G4double zend = posi + fwidth/2.;
715 origparam.Rmin[0] = GetRmin(zstart, fNSegment);
716 origparam.Rmax[0] = GetRmax(zstart, fNSegment);
717 origparam.Rmin[1] = GetRmin(zend, fNSegment);
718 origparam.Rmax[1] = GetRmax(zend, fNSegment);
719 }
720 else
721 {
722 origparam.Z_values[0] = fwidth/2.;
723 origparam.Z_values[1] = -fwidth/2.;
724
725 // The position of the centre of copyNo-th division
726 //
727 G4double posi = fOrigParamMother->Z_values[0]
728 - ( foffset + (2*copyNo + 1) * fwidth/2.);
729
730 // The first and last z sides z values
731 //
732 G4double zstart = posi + fwidth/2.;
733 G4double zend = posi - fwidth/2.;
734 origparam.Rmin[0] = GetRmin(zstart, fNSegment);
735 origparam.Rmax[0] = GetRmax(zstart, fNSegment);
736 origparam.Rmin[1] = GetRmin(zend, fNSegment);
737 origparam.Rmax[1] = GetRmax(zend, fNSegment);
738 }
739
740 // It can happen due to rounding errors
741 //
742 if ( origparam.Rmin[0] < 0.0 ) origparam.Rmin[0] = 0.0;
743 if ( origparam.Rmin[nz-1] < 0.0 ) origparam.Rmin[1] = 0.0;
744 }
745
746 phedra.SetOriginalParameters(&origparam); // copy values & transfer pointers
747 phedra.Reset(); // reset to new solid parameters
748
749#ifdef G4DIVDEBUG
750 if( verbose >= 2 )
751 {
752 G4cout << "G4ParameterisationPolyhedraZ::ComputeDimensions()" << G4endl
753 << "-- Parametrised phedra copy-number: " << copyNo << G4endl;
754 phedra.DumpInfo();
755 }
756#endif
757}
void SetOriginalParameters(G4PolyhedraHistorical *pars)
G4bool Reset()
void DumpInfo() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 611 of file G4ParameterisationPolyhedra.cc.

613{
614 G4double posi;
615 if ( fDivisionType == DivNDIV )
616 {
617 // The position of the centre of copyNo-th mother polycone segment
618
619 posi = ( fOrigParamMother->Z_values[copyNo]
620 + fOrigParamMother->Z_values[copyNo+1])/2;
621 physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
622 }
623
625 {
626 // The position of the centre of copyNo-th division
627
628 posi = fOrigParamMother->Z_values[0];
629
630 if ( !fReflectedSolid )
631 posi += foffset + (2*copyNo + 1) * fwidth/2.;
632 else
633 posi -= foffset + (2*copyNo + 1) * fwidth/2.;
634
635 physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
636 }
637
638 //----- calculate rotation matrix: unit
639
640#ifdef G4DIVDEBUG
641 if( verbose >= 2 )
642 {
643 G4cout << " G4ParameterisationPolyhedraZ - position: " << posi << G4endl
644 << " copyNo: " << copyNo << " - foffset: " << foffset/CLHEP::deg
645 << " - fwidth: " << fwidth/CLHEP::deg << G4endl;
646 }
647#endif
648
649 ChangeRotMatrix( physVol );
650
651#ifdef G4DIVDEBUG
652 if( verbose >= 2 )
653 {
654 G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraZ "
655 << copyNo << G4endl
656 << " Position: (0,0,0) - Width: " << fwidth
657 << " - Axis: " << faxis << G4endl;
658 }
659#endif
660}
CLHEP::Hep3Vector G4ThreeVector
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.0) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaxParameter()

G4double G4ParameterisationPolyhedraZ::GetMaxParameter ( ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 501 of file G4ParameterisationPolyhedra.cc.

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

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