Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ParameterisationPolyhedraZ Class Reference

#include <G4ParameterisationPolyhedra.hh>

Inheritance diagram for G4ParameterisationPolyhedraZ:
Collaboration 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 *)
 
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 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
G4RotationMatrix
fRot = 0
 
static const G4int verbose = 5
 

Detailed Description

Definition at line 188 of file G4ParameterisationPolyhedra.hh.

Constructor & Destructor Documentation

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

Definition at line 434 of file G4ParameterisationPolyhedra.cc.

437  : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType ),
438  fNSegment(0),
439  fOrigParamMother(((G4Polyhedra*)fmotherSolid)->GetOriginalParameters())
440 {
442  SetType( "DivisionPolyhedraZ" );
443 
444  if( divType == DivWIDTH )
445  {
446  fnDiv =
447  CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
448  - fOrigParamMother->Z_values[0] , width, offset );
449  }
450  else if( divType == DivNDIV )
451  {
452  fwidth =
453  CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
454  - fOrigParamMother->Z_values[0] , nDiv, offset );
455  }
456 
457 #ifdef G4DIVDEBUG
458  if( verbose >= 1 )
459  {
460  G4cout << " G4ParameterisationPolyhedraZ - # divisions " << fnDiv << " = "
461  << nDiv << G4endl
462  << " Offset " << foffset << " = " << offset << G4endl
463  << " Width " << fwidth << " = " << width << G4endl;
464  }
465 #endif
466 }
void SetType(const G4String &type)
G4VParameterisationPolyhedra(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
G4GLOB_DLL std::ostream G4cout
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

G4ParameterisationPolyhedraZ::~G4ParameterisationPolyhedraZ ( )

Definition at line 469 of file G4ParameterisationPolyhedra.cc.

470 {
471 }

Member Function Documentation

void G4ParameterisationPolyhedraZ::CheckParametersValidity ( )
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 518 of file G4ParameterisationPolyhedra.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VPVParameterisation.

Definition at line 657 of file G4ParameterisationPolyhedra.cc.

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

Here is the call graph for this function:

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

Implements G4VDivisionParameterisation.

Definition at line 608 of file G4ParameterisationPolyhedra.cc.

609 {
610  if ( fDivisionType == DivNDIV ) {
611  // The position of the centre of copyNo-th mother polycone segment
612  G4double posi = ( fOrigParamMother->Z_values[copyNo]
613  + fOrigParamMother->Z_values[copyNo+1])/2;
614  physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
615  }
616 
618  // The position of the centre of copyNo-th division
619 
620  G4double posi = fOrigParamMother->Z_values[0];
621 
622  if ( ! fReflectedSolid )
623  posi += foffset + (2*copyNo + 1) * fwidth/2.;
624  else
625  posi -= foffset + (2*copyNo + 1) * fwidth/2.;
626 
627  physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
628  }
629 
630  //----- calculate rotation matrix: unit
631 
632 #ifdef G4DIVDEBUG
633  if( verbose >= 2 )
634  {
635  G4cout << " G4ParameterisationPolyhedraZ - position: " << posi << G4endl
636  << " copyNo: " << copyNo << " - foffset: " << foffset/deg
637  << " - fwidth: " << fwidth/deg << G4endl;
638  }
639 #endif
640 
641  ChangeRotMatrix( physVol );
642 
643 #ifdef G4DIVDEBUG
644  if( verbose >= 2 )
645  {
646  G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraZ "
647  << copyNo << G4endl
648  << " Position: " << origin << " - Width: " << fwidth
649  << " - Axis: " << faxis << G4endl;
650  }
651 #endif
652 }
CLHEP::Hep3Vector G4ThreeVector
G4GLOB_DLL std::ostream G4cout
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.) const
void SetTranslation(const G4ThreeVector &v)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static constexpr double deg
Definition: G4SIunits.hh:152

Here is the call graph for this function:

G4double G4ParameterisationPolyhedraZ::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 511 of file G4ParameterisationPolyhedra.cc.

512 {
513  return std::abs (fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
514  -fOrigParamMother->Z_values[0]);
515 }

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