/scratch0/jsnuveri/BDSIM/BDSIMgit/bdsim/src/BDSBeamPipeFactoryElliptical.cc

00001 #include "BDSBeamPipeFactoryBase.hh"
00002 #include "BDSBeamPipeFactoryElliptical.hh"
00003 #include "BDSBeamPipe.hh"
00004 
00005 #include "BDSDebug.hh"
00006 #include "BDSExecOptions.hh"
00007 #include "BDSGlobalConstants.hh"
00008 #include "BDSSDManager.hh"
00009 
00010 #include "globals.hh"                 // geant4 globals / types
00011 #include "G4Box.hh"
00012 #include "G4EllipticalTube.hh"
00013 #include "G4CutTubs.hh"
00014 #include "G4IntersectionSolid.hh"
00015 #include "G4LogicalVolume.hh"
00016 #include "G4Material.hh"
00017 #include "G4SubtractionSolid.hh"
00018 #include "G4ThreeVector.hh"
00019 #include "G4Tubs.hh"
00020 #include "G4VSolid.hh"
00021 
00022 #include <cmath>                           // sin, cos, fabs
00023 #include <utility>                         // for std::pair
00024 
00025 
00026 BDSBeamPipeFactoryElliptical* BDSBeamPipeFactoryElliptical::_instance = 0;
00027 
00028 BDSBeamPipeFactoryElliptical* BDSBeamPipeFactoryElliptical::Instance()
00029 {
00030   if (_instance == 0)
00031     {_instance = new BDSBeamPipeFactoryElliptical();}
00032   return _instance;
00033 }
00034 
00035 BDSBeamPipeFactoryElliptical::BDSBeamPipeFactoryElliptical():BDSBeamPipeFactoryBase()
00036 {
00037 }
00038 
00039 BDSBeamPipeFactoryElliptical::~BDSBeamPipeFactoryElliptical()
00040 {
00041   _instance = 0;
00042 }
00043 
00044 BDSBeamPipe* BDSBeamPipeFactoryElliptical::CreateBeamPipe(G4String    nameIn,              // name
00045                                                           G4double    lengthIn,            // length [mm]
00046                                                           G4double    aper1In,             // aperture parameter 1
00047                                                           G4double    aper2In,             // aperture parameter 2
00048                                                           G4double    /*aper3In*/,         // aperture parameter 3
00049                                                           G4double    /*aper4In*/,         // aperture parameter 4
00050                                                           G4Material* vacuumMaterialIn,    // vacuum material
00051                                                           G4double    beamPipeThicknessIn, // beampipe thickness [mm]
00052                                                           G4Material* beamPipeMaterialIn   // beampipe material
00053                                                           )
00054 {
00055 #ifdef BDSDEBUG
00056   G4cout << __METHOD_NAME__ << G4endl;
00057 #endif
00058   // clean up after last usage
00059   CleanUp();
00060   
00061   // test input parameters - set global options as default if not specified
00062   TestInputParameters(vacuumMaterialIn,beamPipeThicknessIn,beamPipeMaterialIn,aper1In,aper2In);
00063 
00064   // build the solids
00065   vacuumSolid   = new G4EllipticalTube(nameIn + "_vacuum_solid",       // name
00066                                        aper1In,                        // x half width
00067                                        aper2In,                        // y half width
00068                                        (lengthIn*0.5)-2*lengthSafety); // half length
00069 
00070   G4VSolid* beamPipeSolidInner; // construct rectangular beam pipe by subtracting an inner
00071   G4VSolid* beamPipeSolidOuter; // box from an outer one - only way
00072   // beamPipeSolidInner will be the inner edge of the metal beampipe
00073   // therefore it has to be the width of the aperture + lengthSafety
00074   beamPipeSolidInner = new G4EllipticalTube(nameIn + "_pipe_solid_inner",   // name
00075                                             aper1In + lengthSafety,         // x half width - length safety to avoid overlaps
00076                                             aper2In + lengthSafety,         // y half width
00077                                             lengthIn);                      // length - full length fo unambiguous subtraction
00078   // beamPipeSolidOuter will be the outer edge of the metal beampipe
00079   // therefore it has to be the width of the aperture + beampipeThickness
00080   beamPipeSolidOuter = new G4EllipticalTube(nameIn + "_pipe_solid_outer",   // name
00081                                             aper1In + beamPipeThicknessIn,  // x half width
00082                                             aper2In + beamPipeThicknessIn,  // y half width
00083                                             (lengthIn*0.5)-2*lengthSafety); // half length - lengthSafety to fit in container
00084   beamPipeSolid = new G4SubtractionSolid(nameIn + "_pipe_solid",
00085                                          beamPipeSolidOuter,
00086                                          beamPipeSolidInner); // outer minus inner
00087   
00088   G4double containerXHalfWidth = aper1In + beamPipeThicknessIn + lengthSafety;
00089   G4double containerYHalfWidth = aper2In + beamPipeThicknessIn + lengthSafety;
00090   containerSolid = new G4EllipticalTube(nameIn  + "_container_solid",  // name
00091                                         containerXHalfWidth,           // x half width
00092                                         containerYHalfWidth,           // y half width
00093                                         (lengthIn*0.5)-lengthSafety);  // half length
00094                                         
00095   return CommonFinalConstruction(nameIn, vacuumMaterialIn, beamPipeMaterialIn, lengthIn, aper1In, aper2In, beamPipeThicknessIn);
00096 }
00097 
00098 BDSBeamPipe* BDSBeamPipeFactoryElliptical::CreateBeamPipeAngledInOut(G4String    nameIn,              // name
00099                                                                      G4double    lengthIn,            // length [mm]
00100                                                                      G4double    angleInIn,           // the normal angle of the input face
00101                                                                      G4double    angleOutIn,          // the normal angle of the input face
00102                                                                      G4double    aper1In,             // aperture parameter 1
00103                                                                      G4double    aper2In,             // aperture parameter 2
00104                                                                      G4double    /*aper3In*/,         // aperture parameter 3
00105                                                                      G4double    /*aper4In */,        // aperture parameter 4
00106                                                                      G4Material* vacuumMaterialIn,    // vacuum material
00107                                                                      G4double    beamPipeThicknessIn, // beampipe thickness [mm]
00108                                                                       G4Material* beamPipeMaterialIn   // beampipe material
00109                                                                      )
00110 {
00111 #ifdef BDSDEBUG
00112   G4cout << __METHOD_NAME__ << G4endl;
00113 #endif
00114   // clean up after last usage
00115   CleanUp();
00116   
00117   // test input parameters - set global options as default if not specified
00118   TestInputParameters(vacuumMaterialIn,beamPipeThicknessIn,beamPipeMaterialIn,aper1In,aper2In);
00119 
00120   std::pair<G4ThreeVector,G4ThreeVector> faces = CalculateFaces(angleInIn, angleOutIn);
00121   G4ThreeVector inputface  = faces.first;
00122   G4ThreeVector outputface = faces.second;
00123   
00124   CreateGeneralAngledSolids(nameIn, lengthIn, aper1In, aper2In, beamPipeThicknessIn, inputface, outputface);
00125   
00126   return CommonFinalConstruction(nameIn, vacuumMaterialIn, beamPipeMaterialIn, lengthIn, aper1In, aper2In, beamPipeThicknessIn);
00127 }
00128 
00130 
00132 void BDSBeamPipeFactoryElliptical::TestInputParameters(G4Material*&  vacuumMaterialIn,     // reference to a pointer
00133                                                        G4double&     beamPipeThicknessIn,
00134                                                        G4Material*&  beamPipeMaterialIn,
00135                                                        G4double&     aper1In,
00136                                                        G4double&     aper2In)
00137 {
00138   BDSBeamPipeFactoryBase::TestInputParameters(vacuumMaterialIn,beamPipeThicknessIn,beamPipeMaterialIn);
00139 
00140   if (aper1In < 1e-10)
00141     {aper1In = BDSGlobalConstants::Instance()->GetBeamPipeRadius();}
00142 
00143   if (aper2In < 1e-10)
00144     {aper2In = BDSGlobalConstants::Instance()->GetBeamPipeRadius();}
00145 }
00146 
00149 BDSBeamPipe* BDSBeamPipeFactoryElliptical::CommonFinalConstruction(G4String    nameIn,
00150                                                                    G4Material* vacuumMaterialIn,
00151                                                                    G4Material* beamPipeMaterialIn,
00152                                                                    G4double    lengthIn,
00153                                                                    G4double    aper1In,
00154                                                                    G4double    aper2In,
00155                                                                    G4double    beamPipeThicknessIn)
00156 {
00157 #ifdef BDSDEBUG
00158   G4cout << __METHOD_NAME__ << G4endl;
00159 #endif
00160   // prepare a longer container subtraction solid
00161   G4double containerXHalfWidth = aper1In + beamPipeThicknessIn + lengthSafety;
00162   G4double containerYHalfWidth = aper2In + beamPipeThicknessIn + lengthSafety;
00163   // doesn't have to be angled as it's only used for transverse subtraction
00164   containerSubtractionSolid = new G4EllipticalTube(nameIn  + "_container_solid", // name
00165                                                    containerXHalfWidth,          // x half width
00166                                                    containerYHalfWidth,          // y half width
00167                                                    lengthIn);                    // full length for unambiguous subtraction
00168 
00169 
00170   BDSBeamPipeFactoryBase::CommonConstruction(nameIn,
00171                                              vacuumMaterialIn,
00172                                              beamPipeMaterialIn,
00173                                              lengthIn);
00174 
00175   // record extents
00176   std::pair<double,double> extX = std::make_pair(-containerXHalfWidth,containerXHalfWidth);
00177   std::pair<double,double> extY = std::make_pair(-containerYHalfWidth,containerYHalfWidth);
00178   std::pair<double,double> extZ = std::make_pair(-lengthIn*0.5,lengthIn*0.5);
00179   // calculate radius if a tube were to be place around it
00180   G4double containerRadius = std::max(containerXHalfWidth, containerYHalfWidth);
00181   
00182   BDSBeamPipe* aPipe = BuildBeamPipeAndRegisterVolumes(extX,extY,extZ,containerRadius);
00183 
00184   // REGISTER all lvs
00185   aPipe->RegisterLogicalVolume(vacuumLV); //using geometry component base class method
00186   aPipe->RegisterLogicalVolume(beamPipeLV);
00187 
00188   // register sensitive volumes
00189   aPipe->RegisterSensitiveVolume(beamPipeLV);
00190   
00191   return aPipe;
00192 }
00193 
00194 
00197 void BDSBeamPipeFactoryElliptical::CreateGeneralAngledSolids(G4String      nameIn,
00198                                                              G4double      lengthIn,
00199                                                              G4double      aper1In,
00200                                                              G4double      aper2In,
00201                                                              G4double      beamPipeThicknessIn,
00202                                                              G4ThreeVector inputfaceIn,
00203                                                              G4ThreeVector outputfaceIn)
00204 {
00205   // this function will make a longer normal rectangular beampipe and chop it off
00206   // to make angled faces as required
00207   // achieve this using the intersection of the normal beampipe (but a little longer)
00208   // with a large G4CutTubs to get the angled faces.
00209   // note even if one face is flat, we don't save a boolean operation as the intersection
00210   // can be on both sides using a G4CutTubs.  Also, keeping one side flat would require
00211   // shifting the volume from 0 which causes headaches later with SDs.
00212 
00213   // build the solids - vacuum, beampipe and container solids
00214   // extra solids required for booleans
00215   G4VSolid* vacuumSolidLong;
00216   G4VSolid* beamPipeSolidLong;
00217   G4VSolid* angledFaceSolid;
00218   G4VSolid* containerSolidLong;
00219   G4VSolid* angledFaceSolidContainer;
00220 
00221   // build the solid with angled faces for intersection
00222   G4double angledFaceRadius = (std::max(aper1In,aper2In) + beamPipeThicknessIn)*2.0; //huge for unambiguous intersection
00223   angledFaceSolid = new G4CutTubs(nameIn + "_angled_face",       // name
00224                                   0,                             // inner radius
00225                                   angledFaceRadius,              // outer radius
00226                                   (lengthIn*0.5)-2*lengthSafety, // half length - must fit within container
00227                                   0,                             // rotation start angle
00228                                   CLHEP::twopi,                  // rotation finish angle
00229                                   inputfaceIn,                   // input face normal
00230                                   outputfaceIn);                 // output face normal
00231   
00232   vacuumSolidLong = new G4EllipticalTube(nameIn + "_vacuum_solid_long", // name
00233                                          aper1In,                       // x half width
00234                                          aper2In,                       // y half width
00235                                          lengthIn);                     // full length for unambiguous boolean
00236   vacuumSolid     = new G4IntersectionSolid(nameIn + "_vacuum_solid",
00237                                             vacuumSolidLong,
00238                                             angledFaceSolid);
00239   
00240   G4VSolid* beamPipeSolidInner; // construct rectangular beam pipe by subtracting an inner
00241   G4VSolid* beamPipeSolidOuter; // box from an outer one - only way
00242   // beamPipeSolidInner will be the inner edge of the metal beampipe
00243   // therefore it has to be the width of the aperture + lengthSafety
00244   beamPipeSolidInner = new G4EllipticalTube(nameIn + "_pipe_solid_inner",   // name
00245                                             aper1In + lengthSafety,         // x half width - length safety to avoid overlaps
00246                                             aper2In + lengthSafety,         // y half width
00247                                             2*lengthIn);                    // 2*length - full length fo unambiguous subtraction
00248   // beamPipeSolidOuter will be the outer edge of the metal beampipe
00249   // therefore it has to be the width of the aperture + beampipeThickness
00250   beamPipeSolidOuter = new G4EllipticalTube(nameIn + "_pipe_solid_outer",   // name
00251                                             aper1In + beamPipeThicknessIn,  // x half width
00252                                             aper2In + beamPipeThicknessIn,  // y half width
00253                                             lengthIn);                      // full length for unambiguous intersection
00254   beamPipeSolidLong = new G4SubtractionSolid(nameIn + "_pipe_solid_long",
00255                                          beamPipeSolidOuter,
00256                                          beamPipeSolidInner); // outer minus inner
00257   beamPipeSolid = new G4IntersectionSolid(nameIn + "_pipe_solid",
00258                                           beamPipeSolidLong,
00259                                           angledFaceSolid);
00260   
00261   G4double containerXHalfWidth = aper1In + beamPipeThicknessIn + lengthSafety;
00262   G4double containerYHalfWidth = aper2In + beamPipeThicknessIn + lengthSafety;
00263   containerSolidLong = new G4EllipticalTube(nameIn  + "_container_solid_long",// name
00264                                             containerXHalfWidth,              // x half width
00265                                             containerYHalfWidth,              // y half width
00266                                             lengthIn);                        // full length for unambiguous intersection
00267   angledFaceSolidContainer = new G4CutTubs(nameIn + "_angled_face_container",// name
00268                                            0,                                // inner radius
00269                                            angledFaceRadius,                 // outer radius
00270                                            (lengthIn*0.5)-lengthSafety,      // half length - must fit within magnet
00271                                            0,                                // rotation start angle
00272                                            CLHEP::twopi,                     // rotation finish angle
00273                                            inputfaceIn,                      // input face normal
00274                                            outputfaceIn);                    // output face normal
00275   containerSolid = new G4IntersectionSolid(nameIn + "_container_solid",
00276                                            containerSolidLong,
00277                                            angledFaceSolidContainer);
00278 }

Generated on 28 Jun 2015 for BDSIM by  doxygen 1.4.7