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

00001 #include "BDSMagnetOuterFactoryPolesFacetCrop.hh"
00002 
00003 #include "globals.hh"  // geant4 globals / types
00004 
00005 #include "G4Box.hh"
00006 #include "G4Polyhedra.hh"
00007 #include "G4SubtractionSolid.hh"
00008 #include "G4ThreeVector.hh"
00009 #include "G4Tubs.hh"
00010 #include "G4VSolid.hh"
00011 
00012 #include <cmath>
00013 
00014 BDSMagnetOuterFactoryPolesFacetCrop* BDSMagnetOuterFactoryPolesFacetCrop::_instance = 0;
00015 
00016 BDSMagnetOuterFactoryPolesFacetCrop* BDSMagnetOuterFactoryPolesFacetCrop::Instance()
00017 {
00018   if (_instance == 0)
00019     {_instance = new BDSMagnetOuterFactoryPolesFacetCrop();}
00020   return _instance;
00021 }
00022 
00023 BDSMagnetOuterFactoryPolesFacetCrop::~BDSMagnetOuterFactoryPolesFacetCrop()
00024 {
00025   _instance = 0;
00026 }
00027 
00028 // NOTE unforunately can't get this through inheritance as BDSMagnetOuterFactoryPolesFacet
00029 // is a singleton with private constructor - so this method is duplicated
00030 void BDSMagnetOuterFactoryPolesFacetCrop::CreatePoleSolid(G4String     name,
00031                                                           G4double     length,
00032                                                           G4int        order)
00033 {
00034   // use base class to do all the work, then modify the pole by cropping
00035   // it with a box to get the right shape
00036   
00037   // set pole length to be longer so we have unambiguous subtraction
00038   G4double tempPoleFinishRadius = poleFinishRadius; // copy to set back later
00039   poleFinishRadius = 2*poleFinishRadius;
00040   
00041   // base class method will use poleFinishRadius
00042   BDSMagnetOuterFactoryPolesBase::CreatePoleSolid(name,length,order);
00043   
00044   poleFinishRadius = tempPoleFinishRadius; // set it back to what it was
00045 
00046   G4VSolid* baseClassPoleSolid = poleSolid;
00047 
00048   // for box sizes we need something roughly adaptive to the component size
00049   // but just needs to be big - don't want to hard code anything so use
00050   // poleFinishRadius - a bigger component will entail bigger poleFinishRadius
00051   G4VSolid* subtractionBox = new G4Box(name + "_subtraction_box",  // name
00052                                        poleFinishRadius,           // x half width
00053                                        poleFinishRadius,           // y half width
00054                                        length);                    // z half width
00055   // z half width is full length for unambiguous subtraction
00056 
00057   // note translation is to centre of box which is nominal poleFinishRadius + half width
00058   // of box which is poleFinishRadius also - so tranlation is 2*poleFinishRadius
00059   G4ThreeVector boxTranslation(2*poleFinishRadius,0,0);
00060   poleSolid = new G4SubtractionSolid(name + "_pole_solid",        // name
00061                                      baseClassPoleSolid,          // solid 1
00062                                      subtractionBox,              // solid 2 - subtract this one
00063                                      0,                           // rotation
00064                                      boxTranslation);
00065                                      
00066 }
00067 
00068 void BDSMagnetOuterFactoryPolesFacetCrop::CreateYokeAndContainerSolid(G4String      name,
00069                                                                       G4double      length,
00070                                                                       G4int         order)
00071 {
00072   G4double segmentAngle = CLHEP::twopi / (2*order);
00073   
00074   G4double zPlanes[2] = {-length*0.5, length*0.5};
00075   std::vector<G4double> innerRadii;
00076   std::vector<G4double> outerRadii;
00077   for (G4int i = 0; i < order*4; ++i)
00078     {
00079       innerRadii.push_back(yokeStartRadius);
00080       outerRadii.push_back(yokeFinishRadius);
00081     }
00082   
00083   yokeSolid = new G4Polyhedra(name + "_yoke_solid",    // name
00084                               CLHEP::pi*0.5 + segmentAngle*0.25,    // start angle
00085                               CLHEP::twopi,            // sweep angle
00086                               4*order,                 // number of sides
00087                               2,                       // number of z planes
00088                               zPlanes,                 // z plane z coordinates
00089                               innerRadii.data(),
00090                               outerRadii.data());
00091 
00092   // note container must have hole in it for the beampipe to fit in!
00093   // poled geometry doesn't fit tightly to beampipe so can alays use a circular aperture
00094   G4double yokeExtent   = yokeFinishRadius / cos(segmentAngle*0.25);
00095   
00096   containerSolid = new G4Tubs(name + "_container_solid",       // name
00097                               poleStartRadius,                 // start radius
00098                               yokeExtent + lengthSafety,       // finish radius
00099                               length*0.5,                      // z half length
00100                               0,                               // start angle
00101                               CLHEP::twopi);                   // sweep angle
00102 
00103 
00104 }

Generated on 28 Jun 2015 for BDSIM by  doxygen 1.4.7