00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036 #include "G4ParameterisationBox.hh"
00037
00038 #include <iomanip>
00039 #include "G4ThreeVector.hh"
00040 #include "G4Transform3D.hh"
00041 #include "G4RotationMatrix.hh"
00042 #include "G4VPhysicalVolume.hh"
00043 #include "G4ReflectedSolid.hh"
00044 #include "G4Box.hh"
00045
00046
00047 G4VParameterisationBox::
00048 G4VParameterisationBox( EAxis axis, G4int nDiv, G4double width,
00049 G4double offset, G4VSolid* msolid,
00050 DivisionType divType )
00051 : G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid )
00052 {
00053 G4Box* msol = (G4Box*)(msolid);
00054 if (msolid->GetEntityType() == "G4ReflectedSolid")
00055 {
00056
00057 G4VSolid* mConstituentSolid
00058 = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid();
00059 msol = (G4Box*)(mConstituentSolid);
00060 fmotherSolid = msol;
00061 fReflectedSolid = true;
00062 }
00063 }
00064
00065
00066 G4VParameterisationBox::~G4VParameterisationBox()
00067 {
00068 }
00069
00070
00071 G4ParameterisationBoxX::
00072 G4ParameterisationBoxX( EAxis axis, G4int nDiv, G4double width,
00073 G4double offset, G4VSolid* msolid,
00074 DivisionType divType )
00075 : G4VParameterisationBox( axis, nDiv, width, offset, msolid, divType )
00076 {
00077 CheckParametersValidity();
00078 SetType( "DivisionBoxX" );
00079
00080 G4Box* mbox = (G4Box*)(fmotherSolid);
00081 if( divType == DivWIDTH )
00082 {
00083 fnDiv = CalculateNDiv( 2*mbox->GetXHalfLength(), width, offset );
00084 }
00085 else if( divType == DivNDIV )
00086 {
00087 fwidth = CalculateWidth( 2*mbox->GetXHalfLength(), nDiv, offset );
00088 }
00089 #ifdef G4DIVDEBUG
00090 if( verbose >= 1 )
00091 {
00092 G4cout << " G4ParameterisationBoxX - no divisions "
00093 << fnDiv << " = " << nDiv << G4endl
00094 << " Offset " << foffset << " = " << offset << G4endl
00095 << " Width " << fwidth << " = " << width << G4endl;
00096 }
00097 #endif
00098 }
00099
00100
00101 G4ParameterisationBoxX::~G4ParameterisationBoxX()
00102 {
00103 }
00104
00105
00106 G4double G4ParameterisationBoxX::GetMaxParameter() const
00107 {
00108 G4Box* msol = (G4Box*)(fmotherSolid);
00109 return 2*msol->GetXHalfLength();
00110 }
00111
00112
00113 void
00114 G4ParameterisationBoxX::
00115 ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol ) const
00116 {
00117 G4Box* msol = (G4Box*)(fmotherSolid );
00118 G4double mdx = msol->GetXHalfLength( );
00119
00120
00121 G4ThreeVector origin(0.,0.,0.);
00122 G4double posi = -mdx + foffset+(copyNo+0.5)*fwidth;
00123
00124 if( faxis == kXAxis )
00125 {
00126 origin.setX( posi );
00127 }
00128 else
00129 {
00130 std::ostringstream message;
00131 message << "Only axes along X are allowed ! Axis: " << faxis;
00132 G4Exception("G4ParameterisationBoxX::ComputeTransformation()",
00133 "GeomDiv0002", FatalException, message);
00134 }
00135 #ifdef G4DIVDEBUG
00136 if( verbose >= 2 )
00137 {
00138 G4cout << std::setprecision(8) << " G4ParameterisationBoxX: "
00139 << copyNo << G4endl
00140 << " Position " << origin << " Axis " << faxis << G4endl;
00141 }
00142 #endif
00143
00144 physVol->SetTranslation( origin );
00145 }
00146
00147
00148 void
00149 G4ParameterisationBoxX::
00150 ComputeDimensions( G4Box& box, const G4int,
00151 const G4VPhysicalVolume* ) const
00152 {
00153 G4Box* msol = (G4Box*)(fmotherSolid);
00154
00155 G4double pDx = fwidth/2. - fhgap;
00156 G4double pDy = msol->GetYHalfLength();
00157 G4double pDz = msol->GetZHalfLength();
00158
00159 box.SetXHalfLength( pDx );
00160 box.SetYHalfLength( pDy );
00161 box.SetZHalfLength( pDz );
00162
00163 #ifdef G4DIVDEBUG
00164 if( verbose >= 2 )
00165 {
00166 G4cout << " G4ParameterisationBoxX::ComputeDimensions()" << G4endl
00167 << " pDx: " << pDz << G4endl;
00168 box.DumpInfo();
00169 }
00170 #endif
00171 }
00172
00173
00174 G4ParameterisationBoxY::
00175 G4ParameterisationBoxY( EAxis axis, G4int nDiv, G4double width,
00176 G4double offset, G4VSolid* msolid,
00177 DivisionType divType)
00178 : G4VParameterisationBox( axis, nDiv, width, offset, msolid, divType )
00179 {
00180 CheckParametersValidity();
00181 SetType( "DivisionBoxY" );
00182
00183 G4Box* mbox = (G4Box*)(fmotherSolid);
00184 if( divType == DivWIDTH )
00185 {
00186 fnDiv = CalculateNDiv( 2*mbox->GetYHalfLength(), width, offset );
00187 }
00188 else if( divType == DivNDIV )
00189 {
00190 fwidth = CalculateWidth( 2*mbox->GetYHalfLength(), nDiv, offset );
00191 }
00192
00193 #ifdef G4DIVDEBUG
00194 if( verbose >= 1 )
00195 {
00196 G4cout << " G4ParameterisationBoxY - no divisions " << fnDiv << " = "
00197 << nDiv << ". Offset " << foffset << " = " << offset
00198 << ". Width " << fwidth << " = " << width << G4endl;
00199 }
00200 #endif
00201 }
00202
00203
00204 G4ParameterisationBoxY::~G4ParameterisationBoxY()
00205 {
00206 }
00207
00208
00209 G4double G4ParameterisationBoxY::GetMaxParameter() const
00210 {
00211 G4Box* msol = (G4Box*)(fmotherSolid);
00212 return 2*msol->GetYHalfLength();
00213 }
00214
00215
00216 void
00217 G4ParameterisationBoxY::
00218 ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol ) const
00219 {
00220 G4Box* msol = (G4Box*)(fmotherSolid);
00221 G4double mdy = msol->GetYHalfLength();
00222
00223
00224 G4ThreeVector origin(0.,0.,0.);
00225 G4double posi = -mdy + foffset + (copyNo+0.5)*fwidth;
00226 if( faxis == kYAxis )
00227 {
00228 origin.setY( posi );
00229 }
00230 else
00231 {
00232 std::ostringstream message;
00233 message << "Only axes along Y are allowed ! Axis: " << faxis;
00234 G4Exception("G4ParameterisationBoxY::ComputeTransformation()",
00235 "GeomDiv0002", FatalException, message);
00236 }
00237 #ifdef G4DIVDEBUG
00238 if( verbose >= 2 )
00239 {
00240 G4cout << std::setprecision(8) << " G4ParameterisationBoxY: "
00241 << copyNo << G4endl
00242 << " Position " << origin << " Axis " << faxis << G4endl;
00243 }
00244 #endif
00245
00246 physVol->SetTranslation( origin );
00247 }
00248
00249
00250 void
00251 G4ParameterisationBoxY::
00252 ComputeDimensions( G4Box& box, const G4int,
00253 const G4VPhysicalVolume* ) const
00254 {
00255 G4Box* msol = (G4Box*)(fmotherSolid);
00256
00257 G4double pDx = msol->GetXHalfLength();
00258 G4double pDy = fwidth/2. - fhgap;
00259 G4double pDz = msol->GetZHalfLength();
00260
00261 box.SetXHalfLength( pDx );
00262 box.SetYHalfLength( pDy );
00263 box.SetZHalfLength( pDz );
00264
00265 #ifdef G4DIVDEBUG
00266 if( verbose >= 2 )
00267 {
00268 G4cout << " G4ParameterisationBoxY::ComputeDimensions()" << G4endl
00269 << " pDx: " << pDz << G4endl;
00270 box.DumpInfo();
00271 }
00272 #endif
00273 }
00274
00275
00276 G4ParameterisationBoxZ::
00277 G4ParameterisationBoxZ( EAxis axis, G4int nDiv, G4double width,
00278 G4double offset, G4VSolid* msolid,
00279 DivisionType divType )
00280 : G4VParameterisationBox( axis, nDiv, width, offset, msolid, divType )
00281 {
00282 CheckParametersValidity();
00283 SetType( "DivisionBoxZ" );
00284
00285 G4Box* mbox = (G4Box*)(fmotherSolid);
00286 if( divType == DivWIDTH )
00287 {
00288 fnDiv = CalculateNDiv( 2*mbox->GetZHalfLength(), width, offset );
00289 }
00290 else if ( divType == DivNDIV )
00291 {
00292 fwidth = CalculateWidth( 2*mbox->GetZHalfLength(), nDiv, offset );
00293 }
00294 #ifdef G4DIVDEBUG
00295 if( verbose >= 1 )
00296 {
00297 G4cout << " G4ParameterisationBoxZ - no divisions " << fnDiv << " = "
00298 << nDiv << ". Offset " << foffset << " = " << offset
00299 << ". Width " << fwidth << " = " << width << G4endl;
00300 }
00301 #endif
00302 }
00303
00304
00305 G4ParameterisationBoxZ::~G4ParameterisationBoxZ()
00306 {
00307 }
00308
00309
00310 G4double G4ParameterisationBoxZ::GetMaxParameter() const
00311 {
00312 G4Box* msol = (G4Box*)(fmotherSolid);
00313 return 2*msol->GetZHalfLength();
00314 }
00315
00316
00317 void
00318 G4ParameterisationBoxZ::
00319 ComputeTransformation( const G4int copyNo, G4VPhysicalVolume *physVol ) const
00320 {
00321 G4Box* msol = (G4Box*)(fmotherSolid );
00322 G4double mdz = msol->GetZHalfLength();
00323
00324
00325 G4ThreeVector origin(0.,0.,0.);
00326 G4double posi = -mdz + OffsetZ() + (copyNo+0.5)*fwidth;
00327
00328 if( faxis == kZAxis )
00329 {
00330 origin.setZ( posi );
00331 }
00332 else
00333 {
00334 std::ostringstream message;
00335 message << "Only axes along Z are allowed ! Axis: " << faxis;
00336 G4Exception("G4ParameterisationBoxZ::ComputeTransformation()",
00337 "GeomDiv0002", FatalException, message);
00338 }
00339 #ifdef G4DIVDEBUG
00340 if( verbose >= 2 )
00341 {
00342 G4cout << std::setprecision(8) << " G4ParameterisationBoxZ: "
00343 << copyNo << G4endl
00344 << " Position " << origin << " Axis " << faxis << G4endl;
00345 }
00346 #endif
00347
00348 physVol->SetTranslation( origin );
00349 }
00350
00351
00352 void
00353 G4ParameterisationBoxZ::
00354 ComputeDimensions( G4Box& box, const G4int,
00355 const G4VPhysicalVolume* ) const
00356 {
00357 G4Box* msol = (G4Box*)(fmotherSolid);
00358
00359 G4double pDx = msol->GetXHalfLength();
00360 G4double pDy = msol->GetYHalfLength();
00361 G4double pDz = fwidth/2. - fhgap;
00362
00363 box.SetXHalfLength( pDx );
00364 box.SetYHalfLength( pDy );
00365 box.SetZHalfLength( pDz );
00366
00367 #ifdef G4DIVDEBUG
00368 if( verbose >= 2 )
00369 {
00370 G4cout << " G4ParameterisationBoxZ::ComputeDimensions()" << G4endl
00371 << " pDx: " << pDz << G4endl;
00372 box.DumpInfo();
00373 }
00374 #endif
00375 }
00376