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

00001 #include "BDSMagnetOuterFactoryLHC.hh"
00002 
00003 #include "BDSBeamPipe.hh"
00004 #include "BDSBeamPipeType.hh"
00005 #include "BDSBeamPipeFactory.hh"
00006 #include "BDSDebug.hh"
00007 #include "BDSExecOptions.hh"
00008 #include "BDSGeometryComponent.hh"
00009 #include "BDSGlobalConstants.hh"
00010 #include "BDSMagnetColours.hh"
00011 #include "BDSMaterials.hh"
00012 #include "BDSSDManager.hh"
00013 #include "BDSUtilities.hh"                 // for calculateorientation
00014 
00015 #include "globals.hh"                      // geant4 globals / types
00016 #include "G4Box.hh"
00017 #include "G4Colour.hh"
00018 #include "G4CutTubs.hh"
00019 #include "G4IntersectionSolid.hh"
00020 #include "G4LogicalVolume.hh"
00021 #include "G4UnionSolid.hh"
00022 #include "G4Material.hh"
00023 #include "G4PVPlacement.hh"
00024 #include "G4SubtractionSolid.hh"
00025 #include "G4ThreeVector.hh"
00026 #include "G4Tubs.hh"
00027 #include "G4UserLimits.hh"
00028 #include "G4VisAttributes.hh"
00029 #include "G4VSolid.hh"
00030 #include <cmath>
00031 #include <utility>                         // for std::pair
00032 #include <algorithm>                       // for std::max
00033 #include <vector>
00034 
00035 
00036 BDSMagnetOuterFactoryLHC::BDSMagnetOuterFactoryLHC(G4bool isLeftOffsetIn):
00037   isLeftOffset(isLeftOffsetIn)
00038 {
00039   CleanUp();
00040 }
00041 
00042 BDSGeometryComponent* BDSMagnetOuterFactoryLHC::CreateSectorBend(G4String      name,
00043                                                                  G4double      length,
00044                                                                  BDSBeamPipe*  beamPipe,
00045                                                                  G4double      outerDiameter,
00046                                                                  G4double      angle,
00047                                                                  G4Material*   outerMaterial)
00048 
00049 {
00050 #ifdef BDSDEBUG
00051   G4cout << __METHOD_NAME__ << G4endl;
00052 #endif
00053   CleanUp();
00054   
00055   // note this geometry does not respond to outerDiameter - it's hard coded to the
00056   // design of a sector bend for the lhc.  TestInputParameters requires it though
00057   // to be the same check for the other methods
00058 
00059   // test input parameters - set global options as default if not specified
00060   TestInputParameters(beamPipe,outerDiameter,outerMaterial);
00061 
00062   // nominal lhc beampipe parameters for reference
00063   // aper1 = 4.404cm / 2
00064   // aper2 = 3.428cm / 2
00065   // aper3 = 4.404cm / 2
00066   // containerRadius -> 24.599mm for lhc beampipe with these parameters
00067 
00068   // geometrical constants
00069   // [1] LHC design report - Chapter 7, fig 7.3
00070   // [2] LHC design report - Chapter 7, fig 7.1
00071   G4double beamPipeAxisSeparation = 194.00*CLHEP::mm;             // at 1.9K
00072   G4double massShift              = 0.5 * beamPipeAxisSeparation; // central shift to geometry
00073   //G4double collarBoxHalfHeight    = 60*CLHEP::mm;                 // [1] by visual inspection
00074   G4double collarBoxHalfWidth     = 22*CLHEP::mm;                 // fits between outer coils
00075 
00076   // radii
00077   G4double containerInnerRadius   = beamPipe->GetContainerRadius()+1*CLHEP::um;
00078   G4double innerCoilInnerRadius   = containerInnerRadius + lengthSafety;
00079   G4double innerCoilInnerRadiusF  = 24.601*CLHEP::mm; // the fixed radius for the second pipe - F for fixed
00080   G4double innerCoilOuterRadius   = 42*CLHEP::mm;                        // [2] by visual inspection
00081   G4double outerCoilInnerRadius   = innerCoilOuterRadius + lengthSafety;
00082   G4double outerCoilInnerRadiusF  = innerCoilOuterRadius + lengthSafety; // doesn't change
00083   G4double outerCoilOuterRadius   = 60*CLHEP::mm;                        // [2] by visual inspection
00084   G4double collarInnerRadius      = outerCoilOuterRadius + lengthSafety;
00085   G4double collarInnerRadiusF     = outerCoilOuterRadius + lengthSafety;
00086   G4double collarOuterRadius      = 101.18*CLHEP::mm;                    // [1] - at 293K but don't have 1.9K measurement
00087   G4double yokeOuterRadius        = 570.0*0.5*CLHEP::mm;                 // [2] - at 293K but don't have 1.9K measurement
00088 
00089   // angular setup of coils
00090   // these angles were calculated by visually analysing the coil layout graph in [2]
00091   G4double poleInnerFullAngle = 33./180.*2; //33 degrees half angle in radians
00092   G4double poleOuterFullAngle = 100./180.*2; //80 degrees half angle in radians
00093   G4double coilInnerFullAngle = CLHEP::pi - poleInnerFullAngle - 1e-5; // 0.001 is small margin to avoid overlap
00094   G4double coilOuterFullAngle = CLHEP::pi - poleOuterFullAngle - 1e-5;
00095 
00096   // whether to build various components around active beam pipe depending on how wide it is
00097   // these ONLY apply to the components around the active beampipe
00098   G4bool buildInnerCoil           = true;
00099   G4bool buildOuterCoil           = true;
00100   G4bool buildCollar              = true; // collar around the active beam pipe
00101   if (innerCoilInnerRadius > innerCoilOuterRadius)
00102     {buildInnerCoil = false;}
00103   if ((innerCoilInnerRadius > outerCoilInnerRadius) && (innerCoilInnerRadius < outerCoilOuterRadius))
00104     {outerCoilInnerRadius = containerInnerRadius + lengthSafety;}
00105   if (innerCoilInnerRadius > outerCoilOuterRadius)
00106     {buildOuterCoil = false;}
00107   // this still uses the boxHalfWidth but just as good as the collar annulli overlap slightly in the middle
00108   // and this will protect against this
00109   if ((innerCoilInnerRadius > collarInnerRadius) && (innerCoilInnerRadius < (massShift - collarBoxHalfWidth)))
00110     {collarInnerRadius = containerInnerRadius + lengthSafety;}
00111   if (innerCoilInnerRadius > (massShift - collarBoxHalfWidth))
00112     {buildCollar = false;}
00113   if (innerCoilInnerRadius > collarOuterRadius)
00114     {
00115       // pipe is too big to use with this geometry!
00116       G4cerr << __METHOD_NAME__ << "this beam pipe is too big to use with the LHC dipole geometry" << G4endl;
00117       G4cerr << "Please consider using a different magnet geometry for this particular magnet" << G4endl;
00118       G4cerr << "Magnet named: " << name << G4endl;
00119       exit(1);
00120     }
00121 
00122   G4ThreeVector dipolePosition; // translation of whole assembly relative to centre of active pipe
00123   if (isLeftOffset)
00124     {
00125       dipolePosition = G4ThreeVector(massShift,0.,0.);
00126       beamPipeAxisSeparation  *= -1;
00127 #ifdef BDSDEBUG
00128       G4cout << __METHOD_NAME__ << "dipole to the left" << G4endl;
00129 #endif
00130     }
00131   else
00132     {
00133       dipolePosition = G4ThreeVector(-massShift,0.,0.);
00134       //massShift     *= -1;
00135 #ifdef BDSDEBUG
00136       G4cout << __METHOD_NAME__ << "dipole to the right" << G4endl;
00137 #endif
00138     }
00139 
00140   // calculate some geometrical parameters
00141   G4int orientation        = BDS::CalculateOrientation(angle);
00142   G4double zcomponent      = cos(fabs(angle*0.5)); // calculate components of normal vectors (in the end mag(normal) = 1)
00143   G4double xcomponent      = sin(fabs(angle*0.5)); // note full angle here as it's the exit angle
00144   G4ThreeVector inputface  = G4ThreeVector(-orientation*xcomponent, 0.0, -1.0*zcomponent); //-1 as pointing down in z for normal
00145   G4ThreeVector outputface = G4ThreeVector(-orientation*xcomponent, 0.0, zcomponent);   // no output face angle
00146 
00147   // lengths at different points transversely - dependent on left or right geometry as well as angle +ve or -ve
00148   G4double centralHalfLength  = length*0.5 - orientation*0.5*beamPipeAxisSeparation*tan(fabs(angle*0.5)); // central axis of outer cylinder
00149   G4double secondBPHalfLength = length*0.5 - orientation*beamPipeAxisSeparation*tan(fabs(angle*0.5));     // central axis of second beampipe
00150   
00151   // we make a lot of volumes in this class - keep a record for easily attaching properties to all
00152   std::vector<G4LogicalVolume*> allLogicalVolumes;
00153   
00154   if (beamPipe->ContainerIsCircular())
00155     {
00156       //circular beampipe so we can simply use its radius
00157       //container is similar but slightly wider and hollow (to allow placement of beampipe)
00158       //have to do subtraction as cuttubs aperture is central and the beampipe (active one) is not here
00159       G4VSolid* containerSolidOuter = new G4CutTubs(name + "_contiainer_solid_outer",  // name
00160                                                     0,                           // inner radius
00161                                                     yokeOuterRadius,             // outer radius
00162                                                     centralHalfLength,           // half length
00163                                                     0,                           // rotation start angle
00164                                                     CLHEP::twopi,                // rotation finish angle
00165                                                     inputface,                   // input face normal
00166                                                     outputface);                 // output face normal
00167       G4VSolid* containerSolidInner = new G4Tubs(name + "_contiainer_solid_inner",  // name
00168                                                  0,                                 // inner radius
00169                                                  containerInnerRadius,              // outer radius
00170                                                  length,                            // full length for unambiguous subtraction
00171                                                  0,                                 // rotation start angle
00172                                                  CLHEP::twopi);
00173       containerSolid = new G4SubtractionSolid(name + "_container_solid",   // name
00174                                               containerSolidOuter,         // outer bit
00175                                               containerSolidInner,         // subtract this from it
00176                                               0,                           // rotation
00177                                               -dipolePosition);            // translation
00178     }
00179   else
00180     {
00181       //container is similar but slightly wider
00182       G4VSolid* containerSolidOuter = new G4CutTubs(name + "_contiainer_solid_outer",  // name
00183                                                     0,                                 // inner radius
00184                                                     yokeOuterRadius,                   // outer radius
00185                                                     centralHalfLength,                 // half length
00186                                                     0,                                 // rotation start angle
00187                                                     CLHEP::twopi,                      // rotation finish angle
00188                                                     inputface,                         // input face normal
00189                                                     outputface);                       // output face normal
00190       containerSolid = new G4SubtractionSolid(name + "_container_solid",
00191                                               containerSolidOuter,
00192                                               beamPipe->GetContainerSubtractionSolid(),
00193                                               0,                // rotation
00194                                               -dipolePosition); // translation
00195     }
00196 
00197   G4Material* emptyMaterial = BDSMaterials::Instance()->GetMaterial(BDSGlobalConstants::Instance()->GetEmptyMaterial());
00198   G4LogicalVolume* containerLV = new G4LogicalVolume(containerSolid,
00199                                                      emptyMaterial,
00200                                                      name + "_container_lv");
00201     
00202   // coil solids
00203   G4VSolid*        coil1Inner   = NULL;
00204   G4VSolid*        coil1Outer   = NULL;
00205   G4VSolid*        coil2Inner   = NULL;
00206   G4VSolid*        coil2Outer   = NULL;
00207   G4LogicalVolume* coil1InnerLV = NULL;
00208   G4LogicalVolume* coil1OuterLV = NULL;
00209   G4LogicalVolume* coil2InnerLV = NULL;
00210   G4LogicalVolume* coil2OuterLV = NULL;
00211   G4Material*      nbti         = BDSMaterials::Instance()->GetMaterial("NbTi.1");
00212   G4Material* stainlesssteel    = BDSMaterials::Instance()->GetMaterial("stainlesssteel");
00213   G4VisAttributes* coilVisAtt   = new G4VisAttributes(G4Colour(0.9, 0.75, 0.)); //gold-ish colour
00214   coilVisAtt->SetForceLineSegmentsPerCircle(nSegmentsPerCircle);
00215   G4VSolid*        collar1PoleTopInnerSolid     = NULL;
00216   G4VSolid*        collar1PoleBottomInnerSolid  = NULL;
00217   G4VSolid*        collar1PoleTopOuterSolid     = NULL;
00218   G4VSolid*        collar1PoleBottomOuterSolid  = NULL;
00219   G4LogicalVolume* collar1PoleTopInnerLV        = NULL;
00220   G4LogicalVolume* collar1PoleBottomInnerLV     = NULL;
00221   G4LogicalVolume* collar1PoleTopOuterLV        = NULL;
00222   G4LogicalVolume* collar1PoleBottomOuterLV     = NULL;
00223   G4VisAttributes* collarVisAtt = new G4VisAttributes(G4Colour(0.9, 0.9, 0.9)); // grey
00224   collarVisAtt->SetForceLineSegmentsPerCircle(nSegmentsPerCircle);
00225 
00226   if (buildInnerCoil)
00227     {
00228       coil1Inner = new G4CutTubs(name+"_coil1_inner_solid",            // name
00229                                  innerCoilInnerRadius,                 // inner radius
00230                                  innerCoilOuterRadius,                 // outer radius
00231                                  length*0.5-2*lengthSafety,            // length
00232                                  -coilInnerFullAngle*0.5,              // start angle
00233                                  coilInnerFullAngle,                   // sweep angle
00234                                  inputface,                            // input face normal
00235                                  outputface);                          // output face normal
00236       coil2Inner = new G4CutTubs(name+"_coil2_inner_solid",            // name
00237                                  innerCoilInnerRadius,                 // inner radius
00238                                  innerCoilOuterRadius,                 // outer radius
00239                                  length*0.5-2*lengthSafety,            // length
00240                                  CLHEP::pi-coilInnerFullAngle*0.5,     // start angle
00241                                  coilInnerFullAngle,                   // sweep angle
00242                                  inputface,                            // input face normal
00243                                  outputface);                          // output face normal
00244       coil1InnerLV =  new G4LogicalVolume(coil1Inner,
00245                                           nbti,
00246                                           name+"_coil1_Inner_lv");
00247       coil2InnerLV =  new G4LogicalVolume(coil2Inner,
00248                                           nbti,
00249                                           name+"_coil2_Inner_lv");
00250       coil1InnerLV->SetVisAttributes(coilVisAtt);
00251       coil2InnerLV->SetVisAttributes(coilVisAtt);
00252       allLogicalVolumes.push_back(coil1InnerLV);
00253       allLogicalVolumes.push_back(coil2InnerLV);
00254 
00255       new G4PVPlacement(0,                      // rotation
00256                         -dipolePosition,        // position
00257                         coil1InnerLV,           // its logical volume
00258                         name+"_coil1_inner_pv", // its name
00259                         containerLV,            // its mother  volume
00260                         false,                  // no boolean operation
00261                         0, 
00262                         checkOverlaps);
00263       new G4PVPlacement(0,                      // rotation
00264                         -dipolePosition,        // position
00265                         coil2InnerLV,           // its logical volume
00266                         name+"_coil2_inner_pv", // its name
00267                         containerLV,            // its mother  volume
00268                         false,                  // no boolean operation
00269                         0, 
00270                         checkOverlaps);
00271 
00272       collar1PoleTopInnerSolid    = new G4CutTubs(name+"_collar1_pole_inner_top",      // name
00273                                                   innerCoilInnerRadius,                // inner radius
00274                                                   innerCoilOuterRadius,                // outer radius
00275                                                   length*0.5 - 2*lengthSafety,         // length
00276                                                   CLHEP::pi*0.5-poleInnerFullAngle*0.5,// start angle
00277                                                   poleInnerFullAngle,                  // sweep angle
00278                                                   inputface,                           // input face normal
00279                                                   outputface);                         // output face normal
00280       collar1PoleBottomInnerSolid = new G4CutTubs(name+"_collar1_pole_inner_bottom",   // name
00281                                                   innerCoilInnerRadius,                // inner radius
00282                                                   innerCoilOuterRadius,                // outer radius
00283                                                   length*0.5 - 2*lengthSafety,         // length
00284                                                   CLHEP::pi*1.5-poleInnerFullAngle*0.5,// start angle
00285                                                   poleInnerFullAngle,                  // sweep angle
00286                                                   inputface,                           // input face normal
00287                                                   outputface);                         // output face normal
00288       collar1PoleTopInnerLV    = new G4LogicalVolume(collar1PoleTopInnerSolid,
00289                                                      stainlesssteel,
00290                                                      name+"_collar1_pole_top_inner_lv");
00291       collar1PoleBottomInnerLV = new G4LogicalVolume(collar1PoleBottomInnerSolid,
00292                                                      stainlesssteel,
00293                                                      name+"_collar1_pole_bottom_inner_lv");
00294       
00295       collar1PoleTopInnerLV->SetVisAttributes(collarVisAtt);
00296       collar1PoleBottomInnerLV->SetVisAttributes(collarVisAtt);
00297       
00298       allLogicalVolumes.push_back(collar1PoleTopInnerLV); // register locally
00299       allLogicalVolumes.push_back(collar1PoleBottomInnerLV);
00300 
00301       new G4PVPlacement(0,                          // rotation
00302                         -dipolePosition,            // position
00303                         collar1PoleTopInnerLV,      // its logical volume
00304                         name+"_collar1_pole_top_inner_pv",// its name
00305                         containerLV,                // its mother  volume
00306                         false,                      // no boolean operation
00307                         0,
00308                         checkOverlaps);
00309       new G4PVPlacement(0,                          // rotation
00310                         -dipolePosition,            // position
00311                         collar1PoleBottomInnerLV,   // its logical volume
00312                         name+"_collar1_pole_top_inner_pv",// its name
00313                         containerLV,                // its mother  volume
00314                         false,                      // no boolean operation
00315                         0,
00316                         checkOverlaps);
00317     }
00318   
00319   if (buildOuterCoil)
00320     {
00321       coil1Outer = new G4CutTubs(name+"_coil1_outer_solid",            // name
00322                                  outerCoilInnerRadius,                 // inner radius
00323                                  outerCoilOuterRadius,                 // outer radius
00324                                  length*0.5-2*lengthSafety,            // length
00325                                  -coilOuterFullAngle*0.5,              // start angle
00326                                  coilOuterFullAngle,                   // sweep angle
00327                                  inputface,                            // input face normal
00328                                  outputface);                          // output face normal
00329       coil2Outer = new G4CutTubs(name+"_coil2_outer_solid",            // name
00330                                  outerCoilInnerRadius,                 // inner radius
00331                                  outerCoilOuterRadius,                 // outer radius
00332                                  length*0.5-2*lengthSafety,            // length
00333                                  CLHEP::pi-coilOuterFullAngle*0.5,     // start angle
00334                                  coilOuterFullAngle,                   // sweep angle
00335                                  inputface,                            // input face normal
00336                                  outputface);                          // output face normal
00337       coil1OuterLV =  new G4LogicalVolume(coil1Outer,
00338                                           nbti,
00339                                           name+"_coil1_Inner_lv");
00340       coil2OuterLV =  new G4LogicalVolume(coil2Outer,
00341                                           nbti,
00342                                           name+"_coil2_Inner_lv");
00343       coil1OuterLV->SetVisAttributes(coilVisAtt);
00344       coil2OuterLV->SetVisAttributes(coilVisAtt);
00345       allLogicalVolumes.push_back(coil1OuterLV);
00346       allLogicalVolumes.push_back(coil2OuterLV);
00347 
00348       new G4PVPlacement(0,                      // rotation
00349                         -dipolePosition,        // position
00350                         coil1OuterLV,           // its logical volume
00351                         name+"_coil1_outer_pv", // its name
00352                         containerLV,            // its mother  volume
00353                         false,                  // no boolean operation
00354                         0, 
00355                         checkOverlaps);
00356       new G4PVPlacement(0,                      // rotation
00357                         -dipolePosition,        // position
00358                         coil2OuterLV,           // its logical volume
00359                         name+"_coil2_outer_pv", // its name
00360                         containerLV,            // its mother  volume
00361                         false,                  // no boolean operation
00362                         0, 
00363                         checkOverlaps);
00364 
00365       collar1PoleTopOuterSolid    = new G4CutTubs(name+"_collar1_pole_outer_top",      // name
00366                                                   outerCoilInnerRadius,                // inner radius
00367                                                   outerCoilOuterRadius,                // outer radius
00368                                                   length*0.5 - 2*lengthSafety,         // length
00369                                                   CLHEP::pi*0.5-poleOuterFullAngle*0.5,// start angle
00370                                                   poleOuterFullAngle,                  // sweep angle
00371                                                   inputface,                           // input face normal
00372                                                   outputface);                         // output face normal
00373       collar1PoleBottomOuterSolid = new G4CutTubs(name+"_collar1_pole_outer_bottom",   // name
00374                                                   outerCoilInnerRadius,                // inner radius
00375                                                   outerCoilOuterRadius,                // outer radius
00376                                                   length*0.5 - 2*lengthSafety,         // length
00377                                                   CLHEP::pi*1.5-poleOuterFullAngle*0.5,// start angle
00378                                                   poleOuterFullAngle,                  // sweep angle
00379                                                   inputface,                           // input face normal
00380                                                   outputface);                         // output face normal
00381 
00382       collar1PoleTopOuterLV    = new G4LogicalVolume(collar1PoleTopOuterSolid,
00383                                                      stainlesssteel,
00384                                                      name+"_collar1_pole_top_outer_lv");
00385       collar1PoleBottomOuterLV = new G4LogicalVolume(collar1PoleBottomOuterSolid,
00386                                                      stainlesssteel,
00387                                                      name+"_collar1_pole_bottom_outer_lv");
00388       
00389       collar1PoleTopOuterLV->SetVisAttributes(collarVisAtt);
00390       collar1PoleBottomOuterLV->SetVisAttributes(collarVisAtt);
00391       
00392       allLogicalVolumes.push_back(collar1PoleTopOuterLV);
00393       allLogicalVolumes.push_back(collar1PoleBottomOuterLV);
00394 
00395       new G4PVPlacement(0,                                // rotation
00396                         -dipolePosition,                  // position
00397                         collar1PoleTopOuterLV,            // its logical volume
00398                         name+"_collar1_pole_top_inner_pv",// its name
00399                         containerLV,                      // its mother  volume
00400                         false,                            // no boolean operation
00401                         0,
00402                         checkOverlaps);
00403       new G4PVPlacement(0,                                // rotation
00404                         -dipolePosition,                  // position
00405                         collar1PoleBottomOuterLV,         // its logical volume
00406                         name+"_collar1_pole_top_inner_pv",// its name
00407                         containerLV,                      // its mother  volume
00408                         false,                            // no boolean operation
00409                         0,
00410                         checkOverlaps);
00411     }
00412   
00413   // coils on inactive beam pipe - always built
00414   G4VSolid *coil3Inner = new G4CutTubs(name+"_coil3_inner_solid",            // name
00415                                        innerCoilInnerRadiusF,                // inner radius
00416                                        innerCoilOuterRadius,                 // outer radius
00417                                        secondBPHalfLength-2*lengthSafety,    // length
00418                                        -coilInnerFullAngle*0.5,              // start angle
00419                                        coilInnerFullAngle,                   // sweep angle
00420                                        inputface,                            // input face normal
00421                                        outputface);                          // output face normal
00422   G4VSolid *coil3Outer = new G4CutTubs(name+"_coil3_outer_solid",            // name
00423                                        outerCoilInnerRadiusF,                // inner radius
00424                                        outerCoilOuterRadius,                 // outer radius
00425                                        secondBPHalfLength-2*lengthSafety,    // length
00426                                        -coilOuterFullAngle*0.5,              // start angle
00427                                        coilOuterFullAngle,                   // sweep angle
00428                                        inputface,                            // input face normal
00429                                        outputface);                          // output face normal
00430   G4VSolid *coil4Inner = new G4CutTubs(name+"_coil4_inner_solid",            // name
00431                                        innerCoilInnerRadiusF,                // inner radius
00432                                        innerCoilOuterRadius,                 // outer radius
00433                                        secondBPHalfLength-2*lengthSafety,    // length
00434                                        CLHEP::pi-coilInnerFullAngle*0.5,     // start angle
00435                                        coilInnerFullAngle,                   // sweep angle
00436                                        inputface,                            // input face normal
00437                                        outputface);                          // output face normal
00438   G4VSolid *coil4Outer = new G4CutTubs(name+"_coil4_outer_solid",            // name
00439                                        outerCoilInnerRadiusF,                // inner radius
00440                                        outerCoilOuterRadius,                 // outer radius
00441                                        secondBPHalfLength-2*lengthSafety,    // length
00442                                        CLHEP::pi-coilOuterFullAngle*0.5,     // start angle
00443                                        coilOuterFullAngle,                   // sweep angle
00444                                        inputface,                            // input face normal
00445                                        outputface);                          // output face normal
00446   
00447   G4LogicalVolume* coil3InnerLV =  new G4LogicalVolume(coil3Inner,
00448                                                        nbti,
00449                                                        name+"_coil3_Inner_lv");
00450   G4LogicalVolume* coil3OuterLV =  new G4LogicalVolume(coil3Outer,
00451                                                        nbti,
00452                                                        name+"_coil3_Inner_lv");
00453   G4LogicalVolume* coil4InnerLV =  new G4LogicalVolume(coil4Inner,
00454                                                        nbti,
00455                                                        name+"_coil4_Inner_lv");
00456   G4LogicalVolume* coil4OuterLV =  new G4LogicalVolume(coil4Outer,
00457                                                        nbti,
00458                                                        name+"_coil4_Inner_lv");
00459   
00460   coil3InnerLV->SetVisAttributes(coilVisAtt);
00461   coil3OuterLV->SetVisAttributes(coilVisAtt);
00462   coil4InnerLV->SetVisAttributes(coilVisAtt);
00463   coil4OuterLV->SetVisAttributes(coilVisAtt);
00464   
00465   allLogicalVolumes.push_back(coil3InnerLV);
00466   allLogicalVolumes.push_back(coil3OuterLV);
00467   allLogicalVolumes.push_back(coil4InnerLV);
00468   allLogicalVolumes.push_back(coil4OuterLV);
00469 
00470   // coil placement  
00471   new G4PVPlacement(0,                      // rotation
00472                     dipolePosition,         // position
00473                     coil3InnerLV,           // its logical volume
00474                     name+"_coil3_inner_pv", // its name
00475                     containerLV,            // its mother  volume
00476                     false,                  // no boolean operation
00477                     0, 
00478                     checkOverlaps);
00479   new G4PVPlacement(0,                      // rotation
00480                     dipolePosition,         // position
00481                     coil3OuterLV,           // its logical volume
00482                     name+"_coil3_outer_pv", // its name
00483                     containerLV,            // its mother  volume
00484                     false,                  // no boolean operation
00485                     0, 
00486                     checkOverlaps);
00487   new G4PVPlacement(0,                      // rotation
00488                     dipolePosition,         // position
00489                     coil4InnerLV,           // its logical volume
00490                     name+"_coil4_inner_pv", // its name
00491                     containerLV,            // its mother  volume
00492                     false,                  // no boolean operation
00493                     0, 
00494                     checkOverlaps);
00495   new G4PVPlacement(0,                      // rotation
00496                     dipolePosition,         // position
00497                     coil4OuterLV,           // its logical volume
00498                     name+"_coil4_outer_pv", // its name
00499                     containerLV,            // its mother  volume
00500                     false,                  // no boolean operation
00501                     0, 
00502                     checkOverlaps);
00503   
00504   // non-magnetic collars
00505   // collar pole solids  
00506   G4VSolid* collar2PoleTopInnerSolid    = new G4CutTubs(name+"_collar2_pole_inner_top",      // name
00507                                                         innerCoilInnerRadiusF,               // inner radius
00508                                                         innerCoilOuterRadius,                // outer radius
00509                                                         secondBPHalfLength-2*lengthSafety,   // length
00510                                                         CLHEP::pi*0.5-poleInnerFullAngle*0.5,// start angle
00511                                                         poleInnerFullAngle,                  // sweep angle
00512                                                         inputface,                           // input face normal
00513                                                         outputface);                         // output face normal
00514   G4VSolid* collar2PoleTopOuterSolid    = new G4CutTubs(name+"_collar2_pole_outer_top",      // name
00515                                                         outerCoilInnerRadiusF,               // inner radius
00516                                                         outerCoilOuterRadius,                // outer radius
00517                                                         secondBPHalfLength-2*lengthSafety,   // length
00518                                                         CLHEP::pi*0.5-poleOuterFullAngle*0.5,// start angle
00519                                                         poleOuterFullAngle,                  // sweep angle
00520                                                         inputface,                           // input face normal
00521                                                         outputface);                         // output face normal
00522   G4VSolid* collar2PoleBottomInnerSolid = new G4CutTubs(name+"_collar2_pole_inner_bottom",   // name
00523                                                         innerCoilInnerRadiusF,               // inner radius
00524                                                         innerCoilOuterRadius,                // outer radius
00525                                                         secondBPHalfLength-2*lengthSafety,   // length
00526                                                         CLHEP::pi*1.5-poleInnerFullAngle*0.5,// start angle
00527                                                         poleInnerFullAngle,                  // sweep angle
00528                                                         inputface,                           // input face normal
00529                                                         outputface);                         // output face normal
00530   G4VSolid* collar2PoleBottomOuterSolid = new G4CutTubs(name+"_collar2_pole_outer_bottom",   // name
00531                                                         outerCoilInnerRadiusF,               // inner radius
00532                                                         outerCoilOuterRadius,                // outer radius
00533                                                         secondBPHalfLength-2*lengthSafety,   // length
00534                                                         CLHEP::pi*1.5-poleOuterFullAngle*0.5,// start angle
00535                                                         poleOuterFullAngle,                  // sweep angle
00536                                                         inputface,                           // input face normal
00537                                                         outputface);                         // output face normal
00538   
00539   // collar pole logical volumes
00540   G4LogicalVolume* collar2PoleTopInnerLV    = new G4LogicalVolume(collar2PoleTopInnerSolid,
00541                                                                   stainlesssteel,
00542                                                                   name+"_collar2_pole_top_inner_lv");
00543   G4LogicalVolume* collar2PoleTopOuterLV    = new G4LogicalVolume(collar2PoleTopOuterSolid,
00544                                                                   stainlesssteel,
00545                                                                   name+"_collar2_pole_top_outer_lv");
00546   G4LogicalVolume* collar2PoleBottomInnerLV = new G4LogicalVolume(collar2PoleBottomInnerSolid,
00547                                                                   stainlesssteel,
00548                                                                   name+"_collar2_pole_bottom_inner_lv");
00549   G4LogicalVolume* collar2PoleBottomOuterLV = new G4LogicalVolume(collar2PoleBottomOuterSolid,
00550                                                                   stainlesssteel,
00551                                                                   name+"_collar2_pole_bottom_outer_lv");
00552 
00553   // collar pole visualisation
00554   collar2PoleTopInnerLV->SetVisAttributes(collarVisAtt);
00555   collar2PoleTopOuterLV->SetVisAttributes(collarVisAtt);
00556   collar2PoleBottomInnerLV->SetVisAttributes(collarVisAtt);
00557   collar2PoleBottomOuterLV->SetVisAttributes(collarVisAtt);
00558   
00559   allLogicalVolumes.push_back(collar2PoleTopInnerLV);
00560   allLogicalVolumes.push_back(collar2PoleTopOuterLV);
00561   allLogicalVolumes.push_back(collar2PoleBottomInnerLV);
00562   allLogicalVolumes.push_back(collar2PoleBottomOuterLV);
00563   
00564   // collar pole placement
00565   new G4PVPlacement(0,                                // rotation
00566                     dipolePosition,                   // position
00567                     collar2PoleTopInnerLV,            // its logical volume
00568                     name+"_collar2_pole_top_inner_pv",// its name
00569                     containerLV,                      // its mother  volume
00570                     false,                            // no boolean operation
00571                     0,
00572                     checkOverlaps);
00573   new G4PVPlacement(0,                                // rotation
00574                     dipolePosition,                   // position
00575                     collar2PoleTopOuterLV,            // its logical volume
00576                     name+"_collar2_pole_top_inner_pv",// its name
00577                     containerLV,                      // its mother  volume
00578                     false,                            // no boolean operation
00579                     0,
00580                     checkOverlaps);
00581   new G4PVPlacement(0,                                // rotation
00582                     dipolePosition,                   // position
00583                     collar2PoleBottomInnerLV,         // its logical volume
00584                     name+"_collar2_pole_top_inner_pv",// its name
00585                     containerLV,                      // its mother  volume
00586                     false,                            // no boolean operation
00587                     0,
00588                     checkOverlaps);
00589   new G4PVPlacement(0,                                // rotation
00590                     dipolePosition,                   // position
00591                     collar2PoleBottomOuterLV,         // its logical volume
00592                     name+"_collar2_pole_top_inner_pv",// its name
00593                     containerLV,                      // its mother  volume
00594                     false,                            // no boolean operation
00595                     0,
00596                     checkOverlaps);
00597   
00598   // outer annulus of collar - two as slightly different lengths
00599   G4VSolid* collarAnnulus2 = new G4CutTubs(name+"_collar2_annulus_solid",    // name
00600                                            collarInnerRadiusF,               // inner radius
00601                                            collarOuterRadius,                // outer radius
00602                                            secondBPHalfLength-lengthSafety,  // length
00603                                            0,                                // starting angle
00604                                            CLHEP::twopi,                     // angle of sweep
00605                                            inputface,                        // input face normal
00606                                            outputface);                      // output face normal
00607   
00608   // make final solid pointer as collar round active beam pipe optional depending on how big active beam pipe is
00609   G4VSolid* collars = collarAnnulus2;
00610   
00611   if (buildCollar)
00612     {
00613       // building collar around active pipe
00614       G4VSolid* collarAnnulus1 = new G4CutTubs(name+"_collar1_annulus_solid",      // name
00615                                                collarInnerRadius,                  // inner radius
00616                                                collarOuterRadius,                  // outer radius
00617                                                length*0.5-2*lengthSafety,          // length
00618                                                0,                                  // starting angle
00619                                                CLHEP::twopi,                       // angle of sweep
00620                                                inputface,                          // input face normal
00621                                                outputface);                        // output face normal
00622 
00623       collars = new G4UnionSolid(name + "_collars_solid", // name
00624                                  collarAnnulus2,          // solid1
00625                                  collarAnnulus1,          // solid2
00626                                  0,                       // rotation
00627                                  -2*dipolePosition);      // translation
00628     }
00629 
00630   /*
00631   // This part seems to not produce any overlapping volumes and no errors but won't render
00632   // in anything but the raytracer. Requires commented out collarBoxHalfHeight at top of this
00633   // method to be uncommented, plus new subtraction solid for yoke will need to written
00634   // ommitted for now
00635   
00636   G4VSolid* collarBox      = new G4Box(name + "_collar_box_solid",           // name
00637                                        collarBoxHalfWidth,                   // x half width
00638                                        collarBoxHalfHeight,                  // y half width
00639                                        2*centralHalfLength);                 // z half length
00640   G4VSolid* collarBoxFaces = new G4CutTubs(name + "_collar_box_faces_solid", // name
00641                                            0,                                // inner radius
00642                                            50*CLHEP::cm,                     // outer radius
00643                                            centralHalfLength - lengthSafety, // length
00644                                            0,                                // starting angle
00645                                            CLHEP::twopi,                     // sweep angle
00646                                            inputface,                        // input face normal
00647                                            outputface);                      // output face normal
00648   G4VSolid* collarCentralPiece = new G4IntersectionSolid(name + "_collar_central_solid", // name
00649                                                          collarBox,                      // solid 1
00650                                                          collarBoxFaces);                // solid 2   
00651   
00652   G4VSolid* collarTotal = new G4UnionSolid(name + "_collar2_plus_box_solid", // name
00653                                            collars,                          // solid 2
00654                                            collarCentralPiece,               // solid 1
00655                                            0,                                // rotation
00656                                            -dipolePosition);    // translation
00657   */
00658   G4VSolid* collarTotal = collars;
00659   
00660   G4LogicalVolume *collarsLV =  new G4LogicalVolume(collarTotal,
00661                                                     stainlesssteel,
00662                                                     name+"_collars_lv");
00663 
00664   // collar annulus visualisation attributes
00665   collarsLV->SetVisAttributes(collarVisAtt);
00666 
00667   allLogicalVolumes.push_back(collarsLV); // register locally
00668   
00669   new G4PVPlacement(0,                  // rotation
00670                     dipolePosition,     // position
00671                     collarsLV,          // its logical volume
00672                     name+"_collars_pv", // its name
00673                     containerLV,        // its mother  volume
00674                     false,              // no boolean operation
00675                     0,                  // copy number
00676                     checkOverlaps);  
00677   
00678   // outer iron yoke
00679   G4VSolid* yokeCylinder = new G4CutTubs(name+"_yoke_cylinder_solid",     // name
00680                                          0.,                              // inner radius
00681                                          yokeOuterRadius - lengthSafety,  // outer radius
00682                                          centralHalfLength-2*lengthSafety,// length
00683                                          0,                               // starting angle
00684                                          CLHEP::twopi,                    // sweep angle
00685                                          inputface,                       // input face normal
00686                                          outputface);                     // output face normal
00687 
00688   // need to cut hole out for everything inside - note subtraction solid has to be solid
00689   G4VSolid* yokeSubtractionCylinder = new G4Tubs(name + "_yoke_subtraction_cyl_solid", // name
00690                                                  0,                                    // inner radius
00691                                                  collarOuterRadius + lengthSafety,     // outer radius
00692                                                  length,                               // z half length - long for unamibiuous subtraction
00693                                                  0,                                    // start angle
00694                                                  CLHEP::twopi);                        // sweep angle
00695 
00696   G4VSolid* yokeSubtractionSolid = new G4UnionSolid(name + "_yoke_subtraction_solid",  // name
00697                                                     yokeSubtractionCylinder,           // solid 1
00698                                                     yokeSubtractionCylinder,           // solid 2
00699                                                     0,                                 // rotation
00700                                                     2*dipolePosition);                 // translation
00701 
00702   G4VSolid* yoke = new G4SubtractionSolid(name+"_yoke_solid",           // name
00703                                           yokeCylinder,                 // from this
00704                                           yokeSubtractionSolid,         // subtract this
00705                                           0,                            // rotation
00706                                           -dipolePosition);              // translation
00707   
00708   G4LogicalVolume* yokeLV = new G4LogicalVolume(yoke,
00709                                                 BDSMaterials::Instance()->GetMaterial("Iron"),
00710                                                 name+"_yoke_lv");
00711 
00712   // yoke visualisation
00713   G4VisAttributes* LHCblue = new G4VisAttributes(G4Colour(0.0, 0.5, 1.0));
00714   LHCblue->SetForceLineSegmentsPerCircle(nSegmentsPerCircle);
00715   yokeLV->SetVisAttributes(LHCblue);
00716   
00717   allLogicalVolumes.push_back(yokeLV); // register locally
00718 
00719   // yoke placement
00720   new G4PVPlacement((G4RotationMatrix*)0,         // no rotation
00721                     G4ThreeVector(0,0,0),         // position
00722                     yokeLV,                       // lv to be placed
00723                     name + "_yoke_pv",            // name
00724                     containerLV,                  // mother lv to be place in
00725                     false,                        // no boolean operation
00726                     0,                            // copy number
00727                     checkOverlaps);
00728   
00729   G4String defaultMaterialName = BDSGlobalConstants::Instance()->GetBeamPipeMaterialName();
00730   G4Material* beamPipeMaterial = BDSMaterials::Instance()->GetMaterial(defaultMaterialName);
00731   G4Material* vacuumMaterial   = BDSMaterials::Instance()->GetMaterial(BDSGlobalConstants::Instance()->GetVacuumMaterial());
00732 
00733   //use beampipe factories to create another beampipe (note no magnetic field for now...)  
00734   BDSBeamPipe* secondBP = BDSBeamPipeFactory::Instance()->CreateBeamPipeAngledInOut(BDSBeamPipeType::lhcdetailed,
00735                                                                                     name,
00736                                                                                     2*secondBPHalfLength-2*lengthSafety,
00737                                                                                     -angle*0.5,        // entrane angle
00738                                                                                     -angle*0.5,        // exit angle
00739                                                                                     2.202*CLHEP::cm,   // aper1
00740                                                                                     1.714*CLHEP::cm,   // aper2
00741                                                                                     2.202*CLHEP::cm,   // aper3
00742                                                                                     0,                 // aper4
00743                                                                                     vacuumMaterial,    // vacuum material
00744                                                                                     1*CLHEP::mm,       // beampipeThickness
00745                                                                                     beamPipeMaterial); // beampipe material
00746   
00747   G4LogicalVolume* secondBPLV = secondBP->GetContainerLogicalVolume();
00748   allLogicalVolumes.push_back(secondBPLV);
00749   new G4PVPlacement((G4RotationMatrix*)0,         // no rotation
00750                     dipolePosition,               // position
00751                     secondBPLV,                   // lv to be placed
00752                     name + "_second_beampipe_pv", // name
00753                     containerLV,                  // mother lv to be place in
00754                     false,                        // no boolean operation
00755                     0,                            // copy number
00756                     checkOverlaps);
00757   
00758   // visual attributes for container
00759   if (BDSExecOptions::Instance()->GetVisDebug())
00760     {containerLV->SetVisAttributes(BDSGlobalConstants::Instance()->GetVisibleDebugVisAttr());}
00761   else
00762     {containerLV->SetVisAttributes(BDSGlobalConstants::Instance()->GetInvisibleVisAttr());}
00763 
00764   // USER LIMITS for all components
00765 #ifndef NOUSERLIMITS
00766   if (!allLogicalVolumes.empty()) {
00767     G4UserLimits* userLimits = new G4UserLimits("outer_cuts");
00768     userLimits->SetMaxAllowedStep( length * maxStepFactor );
00769     userLimits->SetUserMinEkine(BDSGlobalConstants::Instance()->GetThresholdCutCharged());
00770     userLimits->SetUserMaxTime(BDSGlobalConstants::Instance()->GetMaxTime());
00771     
00772     for (std::vector<G4LogicalVolume*>::iterator i = allLogicalVolumes.begin(); i != allLogicalVolumes.end(); ++i)
00773       {(*i)->SetUserLimits(userLimits);}
00774   }
00775 #endif
00776   
00777   // record extents
00778   // container radius is the same for all methods as all cylindrical
00779   G4double containerRadius = yokeOuterRadius;
00780   // massShift defined at very beginning of this function
00781   std::pair<double,double> extX = std::make_pair(-containerRadius+massShift,containerRadius+massShift); 
00782   std::pair<double,double> extY = std::make_pair(-containerRadius,containerRadius);
00783   std::pair<double,double> extZ = std::make_pair(-length*0.5,length*0.5);
00784   
00785   // build the BDSGeometryComponent instance and return it
00786   BDSGeometryComponent* outer = new BDSGeometryComponent(containerSolid,
00787                                                          containerLV,
00788                                                          extX, extY, extZ,
00789                                                          dipolePosition);
00790   // REGISTER all lvs
00791   outer->RegisterLogicalVolumes(secondBP->GetAllLogicalVolumes());
00792   outer->RegisterLogicalVolumes(allLogicalVolumes);
00793 
00794   // copy sensitive volumes if they exist
00795   outer->RegisterSensitiveVolumes(allLogicalVolumes);
00796   outer->RegisterSensitiveVolumes(secondBP->GetAllSensitiveVolumes());
00797 
00798   // allLogicalVolumes is a local variable and goes out of scope so doesn't
00799   // need to be emptied or reset here
00800   
00801   return outer;
00802 }
00803 
00804 BDSGeometryComponent* BDSMagnetOuterFactoryLHC::CreateRectangularBend(G4String      name,
00805                                                                       G4double      length,
00806                                                                       BDSBeamPipe*  beamPipe,
00807                                                                       G4double      boxSize,
00808                                                                       G4double      /*angle*/,
00809                                                                       G4Material*   outerMaterial)
00810 {
00811 #ifdef BDSDEBUG
00812   G4cout << __METHOD_NAME__ << G4endl;
00813 #endif
00814   CleanUp();
00815   //rectangular bends currently just make a shorter straight volume, so ignore angle for now
00816   CreateCylindricalSolids(name, length, beamPipe, boxSize);
00817   return CommonFinalConstructor(name, length, boxSize, outerMaterial, BDSMagnetColours::Instance()->GetMagnetColour("rectangularbend"));
00818 }
00819 
00820 BDSGeometryComponent* BDSMagnetOuterFactoryLHC::CreateQuadrupole(G4String      name,
00821                                                                  G4double      length,
00822                                                                  BDSBeamPipe*  beamPipe,
00823                                                                  G4double      outerDiameter,
00824                                                                  G4Material*   outerMaterial)
00825 
00826 {
00827 #ifdef BDSDEBUG
00828   G4cout << __METHOD_NAME__ << G4endl;
00829 #endif
00830   CleanUp();
00831 
00832   // test input parameters - set global options as default if not specified
00833   TestInputParameters(beamPipe,outerDiameter,outerMaterial);
00834 
00835   // geometrical constants
00836   // [1] LHC design report - Chapter 7, fig 7.3
00837   // [2] LHC design report - Chapter 7, fig 7.1
00838   G4double beamPipeAxisSeparation = 194.00*CLHEP::mm;               // at 1.9K
00839   G4double massShift              = 0.5 * beamPipeAxisSeparation;   // central shift to geometry
00840   //G4double collarBoxHalfHeight    = 60*CLHEP::mm;                 // [1] by visual inspection
00841   G4double collarBoxHalfWidth     = 22*CLHEP::mm;                   // fits between outer coils
00842 
00843   // radii
00844   G4double containerInnerRadius   = beamPipe->GetContainerRadius()+1*CLHEP::um;
00845   G4double coilInnerRadius        = containerInnerRadius + lengthSafety;
00846   G4double coilInnerRadiusF       = 24.601*CLHEP::mm; // the fixed radius for the second pipe - F for fixed
00847   G4double coilOuterRadius        = 60*CLHEP::mm;                   // [2] by visual inspection
00848   G4double collarInnerRadius      = coilOuterRadius + lengthSafety;
00849   G4double collarInnerRadiusF     = coilOuterRadius + lengthSafety;
00850   G4double collarOuterRadius      = 101.18*CLHEP::mm;               // [1] - at 293K but don't have 1.9K measurement
00851   G4double yokeOuterRadius        = 570.0*0.5*CLHEP::mm;            // [2] - at 293K but don't have 1.9K measurement
00852 
00853   // angular setup of coils
00854   // these angles were calculated by visually analysing the coil layout graph in [2]
00855   G4double poleFullAngle    = CLHEP::pi/6.; //30 degrees angle in radians
00856   G4double coilFullAngle    = CLHEP::pi/2. - poleFullAngle - 1e-5; // 1e-5 to prevent overlaps
00857   G4double coilHalfAngle    = coilFullAngle*0.5;
00858   G4double coilStartAngle   = -coilHalfAngle;
00859   G4double poleStartAngle   = coilHalfAngle;
00860   G4RotationMatrix* coil2rm = new G4RotationMatrix();
00861   G4RotationMatrix* coil3rm = new G4RotationMatrix();
00862   G4RotationMatrix* coil4rm = new G4RotationMatrix();
00863   coil2rm->rotateZ(CLHEP::pi/2.0);
00864   coil3rm->rotateZ(CLHEP::pi);
00865   coil4rm->rotateZ(CLHEP::pi*1.5);
00866 
00867   // whether to build various components around active beam pipe depending on how wide it is
00868   // these ONLY apply to the components around the active beampipe
00869   G4bool buildCoil         = true;
00870   G4bool buildCollar       = true; // collar around the active beam pipe
00871   if (coilInnerRadius > coilOuterRadius)
00872     {buildCoil = false;}
00873   // this still uses the boxHalfWidth but just as good as the collar annulli overlap slightly in the middle
00874   // and this will protect against this
00875   if ((coilInnerRadius > collarInnerRadius) && (coilInnerRadius < (massShift - collarBoxHalfWidth)))
00876     {collarInnerRadius = containerInnerRadius + lengthSafety;}
00877   if (coilInnerRadius > (massShift - collarBoxHalfWidth))
00878     {buildCollar = false;}
00879   if (coilInnerRadius > collarOuterRadius)
00880     {
00881       // pipe is too big to use with this geometry!
00882       G4cerr << __METHOD_NAME__ << "this beam pipe is too big to use with the LHC dipole geometry" << G4endl;
00883       G4cerr << "Please consider using a different magnet geometry for this particular magnet" << G4endl;
00884       G4cerr << "Magnet named: " << name << G4endl;
00885       exit(1);
00886     }
00887 
00888   G4ThreeVector dipolePosition; // translation of whole assembly relative to centre of active pipe
00889   if (isLeftOffset)
00890     {
00891       dipolePosition = G4ThreeVector(massShift,0.,0.);
00892       beamPipeAxisSeparation  *= -1;
00893 #ifdef BDSDEBUG
00894       G4cout << __METHOD_NAME__ << "dipole to the left" << G4endl;
00895 #endif
00896     }
00897   else
00898     {
00899       dipolePosition = G4ThreeVector(-massShift,0.,0.);
00900       //massShift     *= -1;
00901 #ifdef BDSDEBUG
00902       G4cout << __METHOD_NAME__ << "dipole to the right" << G4endl;
00903 #endif
00904     }
00905 
00906   // we make a lot of volumes in this class - keep a record for easily attaching properties to all
00907   std::vector<G4LogicalVolume*> allLogicalVolumes;
00908   
00909   if (beamPipe->ContainerIsCircular())
00910     {
00911       //circular beampipe so we can simply use its radius
00912       //container is similar but slightly wider and hollow (to allow placement of beampipe)
00913       //have to do subtraction as cuttubs aperture is central and the beampipe (active one) is not here
00914       G4VSolid* containerSolidOuter = new G4Tubs(name + "_contiainer_solid_outer",  // name
00915                                                  0,                           // inner radius
00916                                                  yokeOuterRadius,             // outer radius
00917                                                  length*0.5,                  // half length
00918                                                  0,                           // rotation start angle
00919                                                  CLHEP::twopi);               // rotation finish angle
00920                                                     
00921       G4VSolid* containerSolidInner = new G4Tubs(name + "_contiainer_solid_inner",  // name
00922                                                  0,                                 // inner radius
00923                                                  containerInnerRadius,              // outer radius
00924                                                  length,                            // full length for unambiguous subtraction
00925                                                  0,                                 // rotation start angle
00926                                                  CLHEP::twopi);
00927       containerSolid = new G4SubtractionSolid(name + "_container_solid",   // name
00928                                               containerSolidOuter,         // outer bit
00929                                               containerSolidInner,         // subtract this from it
00930                                               0,                           // rotation
00931                                               -dipolePosition);            // translation
00932     }
00933   else
00934     {
00935       //container is similar but slightly wider
00936       G4VSolid* containerSolidOuter = new G4Tubs(name + "_contiainer_solid_outer",  // name
00937                                                  0,                                 // inner radius
00938                                                  yokeOuterRadius,                   // outer radius
00939                                                  length*0.5,                        // half length
00940                                                  0,                                 // rotation start angle
00941                                                  CLHEP::twopi);                     // rotation finish angle
00942                                                  
00943       containerSolid = new G4SubtractionSolid(name + "_container_solid",
00944                                               containerSolidOuter,
00945                                               beamPipe->GetContainerSubtractionSolid(),
00946                                               0,                // rotation
00947                                               -dipolePosition); // translation
00948     }
00949 
00950   G4Material* emptyMaterial = BDSMaterials::Instance()->GetMaterial(BDSGlobalConstants::Instance()->GetEmptyMaterial());
00951   G4LogicalVolume* containerLV = new G4LogicalVolume(containerSolid,
00952                                                      emptyMaterial,
00953                                                      name + "_container_lv");
00954 
00955   // coil solids
00956   // only need one pole & coil per beampipe which can be repeatedly placed
00957   // just inner radius can be different depending on active beam pipe hence two
00958   G4VSolid*        coil1     = NULL;
00959   G4VSolid*        coil2     = NULL;
00960   G4LogicalVolume* coil1LV   = NULL;
00961   G4LogicalVolume* coil2LV   = NULL;
00962 
00963   // pole solids are the ones immediately after coilN clockwise
00964   G4VSolid*        pole1     = NULL;
00965   G4VSolid*        pole2     = NULL;
00966   G4LogicalVolume* pole1LV   = NULL;
00967   G4LogicalVolume* pole2LV   = NULL;
00968 
00969   // collar solid - only need one as can union with itself but offset
00970   G4VSolid*        collar    = NULL;
00971   G4VSolid*        collars   = NULL; // union of two collars
00972   G4LogicalVolume* collarsLV = NULL;
00973 
00974   // materials and visual attributes
00975   G4Material* nbti              = BDSMaterials::Instance()->GetMaterial("NbTi.1");
00976   G4Material* iron              = BDSMaterials::Instance()->GetMaterial("iron");
00977   G4Material* stainlesssteel    = BDSMaterials::Instance()->GetMaterial("stainlesssteel");
00978   G4VisAttributes* coilVisAtt   = new G4VisAttributes(G4Colour(0.9, 0.75, 0.)); //gold-ish colour
00979   coilVisAtt->SetForceLineSegmentsPerCircle(nSegmentsPerCircle);
00980   G4VisAttributes* collarVisAtt = new G4VisAttributes(G4Colour(0.9, 0.9, 0.9)); // grey
00981   collarVisAtt->SetForceLineSegmentsPerCircle(nSegmentsPerCircle);
00982   
00983   if (buildCoil)
00984     {
00985       // coil solid
00986       coil1 = new G4Tubs(name+"_coil1_solid",          // name
00987                          coilInnerRadius,              // inner radius
00988                          coilOuterRadius,              // outer radius
00989                          length*0.5-2*lengthSafety,    // z half length
00990                          coilStartAngle,               // start angle
00991                          coilFullAngle);               // sweep angle
00992       // coil logical volumes
00993       coil1LV = new G4LogicalVolume(coil1,             // solid
00994                                     nbti,              // material
00995                                     name+"_coil1_lv"); // name
00996       coil1LV->SetVisAttributes(coilVisAtt);
00997       allLogicalVolumes.push_back(coil1LV);
00998 
00999       // pole solid
01000       pole1 = new G4Tubs(name+"_pole1_solid",          // name
01001                          coilInnerRadius,              // inner radius
01002                          coilOuterRadius,              // outer radius
01003                          length*0.5 - lengthSafety,    // z half length
01004                          poleStartAngle,               // start angle
01005                          poleFullAngle);               // sweep angle
01006       pole1LV = new G4LogicalVolume(pole1,             // solid
01007                                     stainlesssteel,    // material
01008                                     name+"_pole1_lv"); // name
01009       pole1LV->SetVisAttributes(collarVisAtt);
01010       allLogicalVolumes.push_back(pole1LV);
01011 
01012       // coil placements
01013       new G4PVPlacement(0,                  // rotation
01014                         -dipolePosition,    // position
01015                         coil1LV,            // logical volume
01016                         name + "_coil1_pv", // name
01017                         containerLV,        // mother volume
01018                         false,              // boolean operation
01019                         0,                  // copy number
01020                         checkOverlaps);
01021       new G4PVPlacement(coil2rm,            // rotation
01022                         -dipolePosition,    // position
01023                         coil1LV,            // logical volume
01024                         name + "_coil2_pv", // name
01025                         containerLV,        // mother volume
01026                         false,              // boolean operation
01027                         0,                  // copy number
01028                         checkOverlaps);
01029       new G4PVPlacement(coil3rm,            // rotation
01030                         -dipolePosition,    // position
01031                         coil1LV,            // logical volume
01032                         name + "_coil3_pv", // name
01033                         containerLV,        // mother volume
01034                         false,              // boolean operation
01035                         0,                  // copy number
01036                         checkOverlaps);
01037       new G4PVPlacement(coil4rm,            // rotation
01038                         -dipolePosition,    // position
01039                         coil1LV,            // logical volume
01040                         name + "_coil4_pv", // name
01041                         containerLV,        // mother volume
01042                         false,              // boolean operation
01043                         0,                  // copy number
01044                         checkOverlaps);
01045 
01046       // pole placements
01047       new G4PVPlacement(0,                  // rotation
01048                         -dipolePosition,    // position
01049                         pole1LV,            // logical volume
01050                         name + "_pole1_pv", // name
01051                         containerLV,        // mother volume
01052                         false,              // boolean operation
01053                         0,                  // copy number
01054                         checkOverlaps);
01055       new G4PVPlacement(coil2rm,            // rotation
01056                         -dipolePosition,    // position
01057                         pole1LV,            // logical volume
01058                         name + "_pole2_pv", // name
01059                         containerLV,        // mother volume
01060                         false,              // boolean operation
01061                         0,                  // copy number
01062                         checkOverlaps);
01063       new G4PVPlacement(coil3rm,            // rotation
01064                         -dipolePosition,    // position
01065                         pole1LV,            // logical volume
01066                         name + "_pole3_pv", // name
01067                         containerLV,        // mother volume
01068                         false,              // boolean operation
01069                         0,                  // copy number
01070                         checkOverlaps);
01071       new G4PVPlacement(coil4rm,            // rotation
01072                         -dipolePosition,    // position
01073                         pole1LV,            // logical volume
01074                         name + "_pole4_pv", // name
01075                         containerLV,        // mother volume
01076                         false,              // boolean operation
01077                         0,                  // copy number
01078                         checkOverlaps);
01079     }
01080   
01081   // fixed coil
01082   coil2   = new G4Tubs(name+"_coil2_solid",          // name
01083                        coilInnerRadiusF,             // inner radius
01084                        coilOuterRadius,              // outer radius
01085                        length*0.5-2*lengthSafety,    // length
01086                        coilStartAngle,               // start angle
01087                        coilFullAngle);               // sweep angle
01088   coil2LV = new G4LogicalVolume(coil2,               // solid
01089                                 nbti,                // material
01090                                 name+"_coil2_lv");   // name
01091   coil2LV->SetVisAttributes(coilVisAtt);
01092   allLogicalVolumes.push_back(coil2LV);
01093   
01094   // fixed pole
01095   pole2   = new G4Tubs(name+"_pole2_solid",          // name
01096                        coilInnerRadiusF,             // inner radius
01097                        coilOuterRadius,              // outer radius
01098                        length*0.5 - lengthSafety,    // z half length
01099                        poleStartAngle,               // start angle
01100                        poleFullAngle);               // sweep angle
01101   pole2LV = new G4LogicalVolume(pole2,               // solid
01102                                 stainlesssteel,      // material
01103                                 name+"_pole2_lv");   // name
01104   pole2LV->SetVisAttributes(collarVisAtt);
01105   allLogicalVolumes.push_back(pole2LV);
01106   
01107   // fixed coil placements
01108   new G4PVPlacement(0,                  // rotation
01109                     dipolePosition,    // position
01110                     coil2LV,            // logical volume
01111                     name + "_coil5_pv", // name
01112                     containerLV,        // mother volume
01113                     false,              // boolean operation
01114                     0,                  // copy number
01115                     checkOverlaps);
01116   new G4PVPlacement(coil2rm,            // rotation
01117                     dipolePosition,    // position
01118                     coil2LV,            // logical volume
01119                     name + "_coil6_pv", // name
01120                     containerLV,        // mother volume
01121                     false,              // boolean operation
01122                     0,                  // copy number
01123                     checkOverlaps);
01124   new G4PVPlacement(coil3rm,            // rotation
01125                     dipolePosition,    // position
01126                     coil2LV,            // logical volume
01127                     name + "_coil7_pv", // name
01128                     containerLV,        // mother volume
01129                     false,              // boolean operation
01130                     0,                  // copy number
01131                     checkOverlaps);
01132   new G4PVPlacement(coil4rm,            // rotation
01133                     dipolePosition,    // position
01134                     coil2LV,            // logical volume
01135                     name + "_coil8_pv", // name
01136                     containerLV,        // mother volume
01137                     false,              // boolean operation
01138                     0,                  // copy number
01139                     checkOverlaps);
01140   
01141   // fixed pole placements
01142   new G4PVPlacement(0,                  // rotation
01143                     dipolePosition,    // position
01144                     pole2LV,            // logical volume
01145                     name + "_pole5_pv", // name
01146                     containerLV,        // mother volume
01147                     false,              // boolean operation
01148                     0,                  // copy number
01149                     checkOverlaps);
01150   new G4PVPlacement(coil2rm,            // rotation
01151                     dipolePosition,    // position
01152                     pole2LV,            // logical volume
01153                     name + "_pole6_pv", // name
01154                     containerLV,        // mother volume
01155                     false,              // boolean operation
01156                     0,                  // copy number
01157                     checkOverlaps);
01158   new G4PVPlacement(coil3rm,            // rotation
01159                     dipolePosition,    // position
01160                     pole2LV,            // logical volume
01161                     name + "_pole7_pv", // name
01162                     containerLV,        // mother volume
01163                     false,              // boolean operation
01164                     0,                  // copy number
01165                     checkOverlaps);
01166   new G4PVPlacement(coil4rm,            // rotation
01167                     dipolePosition,    // position
01168                     pole2LV,            // logical volume
01169                     name + "_pole8_pv", // name
01170                     containerLV,        // mother volume
01171                     false,              // boolean operation
01172                     0,                  // copy number
01173                     checkOverlaps);
01174   
01175   // non-magnetic collars
01176   // collar pole solid
01177   collar = new G4Tubs(name+"_collar_solid",        // name
01178                       collarInnerRadiusF,          // inner radius
01179                       collarOuterRadius,           // outer radius
01180                       length*0.5 - 2*lengthSafety, // length
01181                       0,                           // start angle
01182                       CLHEP::twopi);               // sweep angle
01183   collars = collar;
01184   if (buildCollar)
01185     {
01186       if (collarInnerRadius == collarInnerRadiusF)
01187         {
01188           // can save a solid by doing the union with the existing one if they're the same
01189           collars = new G4UnionSolid(name + "_collars_solid",  // name
01190                                      collar,                   // solid 1
01191                                      collar,                   // solid 2
01192                                      0,                        // rotation
01193                                      -2*dipolePosition);       // translation
01194         }
01195       else
01196         {
01197           G4VSolid* collar2 = new G4Tubs(name+"_collar2_solid",      // name
01198                                          collarInnerRadius,          // inner radius
01199                                          collarOuterRadius,          // outer radius
01200                                          length*0.5-2*lengthSafety,  // length
01201                                          0,                          // starting angle
01202                                          CLHEP::twopi);              // angle of sweep
01203           collars = new G4UnionSolid(name + "_collars_solid", // name
01204                                      collar,                  // solid 1
01205                                      collar2,                 // solid 2
01206                                      0,                       // rotation
01207                                      -2*dipolePosition);      // translation
01208         }
01209     }
01210   
01211   collarsLV = new G4LogicalVolume(collars,
01212                                   stainlesssteel,
01213                                   name+"_collars_lv");
01214   collarsLV->SetVisAttributes(collarVisAtt);
01215   allLogicalVolumes.push_back(collarsLV); 
01216 
01217   new G4PVPlacement(0,                  // rotation
01218                     dipolePosition,     // position
01219                     collarsLV,          // its logical volume
01220                     name+"_collars_pv", // its name
01221                     containerLV,        // its mother  volume
01222                     false,              // no boolean operation
01223                     0,                  // copy number
01224                     checkOverlaps);
01225 
01226   // prepare a solid to cut a hole in the outer yoke volume (can just use one twice)
01227   // can't use the existing collar solids as they're not solid - need them to be solid
01228   G4VSolid* collarSubtractionCylinder = new G4Tubs(name+"_collar_subtraction_solid",  // name
01229                                                    0,                                 // inner radius
01230                                                    collarOuterRadius + lengthSafety,  // outer radius
01231                                                    length,                            // double length for unambiguous subtraction
01232                                                    0,                                 // starting angle
01233                                                    CLHEP::twopi);                     // sweep angle
01234 
01235   G4VSolid* collarSubtractionCylinders = new G4UnionSolid(name + "_collar_subtraction_cylinders", // name
01236                                                           collarSubtractionCylinder,              // solid1
01237                                                           collarSubtractionCylinder,              // solid2 (here = solid1)
01238                                                           0,                                      // rotation
01239                                                           2*dipolePosition);                      // translation
01240  
01241   // outer iron yoke
01242   G4VSolid* yokeCylinder = new G4Tubs(name+"_yoke_cylinder_solid",     // name
01243                                       0.,                              // inner radius
01244                                       yokeOuterRadius,                 // outer radius
01245                                       0.5*length-2*lengthSafety,       // length
01246                                       0,                               // starting angle
01247                                       CLHEP::twopi * CLHEP::rad);      // sweep angle
01248 
01249   G4VSolid* yoke = new G4SubtractionSolid(name+"_yoke_solid",             // name
01250                                           yokeCylinder,                   // from this
01251                                           collarSubtractionCylinders,     // subtract this
01252                                           0,
01253                                           -dipolePosition);               
01254   
01255   G4LogicalVolume* yokeLV = new G4LogicalVolume(yoke,
01256                                                 iron,
01257                                                 name+"_yoke_lv");
01258 
01259   // yoke visualisation
01260   G4VisAttributes* LHCred = new G4VisAttributes(G4Colour(1.0, 0., 0.));
01261   LHCred->SetForceLineSegmentsPerCircle(nSegmentsPerCircle);
01262   yokeLV->SetVisAttributes(LHCred);
01263   
01264   allLogicalVolumes.push_back(yokeLV); // register locally
01265 
01266   // yoke placement
01267   new G4PVPlacement((G4RotationMatrix*)0,         // no rotation
01268                     G4ThreeVector(0,0,0),         // position
01269                     yokeLV,                       // lv to be placed
01270                     name + "_yoke_pv",            // name
01271                     containerLV,                  // mother lv to be place in
01272                     false,                        // no boolean operation
01273                     0,                            // copy number
01274                     checkOverlaps);
01275 
01276   G4String defaultMaterialName = BDSGlobalConstants::Instance()->GetBeamPipeMaterialName();
01277   G4Material* beamPipeMaterial = BDSMaterials::Instance()->GetMaterial(defaultMaterialName);
01278   G4Material* vacuumMaterial   = BDSMaterials::Instance()->GetMaterial(BDSGlobalConstants::Instance()->GetVacuumMaterial());
01279 
01280   
01281   //use beampipe factories to create another beampipe (note no magnetic field for now...)
01282   BDSBeamPipe* secondBP = BDSBeamPipeFactory::Instance()->CreateBeamPipe(BDSBeamPipeType::lhcdetailed,
01283                                                                          name,
01284                                                                          length-2*lengthSafety,
01285                                                                          2.202*CLHEP::cm,   // aper1
01286                                                                          1.714*CLHEP::cm,   // aper2
01287                                                                          2.202*CLHEP::cm,   // aper3
01288                                                                          0,                 // aper4
01289                                                                          vacuumMaterial,    // vacuum material
01290                                                                          1*CLHEP::mm,       // beampipeThickness
01291                                                                          beamPipeMaterial); // beampipe material
01292   
01293   G4LogicalVolume* secondBPLV = secondBP->GetContainerLogicalVolume();
01294   allLogicalVolumes.push_back(secondBPLV);
01295   new G4PVPlacement((G4RotationMatrix*)0,         // no rotation
01296                     dipolePosition,               // position
01297                     secondBPLV,                   // lv to be placed
01298                     name + "_second_beampipe_pv", // name
01299                     containerLV,                  // mother lv to be place in
01300                     false,                        // no boolean operation
01301                     0,                            // copy number
01302                     checkOverlaps);
01303   
01304   // visual attributes for container
01305   if (BDSExecOptions::Instance()->GetVisDebug())
01306     {containerLV->SetVisAttributes(BDSGlobalConstants::Instance()->GetVisibleDebugVisAttr());}
01307   else
01308     {containerLV->SetVisAttributes(BDSGlobalConstants::Instance()->GetInvisibleVisAttr());}
01309 
01310   // USER LIMITS and SENSITIVITY for all components
01311 #ifndef NOUSERLIMITS
01312   if (!allLogicalVolumes.empty()) {
01313     G4UserLimits* userLimits = new G4UserLimits("outer_cuts");
01314     userLimits->SetMaxAllowedStep( length * maxStepFactor );
01315     userLimits->SetUserMinEkine(BDSGlobalConstants::Instance()->GetThresholdCutCharged());
01316     userLimits->SetUserMaxTime(BDSGlobalConstants::Instance()->GetMaxTime());
01317     
01318     for (std::vector<G4LogicalVolume*>::iterator i = allLogicalVolumes.begin(); i != allLogicalVolumes.end(); ++i)
01319       {(*i)->SetUserLimits(userLimits);}
01320   }
01321 #endif
01322     
01323   // record extents
01324   // container radius is the same for all methods as all cylindrical
01325   G4double containerRadius = yokeOuterRadius;
01326   // massShift defined at very beginning of this function
01327   std::pair<double,double> extX = std::make_pair(-containerRadius+massShift,containerRadius+massShift); 
01328   std::pair<double,double> extY = std::make_pair(-containerRadius,containerRadius);
01329   std::pair<double,double> extZ = std::make_pair(-length*0.5,length*0.5);
01330   
01331   // build the BDSGeometryComponent instance and return it
01332   BDSGeometryComponent* outer = new BDSGeometryComponent(containerSolid,
01333                                                          containerLV,
01334                                                          extX, extY, extZ,
01335                                                          dipolePosition);
01336   // REGISTER all lvs
01337   outer->RegisterLogicalVolumes(secondBP->GetAllLogicalVolumes());
01338   outer->RegisterLogicalVolumes(allLogicalVolumes);
01339   
01340   return outer;
01341 }
01342 
01343 BDSGeometryComponent* BDSMagnetOuterFactoryLHC::CreateSextupole(G4String      name,
01344                                                                 G4double      length,
01345                                                                 BDSBeamPipe*  beamPipe,
01346                                                                 G4double      boxSize,
01347                                                                 G4Material*   outerMaterial)
01348 {
01349   CreateCylindricalSolids(name, length, beamPipe, boxSize);
01350   return CommonFinalConstructor(name, length, boxSize, outerMaterial, BDSMagnetColours::Instance()->GetMagnetColour("sextupole"));
01351 }
01352 
01353 BDSGeometryComponent* BDSMagnetOuterFactoryLHC::CreateOctupole(G4String      name,
01354                                                                G4double      length,
01355                                                                BDSBeamPipe*  beamPipe,
01356                                                                G4double      boxSize,
01357                                                                G4Material*   outerMaterial)
01358 {
01359 #ifdef BDSDEBUG
01360   G4cout << __METHOD_NAME__ << G4endl;
01361 #endif
01362   CreateCylindricalSolids(name, length, beamPipe, boxSize);
01363   return CommonFinalConstructor(name, length, boxSize, outerMaterial, BDSMagnetColours::Instance()->GetMagnetColour("octupole"));
01364 }
01365 
01366 BDSGeometryComponent* BDSMagnetOuterFactoryLHC::CreateDecapole(G4String      name,
01367                                                                G4double      length,
01368                                                                BDSBeamPipe*  beamPipe,
01369                                                                G4double      boxSize,
01370                                                                G4Material*   outerMaterial)
01371 {
01372 #ifdef BDSDEBUG
01373   G4cout << __METHOD_NAME__ << G4endl;
01374 #endif
01375   CreateCylindricalSolids(name, length, beamPipe, boxSize);
01376   return CommonFinalConstructor(name, length, boxSize, outerMaterial, BDSMagnetColours::Instance()->GetMagnetColour("decapole"));
01377 }
01378 
01379 BDSGeometryComponent* BDSMagnetOuterFactoryLHC::CreateSolenoid(G4String      name,
01380                                                                G4double      length,
01381                                                                BDSBeamPipe*  beamPipe,
01382                                                                G4double      boxSize,
01383                                                                G4Material*   outerMaterial)
01384 {
01385   CreateCylindricalSolids(name, length, beamPipe, boxSize);
01386   return CommonFinalConstructor(name, length, boxSize, outerMaterial, BDSMagnetColours::Instance()->GetMagnetColour("solenoid"));
01387 }
01388 
01389 BDSGeometryComponent* BDSMagnetOuterFactoryLHC::CreateMultipole(G4String      name,
01390                                                                 G4double      length,
01391                                                                 BDSBeamPipe*  beamPipe,
01392                                                                 G4double      boxSize,
01393                                                                 G4Material*   outerMaterial)
01394 {
01395 #ifdef BDSDEBUG
01396   G4cout << __METHOD_NAME__ << G4endl;
01397 #endif
01398   CreateCylindricalSolids(name, length, beamPipe, boxSize);
01399   return CommonFinalConstructor(name, length, boxSize, outerMaterial, BDSMagnetColours::Instance()->GetMagnetColour("multipole"));
01400 }
01401 
01402 BDSGeometryComponent* BDSMagnetOuterFactoryLHC::CreateRfCavity(G4String      name,
01403                                                                G4double      length,
01404                                                                BDSBeamPipe*  beamPipe,
01405                                                                G4double      boxSize,
01406                                                                G4Material*   outerMaterial)
01407 {
01408 #ifdef BDSDEBUG
01409   G4cout << __METHOD_NAME__ << G4endl;
01410 #endif
01411   CreateCylindricalSolids(name, length, beamPipe, boxSize);
01412   return CommonFinalConstructor(name, length, boxSize, outerMaterial, BDSMagnetColours::Instance()->GetMagnetColour("rfcavity"));
01413 }
01414 
01415 BDSGeometryComponent* BDSMagnetOuterFactoryLHC::CreateMuSpoiler(G4String      name,
01416                                                                 G4double      length,
01417                                                                 BDSBeamPipe*  beamPipe,
01418                                                                 G4double      boxSize,
01419                                                                 G4Material*   outerMaterial)
01420 {
01421 #ifdef BDSDEBUG
01422   G4cout << __METHOD_NAME__ << G4endl;
01423 #endif
01424   CreateCylindricalSolids(name, length, beamPipe, boxSize);
01425   return CommonFinalConstructor(name, length, boxSize, outerMaterial, BDSMagnetColours::Instance()->GetMagnetColour("muspoiler"));
01426 }
01427 
01428 BDSGeometryComponent* BDSMagnetOuterFactoryLHC::CreateKicker(G4String      name,
01429                                                              G4double      length,
01430                                                              BDSBeamPipe*  beamPipe,
01431                                                              G4double      boxSize,
01432                                                              G4bool        /*vertical*/,
01433                                                              G4Material*   outerMaterial)
01434 {
01435 #ifdef BDSDEBUG
01436   G4cout << __METHOD_NAME__ << G4endl;
01437 #endif
01438   // in this factory, h and v kickers will look the same so ignore bool vertical
01439   // have to retain it though for virtual base class compatability
01440   CreateCylindricalSolids(name, length, beamPipe, boxSize);
01441   return CommonFinalConstructor(name, length, boxSize, outerMaterial, BDSMagnetColours::Instance()->GetMagnetColour("hkicker"));
01442 }
01443 
01445 
01446 void BDSMagnetOuterFactoryLHC::CreateCylindricalSolids(G4String     name,
01447                                                        G4double     length,
01448                                                        BDSBeamPipe* beamPipe,
01449                                                        G4double     boxSize)
01450 {
01451   if (beamPipe->ContainerIsCircular())
01452     {
01453       //circular beampipe so we can simply use its radius
01454       yokeSolid = new G4Tubs(name + "_yoke_solid",       // name
01455                              beamPipe->GetContainerRadius() + 2*lengthSafety,  // inner radius
01456                              boxSize*0.5,                 // outer radius
01457                              length*0.5-2*lengthSafety,   // half length
01458                              0,                           // rotation start angle
01459                              CLHEP::twopi);               // rotation finish angle
01460 
01461       //container is similar but slightly wider and hollow (to allow placement of beampipe)
01462       containerSolid = new G4Tubs(name + "_contiainer_solid",  // name
01463                                   beamPipe->GetContainerRadius() + lengthSafety, // inner radius
01464                                   boxSize*0.5 + lengthSafety,  // outer radius
01465                                   length*0.5,                  // half length
01466                                   0,                           // rotation start angle
01467                                   CLHEP::twopi);               // rotation finish angle
01468     }
01469   else
01470     {
01471       G4VSolid* yokeSolidCylinder = new G4Tubs(name + "_yoke_solid_cylinder",  // name
01472                                                0,  // inner radius - for unambiguous subtraction
01473                                                boxSize*0.5,                 // outer radius
01474                                                length*0.5-2*lengthSafety,   // half length
01475                                                0,                           // rotation start angle
01476                                                CLHEP::twopi);               // rotation finish angle
01477       yokeSolid = new G4SubtractionSolid(name + "_yoke_solid",
01478                                          yokeSolidCylinder,
01479                                          beamPipe->GetContainerSubtractionSolid());
01480       
01481       //container is similar but slightly wider
01482       G4VSolid* containerSolidCylinder = new G4Tubs(name + "_container_solid_cylinder", // name
01483                                                     0,  // inner radius - for unambiguous subtraction
01484                                                     boxSize*0.5 + lengthSafety,  // outer radius
01485                                                     length*0.5,                  // half length
01486                                                     0,                           // rotation start angle
01487                                                     CLHEP::twopi);               // rotation finish angle
01488       containerSolid = new G4SubtractionSolid(name + "_container_solid",
01489                                               containerSolidCylinder,
01490                                               beamPipe->GetContainerSubtractionSolid());
01491     }
01492 }
01493 
01494 void BDSMagnetOuterFactoryLHC::TestInputParameters(BDSBeamPipe* /*beamPipe*/,
01495                                                    G4double&    outerDiameter,
01496                                                    G4Material*& outerMaterial)// reference to a pointer
01497 {
01498   //function arguments by reference to they can be modified in place
01499   //check outer material is something
01500   if (!outerMaterial)
01501     {outerMaterial = BDSMaterials::Instance()->GetMaterial("stainlesssteel");}
01502 
01503   // ensure outerDiameter is > outerCollarDiameter - hard coded as specific to the lhc design
01504   if (outerDiameter < 202*CLHEP::mm )
01505     {outerDiameter = 202*CLHEP::mm;}
01506 }
01507 
01510 BDSGeometryComponent* BDSMagnetOuterFactoryLHC::CommonFinalConstructor(G4String    name,
01511                                                                        G4double    length,
01512                                                                        G4double    boxSize,
01513                                                                        G4Material* outerMaterial,
01514                                                                        G4Colour*   colour)
01515 {
01516 #ifdef BDSDEBUG
01517   G4cout << __METHOD_NAME__ << G4endl;
01518 #endif
01519   
01520   // build the logical volumes
01521   G4LogicalVolume* yokeLV   = new G4LogicalVolume(yokeSolid,
01522                                                   outerMaterial,
01523                                                   name + "_yoke_lv");
01524 
01525   G4Material* emptyMaterial = BDSMaterials::Instance()->GetMaterial(BDSGlobalConstants::Instance()->GetEmptyMaterial());
01526   G4LogicalVolume* containerLV = new G4LogicalVolume(containerSolid,
01527                                                      emptyMaterial,
01528                                                      name + "_container_lv");
01529   
01530   // VISUAL ATTRIBUTES
01531   // set visual attributes
01532   // outer
01533   G4VisAttributes* outerVisAttr = new G4VisAttributes(*colour);
01534   outerVisAttr->SetVisibility(true);
01535   yokeLV->SetVisAttributes(outerVisAttr);
01536   // container
01537   // visual attributes for container
01538   if (BDSExecOptions::Instance()->GetVisDebug())
01539     {containerLV->SetVisAttributes(BDSGlobalConstants::Instance()->GetVisibleDebugVisAttr());}
01540   else
01541     {containerLV->SetVisAttributes(BDSGlobalConstants::Instance()->GetInvisibleVisAttr());}
01542 
01543   // USER LIMITS - set user limits based on bdsim user specified parameters
01544 #ifndef NOUSERLIMITS
01545   G4UserLimits* outerUserLimits = new G4UserLimits("outer_cuts");
01546   outerUserLimits->SetMaxAllowedStep( length * maxStepFactor );
01547   outerUserLimits->SetUserMinEkine(BDSGlobalConstants::Instance()->GetThresholdCutCharged());
01548   outerUserLimits->SetUserMaxTime(BDSGlobalConstants::Instance()->GetMaxTime());
01549   //attach cuts to volumes
01550   yokeLV->SetUserLimits(outerUserLimits);
01551   containerLV->SetUserLimits(outerUserLimits);
01552 #endif
01553 
01554   // PLACEMENT
01555   // place the components inside the container
01556   // note we don't need the pointer for anything - it's registered upon construction with g4
01557   
01558   new G4PVPlacement((G4RotationMatrix*)0,         // no rotation
01559                     (G4ThreeVector)0,             // position
01560                     yokeLV,                       // lv to be placed
01561                     name + "_outer_pv",           // name
01562                     containerLV,                  // mother lv to be place in
01563                     false,                        // no boolean operation
01564                     0,                            // copy number
01565                     checkOverlaps);               // whether to check overlaps
01566   
01567   // record extents
01568   // container radius is the same for all methods as all cylindrical
01569   G4double containerRadius = boxSize + lengthSafety;
01570   std::pair<double,double> extX = std::make_pair(-containerRadius,containerRadius);
01571   std::pair<double,double> extY = std::make_pair(-containerRadius,containerRadius);
01572   std::pair<double,double> extZ = std::make_pair(-length*0.5,length*0.5);
01573   
01574   // build the BDSGeometryComponent instance and return it
01575   BDSGeometryComponent* outer = new BDSGeometryComponent(containerSolid,
01576                                                          containerLV,
01577                                                          extX, extY, extZ);
01578   // REGISTER all lvs
01579   outer->RegisterLogicalVolume(yokeLV); //using geometry component base class method
01580   
01581   return outer;
01582 }

Generated on 28 Jun 2015 for BDSIM by  doxygen 1.4.7