#include <G4ParameterisationTrd.hh>
Inheritance diagram for G4ParameterisationTrdX:
Public Member Functions | |
G4ParameterisationTrdX (EAxis axis, G4int nCopies, G4double width, G4double offset, G4VSolid *motherSolid, DivisionType divType) | |
~G4ParameterisationTrdX () | |
void | CheckParametersValidity () |
G4double | GetMaxParameter () const |
void | ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const |
void | ComputeDimensions (G4Trd &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const |
void | ComputeDimensions (G4Trap &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const |
G4VSolid * | ComputeSolid (const G4int, G4VPhysicalVolume *) |
Definition at line 81 of file G4ParameterisationTrd.hh.
G4ParameterisationTrdX::G4ParameterisationTrdX | ( | EAxis | axis, | |
G4int | nCopies, | |||
G4double | width, | |||
G4double | offset, | |||
G4VSolid * | motherSolid, | |||
DivisionType | divType | |||
) |
Definition at line 83 of file G4ParameterisationTrd.cc.
References G4VParameterisationTrd::bDivInTrap, G4VDivisionParameterisation::CalculateNDiv(), G4VDivisionParameterisation::CalculateWidth(), CheckParametersValidity(), DivNDIV, DivWIDTH, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fnDiv, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Trd::GetXHalfLength1(), G4Trd::GetXHalfLength2(), G4VDivisionParameterisation::kCarTolerance, G4VDivisionParameterisation::SetType(), and G4VDivisionParameterisation::verbose.
00086 : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType ) 00087 { 00088 CheckParametersValidity(); 00089 SetType( "DivisionTrdX" ); 00090 00091 G4Trd* msol = (G4Trd*)(fmotherSolid); 00092 if( divType == DivWIDTH ) 00093 { 00094 fnDiv = CalculateNDiv( msol->GetXHalfLength1()+msol->GetXHalfLength2(), 00095 width, offset ); 00096 } 00097 else if( divType == DivNDIV ) 00098 { 00099 fwidth = CalculateWidth( msol->GetXHalfLength1()+msol->GetXHalfLength2(), 00100 nDiv, offset ); 00101 } 00102 00103 #ifdef G4DIVDEBUG 00104 if( verbose >= 1 ) 00105 { 00106 G4cout << " G4ParameterisationTrdX - ## divisions " << fnDiv << " = " 00107 << nDiv << G4endl 00108 << " Offset " << foffset << " = " << offset << G4endl 00109 << " Width " << fwidth << " = " << width << G4endl; 00110 } 00111 #endif 00112 00113 G4double mpDx1 = msol->GetXHalfLength1(); 00114 G4double mpDx2 = msol->GetXHalfLength2(); 00115 if( std::fabs(mpDx1 - mpDx2) > kCarTolerance ) 00116 { 00117 bDivInTrap = true; 00118 } 00119 }
G4ParameterisationTrdX::~G4ParameterisationTrdX | ( | ) |
void G4ParameterisationTrdX::CheckParametersValidity | ( | ) | [virtual] |
Reimplemented from G4VDivisionParameterisation.
Definition at line 259 of file G4ParameterisationTrd.cc.
References G4VDivisionParameterisation::CheckParametersValidity().
Referenced by G4ParameterisationTrdX().
00260 { 00261 G4VDivisionParameterisation::CheckParametersValidity(); 00262 /* 00263 G4Trd* msol = (G4Trd*)(fmotherSolid); 00264 00265 G4double mpDx1 = msol->GetXHalfLength1(); 00266 G4double mpDx2 = msol->GetXHalfLength2(); 00267 bDivInTrap = false; 00268 00269 if( std::fabs(mpDx1 - mpDx2) > kCarTolerance ) 00270 { 00271 std::ostringstream message; 00272 message << "Invalid solid specification. NOT supported." << G4endl 00273 << "Making a division of a TRD along axis X," << G4endl 00274 << "while the X half lengths are not equal," << G4endl 00275 << "is not (yet) supported. It will result" << G4endl 00276 << "in non-equal division solids."; 00277 G4Exception("G4ParameterisationTrdX::CheckParametersValidity()", 00278 "GeomDiv0001", FatalException, message); 00279 } 00280 */ 00281 }
void G4ParameterisationTrdX::ComputeDimensions | ( | G4Trap & | trd, | |
const G4int | copyNo, | |||
const G4VPhysicalVolume * | pv | |||
) | const [virtual] |
Reimplemented from G4VPVParameterisation.
Definition at line 220 of file G4ParameterisationTrd.cc.
References G4VSolid::DumpInfo(), G4VDivisionParameterisation::fhgap, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fnDiv, G4VDivisionParameterisation::foffset, G4cout, G4endl, G4Trd::GetXHalfLength1(), G4Trd::GetXHalfLength2(), G4Trd::GetYHalfLength1(), G4Trd::GetYHalfLength2(), G4Trd::GetZHalfLength(), G4Trap::SetAllParameters(), and G4VDivisionParameterisation::verbose.
00222 { 00223 G4Trd* msol = (G4Trd*)(fmotherSolid); 00224 G4double pDy1 = msol->GetYHalfLength1(); 00225 G4double pDy2 = msol->GetYHalfLength2(); 00226 G4double pDz = msol->GetZHalfLength(); 00227 G4double pDx1 = msol->GetXHalfLength1()/fnDiv; 00228 G4double pDx2 = msol->GetXHalfLength2()/fnDiv; 00229 00230 G4double cxy1 = -msol->GetXHalfLength1() + foffset 00231 + (copyNo+0.5)*pDx1*2;// centre of the side at y=-pDy1 00232 G4double cxy2 = -msol->GetXHalfLength2() + foffset 00233 + (copyNo+0.5)*pDx2*2;// centre of the side at y=+pDy1 00234 G4double alp = std::atan( (cxy2-cxy1)/pDz ); 00235 00236 trap.SetAllParameters ( pDz, 00237 0., 00238 0., 00239 pDy1, 00240 pDx1, 00241 pDx2, 00242 alp, 00243 pDy2, 00244 pDx1 - fhgap, 00245 pDx2 - fhgap * pDx2/pDx1, 00246 alp); 00247 00248 #ifdef G4DIVDEBUG 00249 if( verbose >= 2 ) 00250 { 00251 G4cout << " G4ParameterisationTrdX::ComputeDimensions():" 00252 << copyNo << G4endl; 00253 trap.DumpInfo(); 00254 } 00255 #endif 00256 }
void G4ParameterisationTrdX::ComputeDimensions | ( | G4Trd & | trd, | |
const G4int | copyNo, | |||
const G4VPhysicalVolume * | pv | |||
) | const [virtual] |
Reimplemented from G4VPVParameterisation.
Definition at line 183 of file G4ParameterisationTrd.cc.
References G4VSolid::DumpInfo(), G4VDivisionParameterisation::fhgap, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Trd::GetYHalfLength1(), G4Trd::GetYHalfLength2(), G4Trd::GetZHalfLength(), G4Trd::SetAllParameters(), and G4VDivisionParameterisation::verbose.
00184 { 00185 G4Trd* msol = (G4Trd*)(fmotherSolid); 00186 00187 G4double pDy1 = msol->GetYHalfLength1(); 00188 G4double pDy2 = msol->GetYHalfLength2(); 00189 G4double pDz = msol->GetZHalfLength(); 00190 G4double pDx = fwidth/2. - fhgap; 00191 00192 trd.SetAllParameters ( pDx, pDx, pDy1, pDy2, pDz ); 00193 00194 #ifdef G4DIVDEBUG 00195 if( verbose >= 2 ) 00196 { 00197 G4cout << " G4ParameterisationTrdX::ComputeDimensions():" 00198 << copyNo << G4endl; 00199 trd.DumpInfo(); 00200 } 00201 #endif 00202 }
G4VSolid * G4ParameterisationTrdX::ComputeSolid | ( | const | G4int, | |
G4VPhysicalVolume * | ||||
) | [virtual] |
Reimplemented from G4VDivisionParameterisation.
Definition at line 205 of file G4ParameterisationTrd.cc.
References G4VParameterisationTrd::bDivInTrap, G4VDivisionParameterisation::ComputeSolid(), and G4VDivisionParameterisation::fmotherSolid.
00206 { 00207 if( bDivInTrap ) 00208 { 00209 return G4VDivisionParameterisation::ComputeSolid(i, pv); 00210 } 00211 else 00212 { 00213 return fmotherSolid; 00214 } 00215 }
void G4ParameterisationTrdX::ComputeTransformation | ( | const G4int | copyNo, | |
G4VPhysicalVolume * | physVol | |||
) | const [virtual] |
Implements G4VDivisionParameterisation.
Definition at line 136 of file G4ParameterisationTrd.cc.
References G4VParameterisationTrd::bDivInTrap, FatalException, G4VDivisionParameterisation::faxis, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fnDiv, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Exception(), G4Trd::GetXHalfLength1(), G4Trd::GetXHalfLength2(), kXAxis, G4VPhysicalVolume::SetTranslation(), and G4VDivisionParameterisation::verbose.
00138 { 00139 G4Trd* msol = (G4Trd*)(fmotherSolid ); 00140 G4double mdx = ( msol->GetXHalfLength1() + msol->GetXHalfLength2() ) / 2.; 00141 00142 //----- translation 00143 G4ThreeVector origin(0.,0.,0.); 00144 G4double posi; 00145 if( !bDivInTrap ) 00146 { 00147 posi = -mdx + foffset + (copyNo+0.5)*fwidth; 00148 } 00149 else 00150 { 00151 G4double aveHL = (msol->GetXHalfLength1()+msol->GetXHalfLength2())/2.; 00152 posi = - aveHL + foffset + (copyNo+0.5)*aveHL/fnDiv*2; 00153 } 00154 if( faxis == kXAxis ) 00155 { 00156 origin.setX( posi ); 00157 } 00158 else 00159 { 00160 std::ostringstream message; 00161 message << "Only axes along X are allowed ! Axis: " << faxis; 00162 G4Exception("G4ParameterisationTrdX::ComputeTransformation()", 00163 "GeomDiv0002", FatalException, message); 00164 } 00165 00166 #ifdef G4DIVDEBUG 00167 if( verbose >= 2 ) 00168 { 00169 G4cout << std::setprecision(8) 00170 << " G4ParameterisationTrdX::ComputeTransformation() " 00171 << copyNo << G4endl 00172 << " Position: " << origin << " - Axis: " << faxis << G4endl; 00173 } 00174 #endif 00175 00176 //----- set translation 00177 physVol->SetTranslation( origin ); 00178 }
G4double G4ParameterisationTrdX::GetMaxParameter | ( | ) | const [virtual] |
Implements G4VDivisionParameterisation.
Definition at line 127 of file G4ParameterisationTrd.cc.
References G4VDivisionParameterisation::fmotherSolid, G4Trd::GetXHalfLength1(), and G4Trd::GetXHalfLength2().
00128 { 00129 G4Trd* msol = (G4Trd*)(fmotherSolid); 00130 return (msol->GetXHalfLength1()+msol->GetXHalfLength2()); 00131 }