G4ParameterisationTrd.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id: G4ParameterisationTrd.cc 69784 2013-05-15 09:16:06Z gcosmo $
00028 //
00029 // class G4ParameterisationTrd Implementation file
00030 //
00031 // 26.05.03 - P.Arce, Initial version
00032 // 08.04.04 - I.Hrivnacova, Implemented reflection
00033 // 21.04.10 - M.Asai, Added gaps
00034 // --------------------------------------------------------------------
00035 
00036 #include "G4ParameterisationTrd.hh"
00037 
00038 #include <iomanip>
00039 #include "G4ThreeVector.hh"
00040 #include "G4RotationMatrix.hh"
00041 #include "G4VPhysicalVolume.hh"
00042 #include "G4LogicalVolume.hh"
00043 #include "G4ReflectedSolid.hh"
00044 #include "G4Trd.hh"
00045 #include "G4Trap.hh"
00046 
00047 //--------------------------------------------------------------------------
00048 G4VParameterisationTrd::
00049 G4VParameterisationTrd( EAxis axis, G4int nDiv, G4double width,
00050                         G4double offset, G4VSolid* msolid,
00051                         DivisionType divType )
00052   :  G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid ),
00053      bDivInTrap(false)
00054 {
00055   G4Trd* msol = (G4Trd*)(msolid);
00056   if (msolid->GetEntityType() == "G4ReflectedSolid")
00057   {
00058     // Get constituent solid  
00059     G4VSolid* mConstituentSolid 
00060        = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid();
00061     msol = (G4Trd*)(mConstituentSolid);
00062   
00063     // Create a new solid with inversed parameters
00064     G4Trd* newSolid
00065       = new G4Trd(msol->GetName(),
00066                   msol->GetXHalfLength2(), msol->GetXHalfLength1(),
00067                   msol->GetYHalfLength2(), msol->GetYHalfLength1(),
00068                   msol->GetZHalfLength());
00069     msol = newSolid;
00070     fmotherSolid = newSolid;
00071     fReflectedSolid = true;
00072     fDeleteSolid = true;
00073   }    
00074 }
00075 
00076 //------------------------------------------------------------------------
00077 G4VParameterisationTrd::~G4VParameterisationTrd()
00078 {
00079 }
00080 
00081 //------------------------------------------------------------------------
00082 G4ParameterisationTrdX::
00083 G4ParameterisationTrdX( EAxis axis, G4int nDiv,
00084                         G4double width, G4double offset,
00085                         G4VSolid* msolid, DivisionType divType )
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 }
00120 
00121 //------------------------------------------------------------------------
00122 G4ParameterisationTrdX::~G4ParameterisationTrdX()
00123 {
00124 }
00125 
00126 //------------------------------------------------------------------------
00127 G4double G4ParameterisationTrdX::GetMaxParameter() const
00128 {
00129   G4Trd* msol = (G4Trd*)(fmotherSolid);
00130   return (msol->GetXHalfLength1()+msol->GetXHalfLength2());
00131 }
00132 
00133 //------------------------------------------------------------------------
00134 void
00135 G4ParameterisationTrdX::
00136 ComputeTransformation( const G4int copyNo,
00137                        G4VPhysicalVolume *physVol ) const
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 }
00179 
00180 //--------------------------------------------------------------------------
00181 void
00182 G4ParameterisationTrdX::
00183 ComputeDimensions( G4Trd& trd, const G4int, const G4VPhysicalVolume* ) const
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 }
00203 
00204 G4VSolid* G4ParameterisationTrdX::
00205 ComputeSolid(const G4int i, G4VPhysicalVolume * pv)
00206 {
00207   if( bDivInTrap ) 
00208   {
00209     return G4VDivisionParameterisation::ComputeSolid(i, pv);
00210   } 
00211   else 
00212   {
00213     return fmotherSolid;
00214   }
00215 }
00216 
00217 
00218 //--------------------------------------------------------------------------
00219 void
00220 G4ParameterisationTrdX::ComputeDimensions( G4Trap& trap, const G4int copyNo,
00221                                            const G4VPhysicalVolume* ) const
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 }
00257 
00258 //--------------------------------------------------------------------------
00259 void G4ParameterisationTrdX::CheckParametersValidity()
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 }
00282 
00283 //--------------------------------------------------------------------------
00284 void G4ParameterisationTrdX::
00285 ComputeTrapParams()
00286 {
00287 }
00288 
00289 //--------------------------------------------------------------------------
00290 G4ParameterisationTrdY::
00291 G4ParameterisationTrdY( EAxis axis, G4int nDiv,
00292                         G4double width, G4double offset,
00293                         G4VSolid* msolid, DivisionType divType)
00294   : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
00295 {
00296   CheckParametersValidity();
00297   SetType( "DivisionTrdY" );
00298 
00299   G4Trd* msol = (G4Trd*)(fmotherSolid);
00300   if( divType == DivWIDTH )
00301   {
00302     fnDiv = CalculateNDiv( 2*msol->GetYHalfLength1(),
00303                            width, offset );
00304   }
00305   else if( divType == DivNDIV )
00306   {
00307     fwidth = CalculateWidth( 2*msol->GetYHalfLength1(),
00308                              nDiv, offset );
00309   }
00310 
00311 #ifdef G4DIVDEBUG
00312   if( verbose >= 1 )
00313   {
00314      G4cout << " G4ParameterisationTrdY no divisions " << fnDiv
00315             << " = " << nDiv << G4endl
00316             << " Offset " << foffset << " = " << offset << G4endl
00317             << " width " << fwidth << " = " << width << G4endl;
00318   }
00319 #endif
00320 }
00321 
00322 //------------------------------------------------------------------------
00323 G4ParameterisationTrdY::~G4ParameterisationTrdY()
00324 {
00325 }
00326 
00327 //------------------------------------------------------------------------
00328 G4double G4ParameterisationTrdY::GetMaxParameter() const
00329 {
00330   G4Trd* msol = (G4Trd*)(fmotherSolid);
00331   return 2*msol->GetYHalfLength1();
00332 }
00333 
00334 //--------------------------------------------------------------------------
00335 void
00336 G4ParameterisationTrdY::
00337 ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol ) const
00338 {
00339   G4Trd* msol = (G4Trd*)(fmotherSolid );
00340   G4double mdy = msol->GetYHalfLength1();
00341 
00342   //----- translation 
00343   G4ThreeVector origin(0.,0.,0.); 
00344   G4double posi = -mdy + foffset + (copyNo+0.5)*fwidth;
00345 
00346   if( faxis == kYAxis )
00347   {
00348     origin.setY( posi ); 
00349   }
00350   else
00351   { 
00352     std::ostringstream message;
00353     message << "Only axes along Y are allowed !  Axis: " << faxis;
00354     G4Exception("G4ParameterisationTrdY::ComputeTransformation()",
00355                 "GeomDiv0002", FatalException, message);
00356   }
00357 
00358 #ifdef G4DIVDEBUG
00359   if( verbose >= 2 )
00360   {
00361     G4cout << std::setprecision(8)
00362            << " G4ParameterisationTrdY::ComputeTransformation " << copyNo
00363            << " pos " << origin << " rot mat " << " axis " << faxis << G4endl;
00364   }
00365 #endif
00366 
00367    //----- set translation 
00368   physVol->SetTranslation( origin );
00369 }
00370 
00371 //--------------------------------------------------------------------------
00372 void
00373 G4ParameterisationTrdY::
00374 ComputeDimensions(G4Trd& trd, const G4int, const G4VPhysicalVolume*) const
00375 {
00376   //---- The division along Y of a Trd will result a Trd, only 
00377   //--- if Y at -Z and +Z are equal, else use the G4Trap version
00378   G4Trd* msol = (G4Trd*)(fmotherSolid);
00379   
00380   G4double pDx1 = msol->GetXHalfLength1();
00381   G4double pDx2 = msol->GetXHalfLength2();
00382   G4double pDz = msol->GetZHalfLength();
00383   G4double pDy = fwidth/2. - fhgap;
00384  
00385   trd.SetAllParameters ( pDx1, pDx2, pDy, pDy, pDz );
00386 
00387 #ifdef G4DIVDEBUG
00388   if( verbose >= 2 )
00389   {
00390     G4cout << " G4ParameterisationTrdY::ComputeDimensions():" << G4endl;
00391     trd.DumpInfo();
00392   }
00393 #endif
00394 }
00395 
00396 //--------------------------------------------------------------------------
00397 void G4ParameterisationTrdY::CheckParametersValidity()
00398 {
00399   G4VDivisionParameterisation::CheckParametersValidity();
00400 
00401   G4Trd* msol = (G4Trd*)(fmotherSolid);
00402 
00403   G4double mpDy1 = msol->GetYHalfLength1();
00404   G4double mpDy2 = msol->GetYHalfLength2();
00405 
00406   if( std::fabs(mpDy1 - mpDy2) > kCarTolerance )
00407   {
00408     std::ostringstream message;
00409     message << "Invalid solid specification. NOT supported." << G4endl
00410             << "Making a division of a TRD along axis Y while" << G4endl
00411             << "the Y half lengths are not equal is not (yet)" << G4endl
00412             << "supported. It will result in non-equal" << G4endl
00413             << "division solids.";
00414     G4Exception("G4ParameterisationTrdY::CheckParametersValidity()",
00415                 "GeomDiv0001", FatalException, message);
00416   }
00417 }
00418 
00419 //--------------------------------------------------------------------------
00420 G4ParameterisationTrdZ::
00421 G4ParameterisationTrdZ( EAxis axis, G4int nDiv,
00422                         G4double width, G4double offset,
00423                         G4VSolid* msolid, DivisionType divType )
00424   : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
00425 { 
00426   CheckParametersValidity();
00427   SetType( "DivTrdZ" );
00428 
00429   G4Trd* msol = (G4Trd*)(fmotherSolid);
00430   if( divType == DivWIDTH )
00431   {
00432     fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(),
00433                            width, offset );
00434   }
00435   else if( divType == DivNDIV )
00436   {
00437     fwidth = CalculateWidth( 2*msol->GetZHalfLength(),
00438                              nDiv, offset );
00439   }
00440 
00441 #ifdef G4DIVDEBUG
00442   if( verbose >= 1 )
00443   {
00444     G4cout << " G4ParameterisationTrdZ no divisions " << fnDiv
00445            << " = " << nDiv << G4endl
00446            << " Offset " << foffset << " = " << offset << G4endl
00447            << " Width " << fwidth << " = " << width << G4endl;
00448   }
00449 #endif
00450 }
00451 
00452 //------------------------------------------------------------------------
00453 G4ParameterisationTrdZ::~G4ParameterisationTrdZ()
00454 {
00455 }
00456 
00457 //------------------------------------------------------------------------
00458 G4double G4ParameterisationTrdZ::GetMaxParameter() const
00459 {
00460   G4Trd* msol = (G4Trd*)(fmotherSolid);
00461   return 2*msol->GetZHalfLength();
00462 }
00463 
00464 //--------------------------------------------------------------------------
00465 void
00466 G4ParameterisationTrdZ::
00467 ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
00468 {
00469   G4Trd* msol = (G4Trd*)(fmotherSolid );
00470   G4double mdz = msol->GetZHalfLength();
00471 
00472   //----- translation 
00473   G4ThreeVector origin(0.,0.,0.); 
00474   G4double posi = -mdz + OffsetZ() + (copyNo+0.5)*fwidth;
00475   if( faxis == kZAxis )
00476   {
00477     origin.setZ( posi ); 
00478   }
00479   else
00480   { 
00481     std::ostringstream message;
00482     message << "Only axes along Z are allowed !  Axis: " << faxis;
00483     G4Exception("G4ParameterisationTrdZ::ComputeTransformation()",
00484                 "GeomDiv0002", FatalException, message);
00485   }
00486 
00487 #ifdef G4DIVDEBUG
00488   if( verbose >= 1 )
00489   {
00490     G4cout << std::setprecision(8) << " G4ParameterisationTrdZ: "
00491            << copyNo << G4endl
00492            << " Position: " << origin << " - Offset: " << foffset 
00493            << " - Width: " << fwidth << " Axis " << faxis << G4endl;
00494   }
00495 #endif
00496 
00497    //----- set translation 
00498   physVol->SetTranslation( origin );
00499 }
00500 
00501 //--------------------------------------------------------------------------
00502 void
00503 G4ParameterisationTrdZ::
00504 ComputeDimensions(G4Trd& trd, const G4int copyNo,
00505                   const G4VPhysicalVolume*) const
00506 {
00507   //---- The division along Z of a Trd will result a Trd
00508   G4Trd* msol = (G4Trd*)(fmotherSolid);
00509 
00510   G4double pDx1 = msol->GetXHalfLength1();
00511   G4double DDx = (msol->GetXHalfLength2() - msol->GetXHalfLength1() );
00512   G4double pDy1 = msol->GetYHalfLength1();
00513   G4double DDy = (msol->GetYHalfLength2() - msol->GetYHalfLength1() );
00514   G4double pDz = fwidth/2. - fhgap;
00515   G4double zLength = 2*msol->GetZHalfLength();
00516  
00517   trd.SetAllParameters( pDx1+DDx*(OffsetZ()+copyNo*fwidth+fhgap)/zLength,
00518                         pDx1+DDx*(OffsetZ()+(copyNo+1)*fwidth-fhgap)/zLength, 
00519                         pDy1+DDy*(OffsetZ()+copyNo*fwidth+fhgap)/zLength,
00520                         pDy1+DDy*(OffsetZ()+(copyNo+1)*fwidth-fhgap)/zLength,
00521                         pDz );
00522 
00523 #ifdef G4DIVDEBUG
00524   if( verbose >= 1 )
00525   {
00526     G4cout << " G4ParameterisationTrdZ::ComputeDimensions()"
00527            << " - Mother TRD " << G4endl;
00528     msol->DumpInfo();
00529     G4cout << " - Parameterised TRD: "
00530            << copyNo << G4endl;
00531     trd.DumpInfo();
00532   }
00533 #endif
00534 }

Generated on Mon May 27 17:49:15 2013 for Geant4 by  doxygen 1.4.7