Geant4 11.2.2
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 () override
 
void CheckParametersValidity () override
 
G4double GetMaxParameter () const override
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const override
 
void ComputeDimensions (G4Polycone &pcone, const G4int copyNo, const G4VPhysicalVolume *physVol) const override
 
- Public Member Functions inherited from G4VParameterisationPolycone
 G4VParameterisationPolycone (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
 ~G4VParameterisationPolycone () 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 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 345 of file G4ParameterisationPolycone.cc.

349 : G4VParameterisationPolycone( axis, nDiv, width, offset, msolid, divType ),
350 fOrigParamMother(((G4Polycone*)fmotherSolid)->GetOriginalParameters())
351{
352
354 SetType( "DivisionPolyconeZ" );
355
356 if( divType == DivWIDTH )
357 {
358 fnDiv =
359 CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
360 - fOrigParamMother->Z_values[0] , width, offset );
361 }
362 else if( divType == DivNDIV )
363 {
364 fwidth =
365 CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
366 - fOrigParamMother->Z_values[0] , nDiv, offset );
367 }
368
369#ifdef G4DIVDEBUG
370 if( verbose >= 1 )
371 {
372 G4cout << " G4ParameterisationPolyconeZ - # divisions " << fnDiv << " = "
373 << nDiv << G4endl
374 << " Offset " << foffset << " = " << offset << G4endl
375 << " Width " << fwidth << " = " << width << G4endl;
376 }
377#endif
378}
#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
G4VParameterisationPolycone(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)

◆ ~G4ParameterisationPolyconeZ()

G4ParameterisationPolyconeZ::~G4ParameterisationPolyconeZ ( )
overridedefault

Member Function Documentation

◆ CheckParametersValidity()

void G4ParameterisationPolyconeZ::CheckParametersValidity ( )
overridevirtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 428 of file G4ParameterisationPolycone.cc.

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

◆ ComputeDimensions()

void G4ParameterisationPolyconeZ::ComputeDimensions ( G4Polycone & pcone,
const G4int copyNo,
const G4VPhysicalVolume * physVol ) const
overridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 587 of file G4ParameterisationPolycone.cc.

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

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 534 of file G4ParameterisationPolycone.cc.

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

◆ GetMaxParameter()

G4double G4ParameterisationPolyconeZ::GetMaxParameter ( ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 421 of file G4ParameterisationPolycone.cc.

422{
423 return std::abs (fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
424 -fOrigParamMother->Z_values[0]);
425}

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