#include <G4ParameterisationPolyhedra.hh>
Inheritance diagram for G4ParameterisationPolyhedraRho:
Public Member Functions | |
G4ParameterisationPolyhedraRho (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType) | |
~G4ParameterisationPolyhedraRho () | |
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 85 of file G4ParameterisationPolyhedra.hh.
G4ParameterisationPolyhedraRho::G4ParameterisationPolyhedraRho | ( | EAxis | axis, | |
G4int | nCopies, | |||
G4double | offset, | |||
G4double | step, | |||
G4VSolid * | motherSolid, | |||
DivisionType | divType | |||
) |
Definition at line 129 of file G4ParameterisationPolyhedra.cc.
References G4VDivisionParameterisation::CalculateNDiv(), G4VDivisionParameterisation::CalculateWidth(), CheckParametersValidity(), DivNDIV, DivWIDTH, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fnDiv, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Polyhedra::GetOriginalParameters(), G4PolyhedraHistorical::Rmax, G4PolyhedraHistorical::Rmin, G4VDivisionParameterisation::SetType(), and G4VDivisionParameterisation::verbose.
00132 : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType ) 00133 { 00134 00135 CheckParametersValidity(); 00136 SetType( "DivisionPolyhedraRho" ); 00137 00138 G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid); 00139 G4PolyhedraHistorical* original_pars = msol->GetOriginalParameters(); 00140 00141 if( divType == DivWIDTH ) 00142 { 00143 fnDiv = CalculateNDiv( original_pars->Rmax[0] 00144 - original_pars->Rmin[0], width, offset ); 00145 } 00146 else if( divType == DivNDIV ) 00147 { 00148 fwidth = CalculateWidth( original_pars->Rmax[0] 00149 - original_pars->Rmin[0], nDiv, offset ); 00150 } 00151 00152 #ifdef G4DIVDEBUG 00153 if( verbose >= 1 ) 00154 { 00155 G4cout << " G4ParameterisationPolyhedraRho - # divisions " << fnDiv 00156 << " = " << nDiv << G4endl 00157 << " Offset " << foffset << " = " << offset << G4endl 00158 << " Width " << fwidth << " = " << width << G4endl; 00159 } 00160 #endif 00161 }
G4ParameterisationPolyhedraRho::~G4ParameterisationPolyhedraRho | ( | ) |
void G4ParameterisationPolyhedraRho::CheckParametersValidity | ( | ) | [virtual] |
Reimplemented from G4VDivisionParameterisation.
Definition at line 169 of file G4ParameterisationPolyhedra.cc.
References G4VDivisionParameterisation::CheckParametersValidity(), DivNDIVandWIDTH, DivWIDTH, G4VDivisionParameterisation::fDivisionType, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::foffset, G4endl, G4Exception(), G4VSolid::GetName(), and JustWarning.
Referenced by G4ParameterisationPolyhedraRho().
00170 { 00171 G4VDivisionParameterisation::CheckParametersValidity(); 00172 00173 G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid); 00174 00175 if( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH ) 00176 { 00177 std::ostringstream message; 00178 message << "In solid " << msol->GetName() << G4endl 00179 << "Division along R will be done with a width " 00180 << "different for each solid section." << G4endl 00181 << "WIDTH will not be used !"; 00182 G4Exception("G4ParameterisationPolyhedraRho::CheckParametersValidity()", 00183 "GeomDiv1001", JustWarning, message); 00184 } 00185 if( foffset != 0. ) 00186 { 00187 std::ostringstream message; 00188 message << "In solid " << msol->GetName() << G4endl 00189 << "Division along R will be done with a width " 00190 << "different for each solid section." << G4endl 00191 << "OFFSET will not be used !"; 00192 G4Exception("G4ParameterisationPolyhedraRho::CheckParametersValidity()", 00193 "GeomDiv1001", JustWarning, message); 00194 } 00195 }
void G4ParameterisationPolyhedraRho::ComputeDimensions | ( | G4Polyhedra & | phedra, | |
const G4int | copyNo, | |||
const G4VPhysicalVolume * | physVol | |||
) | const [virtual] |
Reimplemented from G4VPVParameterisation.
Definition at line 244 of file G4ParameterisationPolyhedra.cc.
References G4VDivisionParameterisation::CalculateWidth(), G4VSolid::DumpInfo(), G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fnDiv, G4VDivisionParameterisation::foffset, G4cout, G4endl, G4Polyhedra::GetOriginalParameters(), G4PolyhedraHistorical::Num_z_planes, G4Polyhedra::Reset(), G4PolyhedraHistorical::Rmax, G4PolyhedraHistorical::Rmin, G4Polyhedra::SetOriginalParameters(), and G4VDivisionParameterisation::verbose.
00246 { 00247 G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid); 00248 00249 G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters(); 00250 G4PolyhedraHistorical origparam( *origparamMother ); 00251 G4int nZplanes = origparamMother->Num_z_planes; 00252 00253 G4double width = 0.; 00254 for( G4int ii = 0; ii < nZplanes; ii++ ) 00255 { 00256 width = CalculateWidth( origparamMother->Rmax[ii] 00257 - origparamMother->Rmin[ii], fnDiv, foffset ); 00258 origparam.Rmin[ii] = origparamMother->Rmin[ii]+foffset+width*copyNo; 00259 origparam.Rmax[ii] = origparamMother->Rmin[ii]+foffset+width*(copyNo+1); 00260 } 00261 00262 phedra.SetOriginalParameters(&origparam); // copy values & transfer pointers 00263 phedra.Reset(); // reset to new solid parameters 00264 00265 #ifdef G4DIVDEBUG 00266 if( verbose >= -2 ) 00267 { 00268 G4cout << "G4ParameterisationPolyhedraRho::ComputeDimensions()" << G4endl 00269 << "-- Parametrised phedra copy-number: " << copyNo << G4endl; 00270 phedra.DumpInfo(); 00271 } 00272 #endif 00273 }
void G4ParameterisationPolyhedraRho::ComputeTransformation | ( | const G4int | copyNo, | |
G4VPhysicalVolume * | physVol | |||
) | const [virtual] |
Implements G4VDivisionParameterisation.
Definition at line 208 of file G4ParameterisationPolyhedra.cc.
References G4VDivisionParameterisation::ChangeRotMatrix(), G4VDivisionParameterisation::faxis, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4VPhysicalVolume::SetTranslation(), and G4VDivisionParameterisation::verbose.
00209 { 00210 //----- translation 00211 G4ThreeVector origin(0.,0.,0.); 00212 00213 //----- set translation 00214 physVol->SetTranslation( origin ); 00215 00216 //----- calculate rotation matrix: unit 00217 00218 #ifdef G4DIVDEBUG 00219 if( verbose >= 2 ) 00220 { 00221 G4cout << " G4ParameterisationPolyhedraRho " << G4endl 00222 << " foffset: " << foffset/deg 00223 << " - fwidth: " << fwidth/deg << G4endl; 00224 } 00225 #endif 00226 00227 ChangeRotMatrix( physVol ); 00228 00229 #ifdef G4DIVDEBUG 00230 if( verbose >= 2 ) 00231 { 00232 G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraRho " 00233 << G4endl 00234 << " Position: " << origin 00235 << " - Width: " << fwidth 00236 << " - Axis: " << faxis << G4endl; 00237 } 00238 #endif 00239 }
G4double G4ParameterisationPolyhedraRho::GetMaxParameter | ( | ) | const [virtual] |
Implements G4VDivisionParameterisation.
Definition at line 198 of file G4ParameterisationPolyhedra.cc.
References G4VDivisionParameterisation::fmotherSolid, G4Polyhedra::GetOriginalParameters(), G4PolyhedraHistorical::Rmax, and G4PolyhedraHistorical::Rmin.
00199 { 00200 G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid); 00201 G4PolyhedraHistorical* original_pars = msol->GetOriginalParameters(); 00202 return original_pars->Rmax[0] - original_pars->Rmin[0]; 00203 }