#include <G4ParameterisationTubs.hh>
Inheritance diagram for G4ParameterisationTubsPhi:
Public Member Functions | |
G4ParameterisationTubsPhi (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType) | |
~G4ParameterisationTubsPhi () | |
G4double | GetMaxParameter () const |
void | ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const |
void | ComputeDimensions (G4Tubs &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const |
Definition at line 117 of file G4ParameterisationTubs.hh.
G4ParameterisationTubsPhi::G4ParameterisationTubsPhi | ( | EAxis | axis, | |
G4int | nCopies, | |||
G4double | offset, | |||
G4double | step, | |||
G4VSolid * | motherSolid, | |||
DivisionType | divType | |||
) |
Definition at line 183 of file G4ParameterisationTubs.cc.
References G4VDivisionParameterisation::CalculateNDiv(), G4VDivisionParameterisation::CalculateWidth(), G4VDivisionParameterisation::CheckParametersValidity(), DivNDIV, DivWIDTH, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fnDiv, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Tubs::GetDeltaPhiAngle(), G4VDivisionParameterisation::SetType(), and G4VDivisionParameterisation::verbose.
00186 : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType ) 00187 { 00188 CheckParametersValidity(); 00189 SetType( "DivisionTubsPhi" ); 00190 00191 G4Tubs* msol = (G4Tubs*)(fmotherSolid); 00192 if( divType == DivWIDTH ) 00193 { 00194 fnDiv = CalculateNDiv( msol->GetDeltaPhiAngle(), width, offset ); 00195 } 00196 else if( divType == DivNDIV ) 00197 { 00198 fwidth = CalculateWidth( msol->GetDeltaPhiAngle(), nDiv, offset ); 00199 } 00200 00201 #ifdef G4DIVDEBUG 00202 if( verbose >= 1 ) 00203 { 00204 G4cout << " G4ParameterisationTubsPhi no divisions " << fnDiv << " = " 00205 << nDiv << G4endl 00206 << " Offset " << foffset << " = " << offset << G4endl 00207 << " Width " << fwidth << " = " << width << G4endl; 00208 } 00209 #endif 00210 }
G4ParameterisationTubsPhi::~G4ParameterisationTubsPhi | ( | ) |
void G4ParameterisationTubsPhi::ComputeDimensions | ( | G4Tubs & | tubs, | |
const G4int | copyNo, | |||
const G4VPhysicalVolume * | physVol | |||
) | const [virtual] |
Reimplemented from G4VPVParameterisation.
Definition at line 262 of file G4ParameterisationTubs.cc.
References G4VSolid::DumpInfo(), G4VDivisionParameterisation::fhgap, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Tubs::GetInnerRadius(), G4Tubs::GetOuterRadius(), G4Tubs::GetStartPhiAngle(), G4Tubs::GetZHalfLength(), G4Tubs::SetDeltaPhiAngle(), G4Tubs::SetInnerRadius(), G4Tubs::SetOuterRadius(), G4Tubs::SetStartPhiAngle(), G4Tubs::SetZHalfLength(), and G4VDivisionParameterisation::verbose.
00264 { 00265 G4Tubs* msol = (G4Tubs*)(fmotherSolid); 00266 00267 G4double pRMin = msol->GetInnerRadius(); 00268 G4double pRMax = msol->GetOuterRadius(); 00269 G4double pDz = msol->GetZHalfLength(); 00270 //----- already rotated in 'ComputeTransformation' 00271 G4double pSPhi = msol->GetStartPhiAngle() + fhgap; 00272 G4double pDPhi = fwidth - 2.*fhgap; 00273 00274 tubs.SetInnerRadius( pRMin ); 00275 tubs.SetOuterRadius( pRMax ); 00276 tubs.SetZHalfLength( pDz ); 00277 tubs.SetStartPhiAngle( pSPhi, false ); 00278 tubs.SetDeltaPhiAngle( pDPhi ); 00279 00280 #ifdef G4DIVDEBUG 00281 if( verbose >= 2 ) 00282 { 00283 G4cout << " G4ParameterisationTubsPhi::ComputeDimensions" << G4endl 00284 << " pSPhi: " << pSPhi << " - pDPhi: " << pDPhi << G4endl; 00285 tubs.DumpInfo(); 00286 } 00287 #endif 00288 }
void G4ParameterisationTubsPhi::ComputeTransformation | ( | const G4int | copyNo, | |
G4VPhysicalVolume * | physVol | |||
) | const [virtual] |
Implements G4VDivisionParameterisation.
Definition at line 227 of file G4ParameterisationTubs.cc.
References G4VDivisionParameterisation::ChangeRotMatrix(), G4VDivisionParameterisation::faxis, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4VPhysicalVolume::SetTranslation(), and G4VDivisionParameterisation::verbose.
00228 { 00229 //----- translation 00230 G4ThreeVector origin(0.,0.,0.); 00231 //----- set translation 00232 physVol->SetTranslation( origin ); 00233 00234 //----- calculate rotation matrix (so that all volumes point to the centre) 00235 G4double posi = foffset + copyNo*fwidth; 00236 00237 #ifdef G4DIVDEBUG 00238 if( verbose >= 2 ) 00239 { 00240 G4cout << " G4ParameterisationTubsPhi - position: " << posi/deg << G4endl 00241 << " copyNo: " << copyNo << " - foffset: " << foffset/deg 00242 << " - fwidth: " << fwidth/deg << G4endl; 00243 } 00244 #endif 00245 00246 ChangeRotMatrix( physVol, -posi ); 00247 00248 #ifdef G4DIVDEBUG 00249 if( verbose >= 2 ) 00250 { 00251 G4cout << std::setprecision(8) << " G4ParameterisationTubsPhi " << copyNo 00252 << G4endl 00253 << " Position: " << origin << " - Width: " << fwidth 00254 << " - Axis: " << faxis << G4endl; 00255 } 00256 #endif 00257 }
G4double G4ParameterisationTubsPhi::GetMaxParameter | ( | ) | const [virtual] |
Implements G4VDivisionParameterisation.
Definition at line 218 of file G4ParameterisationTubs.cc.
References G4VDivisionParameterisation::fmotherSolid, and G4Tubs::GetDeltaPhiAngle().
00219 { 00220 G4Tubs* msol = (G4Tubs*)(fmotherSolid); 00221 return msol->GetDeltaPhiAngle(); 00222 }