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

00001 /* BDSIM code.    Version 1.0
00002    Author: Grahame A. Blair, Royal Holloway, Univ. of London.
00003    Last modified 24.7.2002
00004    Copyright (c) 2002 by G.A.Blair.  ALL RIGHTS RESERVED. 
00005 
00006    Modified 22.03.05 by J.C.Carter, Royal Holloway, Univ. of London.
00007    Added extra parameter to BuildLogicalVolume so that it is 
00008      possible to set the material as either Iron or Vacuum
00009    Removed StringFromInt function
00010 */
00011 
00012 #include "BDSExecOptions.hh"
00013 #include "BDSGlobalConstants.hh" 
00014 #include "BDSDebug.hh"
00015 
00016 #include <cstdlib>
00017 #include <cstddef>
00018 #include <cmath>
00019 #include <string>
00020 #include <algorithm> // for std::max
00021 
00022 #include "G4Box.hh"
00023 #include "G4CutTubs.hh"
00024 #include "G4LogicalVolume.hh"
00025 #include "G4MagIntegratorStepper.hh"
00026 #include "G4MagneticField.hh"
00027 #include "G4PVPlacement.hh"
00028 #include "G4UserLimits.hh"
00029 #include "G4VisAttributes.hh"
00030 #include "G4VPhysicalVolume.hh"
00031 
00032 #include "BDSBeamPipe.hh"
00033 #include "BDSBeamPipeFactory.hh"
00034 #include "BDSBeamPipeType.hh"
00035 #include "BDSBeamPipeInfo.hh"
00036 #include "BDSDebug.hh"
00037 #include "BDSGlobalConstants.hh"
00038 #include "BDSMaterials.hh"
00039 #include "BDSMagnetOuterFactory.hh"
00040 #include "BDSMagnetType.hh"
00041 #include "BDSMagnet.hh"
00042 #include "BDSMultipoleOuterMagField.hh"
00043 #include "BDSUtilities.hh"
00044 
00045 class BDSTiltOffset;
00046 
00047 BDSMagnet::BDSMagnet(BDSMagnetType      type,
00048                      G4String           name,
00049                      G4double           length,
00050                      BDSBeamPipeInfo*   beamPipeInfoIn,
00051                      BDSMagnetOuterInfo magnetOuterInfo,
00052                      BDSTiltOffset      tiltOffset):
00053   BDSAcceleratorComponent(name, length, 0, "magnet", tiltOffset),
00054   itsType(type),
00055   beamPipeInfo(beamPipeInfoIn),
00056   outerDiameter(magnetOuterInfo.outerDiameter),
00057   itsMagnetOuterInfo(magnetOuterInfo)
00058 {
00059   inputface  = G4ThreeVector(0,0,0);
00060   outputface = G4ThreeVector(0,0,0);
00061   itsK1 = 0.0; itsK2 = 0.0; itsK3 = 0.0;
00062   itsStepper=NULL;
00063   itsMagField=NULL;
00064   itsEqRhs=NULL;
00065 
00066   itsBeampipeLogicalVolume=NULL;
00067   itsInnerBPLogicalVolume=NULL;
00068 
00069   itsBeampipeUserLimits=NULL;
00070   itsPhysiComp=NULL; 
00071   itsPhysiInner=NULL;
00072   itsBPFieldMgr=NULL;
00073   itsOuterFieldMgr=NULL;
00074   
00075   itsChordFinder=NULL;
00076   itsOuterMagField=NULL;
00077 
00078   beampipe = NULL;
00079   outer    = NULL;
00080 }
00081 
00082 void BDSMagnet::Build()
00083 {
00084 #ifdef BDSDEBUG
00085   G4cout << __METHOD_NAME__ << G4endl;
00086 #endif
00087   BuildBPFieldAndStepper();
00088   BuildBPFieldMgr(itsStepper, itsMagField);
00089 
00090   BDSAcceleratorComponent::Build(); //builds marker logical volume & itVisAttributes
00091   BuildBeampipe();
00092   BuildOuterVolume();
00093 }
00094 
00095 void BDSMagnet::BuildBeampipe()
00096 {
00097 #ifdef BDSDEBUG
00098   G4cout << __METHOD_NAME__ << G4endl;
00099 #endif
00100   
00101   beampipe = BDSBeamPipeFactory::Instance()->CreateBeamPipe(name,
00102                                                             chordLength,
00103                                                             beamPipeInfo);
00104   BeamPipeCommonTasks();
00105 }
00106 
00107 void BDSMagnet::BeamPipeCommonTasks()
00108 {
00109   //actual beam pipe ->SetFieldManager(BDSGlobalConstants::Instance()->GetZeroFieldManager(),false);
00110   // SET FIELD
00111   if(itsBPFieldMgr)
00112     {beampipe->GetVacuumLogicalVolume()->SetFieldManager(itsBPFieldMgr,false);}
00113 
00114   // now protect the fields inside the marker volume by giving the
00115   // marker a null magnetic field (otherwise G4VPlacement can
00116   // over-ride the already-created fields, by calling 
00117   // G4LogicalVolume::AddDaughter, which calls 
00118   // pDaughterLogical->SetFieldManager(fFieldManager, true) - the
00119   // latter 'true' over-writes all the other fields
00120   
00121   containerLogicalVolume->
00122     SetFieldManager(BDSGlobalConstants::Instance()->GetZeroFieldManager(),false);
00123   
00124 
00125   // register logical volumes using geometry component base class
00126   RegisterLogicalVolumes(beampipe->GetAllLogicalVolumes());
00127 
00128   // register components as sensitive if required
00129   if(BDSGlobalConstants::Instance()->GetSensitiveBeamPipe())
00130     {RegisterSensitiveVolumes(beampipe->GetAllSensitiveVolumes());}
00131 
00132   // place beampipe
00133   itsPhysiComp = new G4PVPlacement(0,                         // rotation
00134                                    (G4ThreeVector)0,          // at (0,0,0)
00135                                    beampipe->GetContainerLogicalVolume(),  // its logical volume
00136                                    name + "_beampipe_pv",     // its name
00137                                    containerLogicalVolume,    // its mother  volume
00138                                    false,                     // no boolean operation
00139                                    0,                         // copy number
00140                                    BDSGlobalConstants::Instance()->GetCheckOverlaps());
00141 }
00142 
00143 void BDSMagnet::BuildBPFieldMgr(G4MagIntegratorStepper* aStepper,
00144                                    G4MagneticField* aField)
00145 {
00146   itsChordFinder= 
00147     new G4ChordFinder(aField,
00148                       BDSGlobalConstants::Instance()->GetChordStepMinimum(),
00149                       aStepper);
00150 
00151   itsChordFinder->SetDeltaChord(BDSGlobalConstants::Instance()->GetDeltaChord());
00152   itsBPFieldMgr= new G4FieldManager();
00153   itsBPFieldMgr->SetDetectorField(aField);
00154   itsBPFieldMgr->SetChordFinder(itsChordFinder);
00155   if(BDSGlobalConstants::Instance()->GetDeltaIntersection()>0){
00156     itsBPFieldMgr->SetDeltaIntersection(BDSGlobalConstants::Instance()->GetDeltaIntersection());
00157   }
00158   if(BDSGlobalConstants::Instance()->GetMinimumEpsilonStep()>0)
00159     itsBPFieldMgr->SetMinimumEpsilonStep(BDSGlobalConstants::Instance()->GetMinimumEpsilonStep());
00160   if(BDSGlobalConstants::Instance()->GetMaximumEpsilonStep()>0)
00161     itsBPFieldMgr->SetMaximumEpsilonStep(BDSGlobalConstants::Instance()->GetMaximumEpsilonStep());
00162   if(BDSGlobalConstants::Instance()->GetDeltaOneStep()>0)
00163     itsBPFieldMgr->SetDeltaOneStep(BDSGlobalConstants::Instance()->GetDeltaOneStep());
00164 }
00165 
00166 void BDSMagnet::BuildOuterVolume()
00167 {
00168   G4Material* outerMaterial          = itsMagnetOuterInfo.outerMaterial;
00169   BDSMagnetGeometryType geometryType = itsMagnetOuterInfo.geometryType;
00170   G4double outerDiameter             = itsMagnetOuterInfo.outerDiameter;
00171 
00172   //build the right thing depending on the magnet type
00173   //saves basically the same funciton in each derived class
00174   BDSMagnetOuterFactory* theFactory = BDSMagnetOuterFactory::Instance();
00175   switch(itsType.underlying()){
00176   case BDSMagnetType::decapole:
00177     outer = theFactory->CreateDecapole(geometryType,name,chordLength,beampipe,
00178                                        outerDiameter,outerMaterial);
00179     break;
00180   case BDSMagnetType::vkicker:
00181     outer = theFactory->CreateKicker(geometryType,name,chordLength,beampipe,
00182                                      outerDiameter,true,outerMaterial);
00183     break;
00184   case BDSMagnetType::hkicker:
00185     outer = theFactory->CreateKicker(geometryType,name,chordLength,beampipe,
00186                                      outerDiameter,false,outerMaterial);
00187     break;
00188   case BDSMagnetType::muspoiler:
00189     outer = theFactory->CreateMuSpoiler(geometryType,name,chordLength,beampipe,
00190                                         outerDiameter,outerMaterial);
00191     break;
00192   case BDSMagnetType::octupole:
00193     outer = theFactory->CreateOctupole(geometryType,name,chordLength,beampipe,
00194                                        outerDiameter,outerMaterial);
00195     break;
00196   case BDSMagnetType::quadrupole:
00197     outer = theFactory->CreateQuadrupole(geometryType,name,chordLength,beampipe,
00198                                          outerDiameter,outerMaterial);
00199     break;
00200   case BDSMagnetType::rectangularbend:
00201     outer = theFactory->CreateRectangularBend(geometryType,name,chordLength,beampipe,
00202                                               outerDiameter,angle,outerMaterial);
00203     break;
00204   case BDSMagnetType::rfcavity:
00205     outer = theFactory->CreateRfCavity(geometryType,name,chordLength,beampipe,
00206                                        outerDiameter,outerMaterial);
00207     break;
00208   case BDSMagnetType::sectorbend:
00209     outer = theFactory->CreateSectorBend(geometryType,name,chordLength,beampipe,
00210                                          outerDiameter,angle,outerMaterial);
00211     break;
00212   case BDSMagnetType::sextupole:
00213     outer = theFactory->CreateSextupole(geometryType,name,chordLength,beampipe,
00214                                         outerDiameter,outerMaterial);
00215     break;
00216   case BDSMagnetType::solenoid:
00217     outer = theFactory->CreateSolenoid(geometryType,name,chordLength,beampipe,
00218                                        outerDiameter,outerMaterial);
00219     break;
00220   case BDSMagnetType::multipole:
00221     outer = theFactory->CreateMultipole(geometryType,name,chordLength,beampipe,
00222                                         outerDiameter,outerMaterial);
00223     break;
00224   default:
00225     G4cout << __METHOD_NAME__ << "unknown magnet type - no outer volume built" << G4endl;
00226     outer = NULL;
00227     break;
00228   }
00229 
00230   if(outer)
00231     {
00232       // register logical volumes using geometry component base class
00233       RegisterLogicalVolumes(outer->GetAllLogicalVolumes());
00234 
00235       // register components as sensitive if required
00236       if(BDSGlobalConstants::Instance()->GetSensitiveComponents())
00237         {RegisterSensitiveVolumes(outer->GetAllLogicalVolumes());}
00238       
00239       // place outer volume
00240       itsPhysiComp = new G4PVPlacement(0,                           // rotation
00241                                        outer->GetPlacementOffset(), // at normally (0,0,0)
00242                                        outer->GetContainerLogicalVolume(), // its logical volume
00243                                        name+"_outer_phys",          // its name
00244                                        containerLogicalVolume,      // its mother  volume
00245                                        false,                       // no boolean operation
00246                                        0,                           // copy number
00247                                        BDSGlobalConstants::Instance()->GetCheckOverlaps());
00248 
00249       //update extents
00250       SetExtentX(outer->GetExtentX());
00251       SetExtentY(outer->GetExtentY());
00252       SetExtentZ(outer->GetExtentZ());
00253     }
00254 }
00255 
00256 void BDSMagnet::BuildContainerLogicalVolume()
00257 {
00258 #ifdef BDSDEBUG
00259   G4cout << __METHOD_NAME__ << G4endl;
00260 #endif
00261   if (BDS::IsFinite(angle))
00262     {
00263       containerSolid = new G4CutTubs(name+"_container_solid",          // name
00264                                      0.0,                              // minimum radius - solid here
00265                                      outerDiameter*0.5 + lengthSafety, // radius - determined above
00266                                      chordLength*0.5,                  // length about centre point
00267                                      0.0,                              // starting angle
00268                                      2.0*CLHEP::pi,                    // sweep angle
00269                                      inputface,                        // input face normal vector
00270                                      outputface);                      // output face normal vector
00271     }
00272   else
00273     {
00274       containerSolid = new G4Box(name + "_container_solid",
00275                                  outerDiameter*0.5 + lengthSafety,
00276                                  outerDiameter*0.5 + lengthSafety,
00277                                  chordLength/2);
00278     }
00279     
00280   containerLogicalVolume = new G4LogicalVolume(containerSolid,
00281                                                emptyMaterial,
00282                                                name + "_container_lv");
00283 
00284   // taken from FinaliseBeamPipe method - supposed to protect against fields being overridden
00285   containerLogicalVolume->
00286     SetFieldManager(BDSGlobalConstants::Instance()->GetZeroFieldManager(),false);
00287 
00288   // USER LIMITS
00289 #ifndef NOUSERLIMITS
00290   G4double maxStepFactor=0.5;
00291   G4UserLimits* itsMarkerUserLimits =  new G4UserLimits();
00292   itsMarkerUserLimits->SetMaxAllowedStep(chordLength*maxStepFactor);
00293   itsMarkerUserLimits->SetUserMinEkine(BDSGlobalConstants::Instance()->GetThresholdCutCharged());
00294   containerLogicalVolume->SetUserLimits(itsMarkerUserLimits);
00295 #endif
00296 
00297   // VIS ATTR
00298   if (BDSExecOptions::Instance()->GetVisDebug())
00299     {containerLogicalVolume->SetVisAttributes(BDSGlobalConstants::Instance()->GetVisibleDebugVisAttr());}
00300   else
00301     {containerLogicalVolume->SetVisAttributes(BDSGlobalConstants::Instance()->GetInvisibleVisAttr());}
00302 
00303   G4double transverseSize = outerDiameter*0.5 + lengthSafety;
00304   SetExtentX(-transverseSize, transverseSize);
00305   SetExtentY(-transverseSize, transverseSize);
00306   SetExtentZ(-chordLength*0.5, chordLength*0.5);
00307 }
00308 
00309 void BDSMagnet::BuildOuterFieldManager(G4int nPoles, G4double poleField,
00310                                           G4double phiOffset)
00311 {
00312   if(nPoles<=0 || nPoles>10 || nPoles%2 !=0)
00313     G4Exception("BDSMagnet: Invalid number of poles", "-1", FatalException, "");
00314   itsOuterFieldMgr=NULL;
00315   if (poleField==0) return;
00316 
00317   itsOuterMagField=new BDSMultipoleOuterMagField(nPoles,poleField,phiOffset);
00318 
00319   itsOuterFieldMgr=new G4FieldManager(itsOuterMagField);
00320   if(BDSGlobalConstants::Instance()->GetDeltaIntersection()>0){
00321     itsOuterFieldMgr->SetDeltaIntersection(BDSGlobalConstants::Instance()->GetDeltaIntersection());
00322   }
00323   if(BDSGlobalConstants::Instance()->GetMinimumEpsilonStep()>0)
00324     itsOuterFieldMgr->SetMinimumEpsilonStep(BDSGlobalConstants::Instance()->GetMinimumEpsilonStep());
00325   if(BDSGlobalConstants::Instance()->GetMaximumEpsilonStep()>0)
00326     itsOuterFieldMgr->SetMaximumEpsilonStep(BDSGlobalConstants::Instance()->GetMaximumEpsilonStep());
00327   if(BDSGlobalConstants::Instance()->GetDeltaOneStep()>0)
00328     itsOuterFieldMgr->SetDeltaOneStep(BDSGlobalConstants::Instance()->GetDeltaOneStep());
00329   outer->GetContainerLogicalVolume()->SetFieldManager(itsOuterFieldMgr,false);
00330 }
00331 
00332 BDSMagnet::~BDSMagnet()
00333 {
00334   delete itsBPFieldMgr;
00335   delete itsChordFinder;
00336 #ifndef NOUSERLIMITS
00337   delete itsBeampipeUserLimits;
00338 #endif
00339   delete itsMagField;
00340   delete itsEqRhs;
00341   delete itsStepper;
00342 }

Generated on 28 Jun 2015 for BDSIM by  doxygen 1.4.7