#include <G4ParameterisationPolyhedra.hh>
Inheritance diagram for G4ParameterisationPolyhedraPhi:
Public Member Functions | |
G4ParameterisationPolyhedraPhi (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType) | |
~G4ParameterisationPolyhedraPhi () | |
void | CheckParametersValidity () |
G4double | GetMaxParameter () const |
void | ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const |
void | ComputeDimensions (G4Polyhedra &phedra, const G4int copyNo, const G4VPhysicalVolume *physVol) const |
Definition at line 134 of file G4ParameterisationPolyhedra.hh.
G4ParameterisationPolyhedraPhi::G4ParameterisationPolyhedraPhi | ( | EAxis | axis, | |
G4int | nCopies, | |||
G4double | offset, | |||
G4double | step, | |||
G4VSolid * | motherSolid, | |||
DivisionType | divType | |||
) |
Definition at line 277 of file G4ParameterisationPolyhedra.cc.
References G4VDivisionParameterisation::CalculateWidth(), CheckParametersValidity(), DivWIDTH, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fnDiv, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Polyhedra::GetEndPhi(), G4Polyhedra::GetNumSide(), G4Polyhedra::GetStartPhi(), G4VDivisionParameterisation::SetType(), and G4VDivisionParameterisation::verbose.
00280 : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType ) 00281 { 00282 CheckParametersValidity(); 00283 SetType( "DivisionPolyhedraPhi" ); 00284 00285 G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid); 00286 G4double deltaPhi = msol->GetEndPhi() - msol->GetStartPhi(); 00287 00288 if( divType == DivWIDTH ) 00289 { 00290 fnDiv = msol->GetNumSide(); 00291 } 00292 00293 fwidth = CalculateWidth( deltaPhi, fnDiv, 0.0 ); 00294 00295 #ifdef G4DIVDEBUG 00296 if( verbose >= 1 ) 00297 { 00298 G4cout << " G4ParameterisationPolyhedraPhi - # divisions " << fnDiv 00299 << " = " << nDiv << G4endl 00300 << " Offset " << foffset << " = " << offset << G4endl 00301 << " Width " << fwidth << " = " << width << G4endl; 00302 } 00303 #endif 00304 }
G4ParameterisationPolyhedraPhi::~G4ParameterisationPolyhedraPhi | ( | ) |
void G4ParameterisationPolyhedraPhi::CheckParametersValidity | ( | ) | [virtual] |
Reimplemented from G4VDivisionParameterisation.
Definition at line 319 of file G4ParameterisationPolyhedra.cc.
References G4VDivisionParameterisation::CheckParametersValidity(), DivNDIVandWIDTH, DivWIDTH, FatalException, G4VDivisionParameterisation::fDivisionType, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fnDiv, G4VDivisionParameterisation::foffset, G4endl, G4Exception(), G4VSolid::GetName(), G4Polyhedra::GetOriginalParameters(), JustWarning, and G4PolyhedraHistorical::numSide.
Referenced by G4ParameterisationPolyhedraPhi().
00320 { 00321 G4VDivisionParameterisation::CheckParametersValidity(); 00322 00323 G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid); 00324 00325 if( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH ) 00326 { 00327 std::ostringstream message; 00328 message << "In solid " << msol->GetName() << G4endl 00329 << " Division along PHI will be done splitting " 00330 << "in the defined numSide." << G4endl 00331 << "WIDTH will not be used !"; 00332 G4Exception("G4ParameterisationPolyhedraPhi::CheckParametersValidity()", 00333 "GeomDiv1001", JustWarning, message); 00334 } 00335 if( foffset != 0. ) 00336 { 00337 std::ostringstream message; 00338 message << "In solid " << msol->GetName() << G4endl 00339 << "Division along PHI will be done splitting " 00340 << "in the defined numSide." << G4endl 00341 << "OFFSET will not be used !"; 00342 G4Exception("G4ParameterisationPolyhedraPhi::CheckParametersValidity()", 00343 "GeomDiv1001", JustWarning, message); 00344 } 00345 00346 G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters(); 00347 00348 if( origparamMother->numSide != fnDiv && fDivisionType != DivWIDTH) 00349 { 00350 std::ostringstream message; 00351 message << "Configuration not supported." << G4endl 00352 << "Division along PHI will be done splitting in the defined" 00353 << G4endl 00354 << "numSide, i.e, the number of division would be :" 00355 << origparamMother->numSide << " instead of " << fnDiv << " !"; 00356 G4Exception("G4ParameterisationPolyhedraPhi::CheckParametersValidity()", 00357 "GeomDiv0001", FatalException, message); 00358 } 00359 }
void G4ParameterisationPolyhedraPhi::ComputeDimensions | ( | G4Polyhedra & | phedra, | |
const G4int | copyNo, | |||
const G4VPhysicalVolume * | physVol | |||
) | const [virtual] |
Reimplemented from G4VPVParameterisation.
Definition at line 400 of file G4ParameterisationPolyhedra.cc.
References G4VSolid::DumpInfo(), G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Polyhedra::GetOriginalParameters(), G4PolyhedraHistorical::numSide, G4PolyhedraHistorical::Opening_angle, G4Polyhedra::Reset(), G4Polyhedra::SetOriginalParameters(), G4PolyhedraHistorical::Start_angle, and G4VDivisionParameterisation::verbose.
00402 { 00403 G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid); 00404 00405 G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters(); 00406 G4PolyhedraHistorical origparam( *origparamMother ); 00407 00408 origparam.numSide = 1; 00409 origparam.Start_angle = origparamMother->Start_angle; 00410 origparam.Opening_angle = fwidth; 00411 00412 phedra.SetOriginalParameters(&origparam); // copy values & transfer pointers 00413 phedra.Reset(); // reset to new solid parameters 00414 00415 #ifdef G4DIVDEBUG 00416 if( verbose >= 2 ) 00417 { 00418 G4cout << "G4ParameterisationPolyhedraPhi::ComputeDimensions():" << G4endl; 00419 phedra.DumpInfo(); 00420 } 00421 #endif 00422 }
void G4ParameterisationPolyhedraPhi::ComputeTransformation | ( | const G4int | copyNo, | |
G4VPhysicalVolume * | physVol | |||
) | const [virtual] |
Implements G4VDivisionParameterisation.
Definition at line 364 of file G4ParameterisationPolyhedra.cc.
References G4VDivisionParameterisation::ChangeRotMatrix(), G4VDivisionParameterisation::faxis, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4VPhysicalVolume::SetTranslation(), and G4VDivisionParameterisation::verbose.
00365 { 00366 //----- translation 00367 G4ThreeVector origin(0.,0.,0.); 00368 //----- set translation 00369 physVol->SetTranslation( origin ); 00370 00371 //----- calculate rotation matrix (so that all volumes point to the centre) 00372 G4double posi = copyNo*fwidth; 00373 00374 #ifdef G4DIVDEBUG 00375 if( verbose >= 2 ) 00376 { 00377 G4cout << " G4ParameterisationPolyhedraPhi - position: " << posi/deg 00378 << G4endl 00379 << " copyNo: " << copyNo 00380 << " - fwidth: " << fwidth/deg << G4endl; 00381 } 00382 #endif 00383 00384 ChangeRotMatrix( physVol, -posi ); 00385 00386 #ifdef G4DIVDEBUG 00387 if( verbose >= 2 ) 00388 { 00389 G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraPhi " << copyNo 00390 << G4endl 00391 << " Position: " << origin << " - Width: " << fwidth 00392 << " - Axis: " << faxis << G4endl; 00393 } 00394 #endif 00395 }
G4double G4ParameterisationPolyhedraPhi::GetMaxParameter | ( | ) | const [virtual] |
Implements G4VDivisionParameterisation.
Definition at line 312 of file G4ParameterisationPolyhedra.cc.
References G4VDivisionParameterisation::fmotherSolid, G4Polyhedra::GetEndPhi(), and G4Polyhedra::GetStartPhi().
00313 { 00314 G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid); 00315 return msol->GetEndPhi() - msol->GetStartPhi(); 00316 }