#include <G4ParameterisationCons.hh>
Inheritance diagram for G4ParameterisationConsPhi:
Public Member Functions | |
G4ParameterisationConsPhi (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType) | |
~G4ParameterisationConsPhi () | |
G4double | GetMaxParameter () const |
void | ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const |
void | ComputeDimensions (G4Cons &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const |
Definition at line 118 of file G4ParameterisationCons.hh.
G4ParameterisationConsPhi::G4ParameterisationConsPhi | ( | EAxis | axis, | |
G4int | nCopies, | |||
G4double | offset, | |||
G4double | step, | |||
G4VSolid * | motherSolid, | |||
DivisionType | divType | |||
) |
Definition at line 216 of file G4ParameterisationCons.cc.
References G4VDivisionParameterisation::CalculateNDiv(), G4VDivisionParameterisation::CalculateWidth(), G4VDivisionParameterisation::CheckParametersValidity(), DivNDIV, DivWIDTH, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fnDiv, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Cons::GetDeltaPhiAngle(), G4VDivisionParameterisation::SetType(), and G4VDivisionParameterisation::verbose.
00219 : G4VParameterisationCons( axis, nDiv, width, offset, msolid, divType ) 00220 { 00221 CheckParametersValidity(); 00222 SetType( "DivisionConsPhi" ); 00223 00224 G4Cons* msol = (G4Cons*)(fmotherSolid); 00225 if( divType == DivWIDTH ) 00226 { 00227 fnDiv = CalculateNDiv( msol->GetDeltaPhiAngle(), width, offset ); 00228 } 00229 else if( divType == DivNDIV ) 00230 { 00231 fwidth = CalculateWidth( msol->GetDeltaPhiAngle(), nDiv, offset ); 00232 } 00233 00234 #ifdef G4DIVDEBUG 00235 if( verbose >= 1 ) 00236 { 00237 G4cout << " G4ParameterisationConsPhi no divisions " << fnDiv << " = " 00238 << nDiv << G4endl 00239 << " Offset " << foffset << " = " << offset << G4endl 00240 << " Width " << fwidth << " = " << width << G4endl; 00241 } 00242 #endif 00243 }
G4ParameterisationConsPhi::~G4ParameterisationConsPhi | ( | ) |
void G4ParameterisationConsPhi::ComputeDimensions | ( | G4Cons & | tubs, | |
const G4int | copyNo, | |||
const G4VPhysicalVolume * | physVol | |||
) | const [virtual] |
Reimplemented from G4VPVParameterisation.
Definition at line 287 of file G4ParameterisationCons.cc.
References G4VSolid::DumpInfo(), G4VDivisionParameterisation::fhgap, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Cons::GetInnerRadiusMinusZ(), G4Cons::GetInnerRadiusPlusZ(), G4Cons::GetOuterRadiusMinusZ(), G4Cons::GetOuterRadiusPlusZ(), G4Cons::GetStartPhiAngle(), G4Cons::GetZHalfLength(), G4Cons::SetDeltaPhiAngle(), G4Cons::SetInnerRadiusMinusZ(), G4Cons::SetInnerRadiusPlusZ(), G4Cons::SetOuterRadiusMinusZ(), G4Cons::SetOuterRadiusPlusZ(), G4Cons::SetStartPhiAngle(), G4Cons::SetZHalfLength(), and G4VDivisionParameterisation::verbose.
00289 { 00290 G4Cons* msol = (G4Cons*)(fmotherSolid); 00291 00292 G4double pRMin1 = msol->GetInnerRadiusMinusZ(); 00293 G4double pRMax1 = msol->GetOuterRadiusMinusZ(); 00294 G4double pRMin2 = msol->GetInnerRadiusPlusZ(); 00295 G4double pRMax2 = msol->GetOuterRadiusPlusZ(); 00296 G4double pDz = msol->GetZHalfLength(); 00297 00298 //- already rotated double pSPhi = foffset + copyNo*fwidth; 00299 G4double pSPhi = foffset + msol->GetStartPhiAngle() + fhgap; 00300 G4double pDPhi = fwidth - 2.*fhgap; 00301 00302 cons.SetInnerRadiusMinusZ( pRMin1 ); 00303 cons.SetOuterRadiusMinusZ( pRMax1 ); 00304 cons.SetInnerRadiusPlusZ( pRMin2 ); 00305 cons.SetOuterRadiusPlusZ( pRMax2 ); 00306 cons.SetZHalfLength( pDz ); 00307 cons.SetStartPhiAngle( pSPhi, false ); 00308 cons.SetDeltaPhiAngle( pDPhi ); 00309 00310 #ifdef G4DIVDEBUG 00311 if( verbose >= 2 ) 00312 { 00313 G4cout << " G4ParameterisationConsPhi::ComputeDimensions" << G4endl 00314 << " pSPhi: " << pSPhi << " - pDPhi: " << pDPhi << G4endl; 00315 if( verbose >= 4 ) cons.DumpInfo(); 00316 } 00317 #endif 00318 }
void G4ParameterisationConsPhi::ComputeTransformation | ( | const G4int | copyNo, | |
G4VPhysicalVolume * | physVol | |||
) | const [virtual] |
Implements G4VDivisionParameterisation.
Definition at line 260 of file G4ParameterisationCons.cc.
References G4VDivisionParameterisation::ChangeRotMatrix(), G4VDivisionParameterisation::faxis, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4VPhysicalVolume::SetTranslation(), and G4VDivisionParameterisation::verbose.
00261 { 00262 //----- translation 00263 G4ThreeVector origin(0.,0.,0.); 00264 //----- set translation 00265 physVol->SetTranslation( origin ); 00266 00267 //----- calculate rotation matrix (so that all volumes point to the centre) 00268 G4double posi = foffset + copyNo*fwidth; 00269 00270 #ifdef G4DIVDEBUG 00271 if( verbose >= 2 ) 00272 { 00273 G4cout << " G4ParameterisationConsPhi - position: " << posi/deg << G4endl 00274 << " Origin: " << origin << " copyNo: " << copyNo 00275 << " - foffset: " << foffset/deg 00276 << " - fwidth: " << fwidth/deg << G4endl 00277 << " - Axis: " << faxis << G4endl; 00278 } 00279 #endif 00280 00281 ChangeRotMatrix( physVol, -posi ); 00282 }
G4double G4ParameterisationConsPhi::GetMaxParameter | ( | ) | const [virtual] |
Implements G4VDivisionParameterisation.
Definition at line 251 of file G4ParameterisationCons.cc.
References G4VDivisionParameterisation::fmotherSolid, and G4Cons::GetDeltaPhiAngle().
00252 { 00253 G4Cons* msol = (G4Cons*)(fmotherSolid); 00254 return msol->GetDeltaPhiAngle(); 00255 }