src/BDSPCLDrift.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 #include "BDSGlobalConstants.hh" 
00007 
00008 #include "BDSPCLDrift.hh"
00009 #include "BDSMagField.hh"
00010 #include "BDSDriftStepper.hh"
00011 #include "BDSPCLTube.hh"
00012 
00013 #include "G4Box.hh"
00014 #include "G4Tubs.hh"
00015 #include "G4VisAttributes.hh"
00016 #include "G4LogicalVolume.hh"
00017 #include "G4VPhysicalVolume.hh"
00018 #include "G4UserLimits.hh"
00019 #include "G4TransportationManager.hh"
00020 #include "G4CashKarpRKF45.hh"
00021 #include "G4ExactHelixStepper.hh"
00022 #include "G4UnionSolid.hh"
00023 
00024 #include <map>
00025 
00026 //============================================================
00027 
00028 typedef std::map<G4String,int> LogVolCountMap;
00029 extern LogVolCountMap* LogVolCount;
00030 
00031 typedef std::map<G4String,G4LogicalVolume*> LogVolMap;
00032 extern LogVolMap* LogVol;
00033 
00034 //============================================================
00035 
00036 BDSPCLDrift::BDSPCLDrift (G4String aName, G4double aLength, 
00037                           std::list<G4double> blmLocZ, std::list<G4double> blmLocTheta, G4double aperX, G4double aperYUp, G4double aperYDown, G4double aperDy, G4String tunnelMaterial, G4double aper, G4double tunnelRadius, G4double tunnelOffsetX, G4String aMaterial):
00038   BDSMultipole(aName, aLength, aper, aper, SetVisAttributes(),  blmLocZ, blmLocTheta, tunnelMaterial, aMaterial, aper, aper, 0, tunnelRadius, tunnelOffsetX),
00039   itsYAperUp(aperYUp), itsYAperDown(aperYDown), itsDyAper(aperDy),
00040   outer_solid(NULL),inner_solid(NULL),itsOuterBeamPipeLogicalVolume(NULL),
00041   itsInnerBeamPipeLogicalVolume(NULL),itsPhysiInner(NULL),itsPhysiOuter(NULL),
00042   itsBeampipeVisAtt(NULL),itsInnerBeampipeVisAtt(NULL),
00043   itsStepper(NULL),itsMagField(NULL),itsEqRhs(NULL)
00044 {
00045   itsType="pcldrift";
00046   itsXAper=aperX;
00047 
00048   if (!(*LogVolCount)[itsName])
00049     {
00050       //
00051       // build external volume
00052       // 
00053       BuildDefaultMarkerLogicalVolume();
00054 
00055       //
00056       // build beampipe (geometry + magnetic field)
00057       //
00058       if(BDSGlobalConstants::Instance()->GetBuildTunnel()){
00059         BuildTunnel();
00060       }
00061       BuildBpFieldAndStepper();
00062       BuildBPFieldMgr(itsStepper, itsMagField);
00063       BuildBeampipe(itsMaterial);
00064       BuildBLMs();
00065   
00066       //
00067       // define sensitive volumes for hit generation
00068       //
00069       if(BDSGlobalConstants::Instance()->GetSensitiveBeamPipe()){
00070         SetMultipleSensitiveVolumes(itsOuterBeamPipeLogicalVolume);     
00071       }
00072       
00073       //
00074       // append marker logical volume to volume map
00075       //
00076       (*LogVolCount)[itsName]=1;
00077       (*LogVol)[itsName]=itsMarkerLogicalVolume;
00078     }
00079   else
00080     {
00081       (*LogVolCount)[itsName]++;
00082       
00083       //
00084       // use already defined marker volume
00085       //
00086       itsMarkerLogicalVolume=(*LogVol)[itsName];
00087     }
00088 }
00089 
00090 void BDSPCLDrift::BuildBeampipe(G4String materialName){
00091   G4Material *material;
00092   if(materialName != ""){
00093     material = BDSMaterials::Instance()->GetMaterial( materialName );
00094   } else {
00095     material = BDSMaterials::Instance()->GetMaterial( BDSGlobalConstants::Instance()->GetPipeMaterialName());
00096   }
00097   
00098   // build beampipe
00099 
00100 #ifdef DEBUG 
00101   G4cout << "Outer pipe :"
00102          << " r= " << itsBpRadius/m << " m"
00103          << " l= " << itsLength/(2.)/m << " m"
00104          << G4endl;
00105   G4cout << "PCLDrift aperX: " << itsXAper/m << " m" << G4endl;
00106   G4cout << "PCLDrift aperYUp: " << itsYAperUp/m << " m" << G4endl;
00107   G4cout << "PCLDrift aperYDown: " << itsYAperDown/m << " m" << G4endl;
00108   G4cout << "PCLDrift Dy: " << itsDyAper/m << " m" << G4endl;
00109 #endif
00110 
00111 
00112 
00113   G4double ts = BDSGlobalConstants::Instance()->GetLengthSafety()+itsBeampipeThickness/2;
00114 
00115   BDSPCLTube* innerTube = new BDSPCLTube(itsXAper-ts, itsYAperUp-ts, itsYAperDown-ts, itsDyAper, -1, itsLength, itsName+"_inner");
00116 
00117   BDSPCLTube* outerTube = new BDSPCLTube(itsXAper, itsYAperUp, itsYAperDown, itsDyAper, itsBeampipeThickness, itsLength, itsName+"_outer");
00118   
00119   //The inner beam pipe solids need to be fused together to form one seamless shape in order to create a vacuum and field without gaps in it.
00120 
00121   inner_solid = innerTube->GetSolid();
00122   outer_solid = outerTube->GetSolid();
00123 
00124   delete innerTube;
00125   delete outerTube;
00126 
00127 #ifdef DEBUG
00128   G4cout << "BDSPCLDrift.cc: Making logical..." << G4endl;
00129 #endif
00130 
00131   itsInnerBeamPipeLogicalVolume=        
00132     new G4LogicalVolume(inner_solid,
00133                         BDSMaterials::Instance()->GetMaterial(BDSGlobalConstants::Instance()->GetVacuumMaterial()),
00134                         itsName+"_inner_bmp_log");
00135 
00136   itsOuterBeamPipeLogicalVolume=        
00137     new G4LogicalVolume(outer_solid,
00138                         material,
00139                         itsName+"_outer_bmp_log");
00140   
00141   //
00142   // set visualization attributes
00143   //
00144   itsInnerBeampipeVisAtt =  new G4VisAttributes(G4Colour(0.5, 0.5, 0.5));
00145   itsInnerBeampipeVisAtt->SetVisibility(true);
00146   itsInnerBeampipeVisAtt->SetForceSolid(true);
00147 
00148   itsBeampipeVisAtt = new G4VisAttributes(G4Colour(0.4, 0.4, 0.4));
00149   itsBeampipeVisAtt->SetForceWireframe(true);
00150   itsBeampipeVisAtt->SetVisibility(false);
00151 
00152   itsOuterBeamPipeLogicalVolume->SetVisAttributes(itsBeampipeVisAtt);
00153   itsInnerBeamPipeLogicalVolume->SetVisAttributes(itsInnerBeampipeVisAtt);
00154 
00155 #ifdef DEBUG  
00156   G4cout << "BDSPCLDrift.cc: Placing..." << G4endl;
00157 #endif
00158 
00159   G4ThreeVector threeVector1;
00160   threeVector1.setY(0);
00161   itsPhysiInner = new G4PVPlacement(
00162                                     (G4RotationMatrix*)0,                       // no rotation
00163                                     threeVector1,                      
00164                                     itsInnerBeamPipeLogicalVolume,  // its logical volume
00165                                     itsName+"_inner_bmp_phys",// its name
00166                                     itsMarkerLogicalVolume,   // its mother  volume
00167                                     false,                      // no boolean operation
00168                                     0, BDSGlobalConstants::Instance()->GetCheckOverlaps());                     // copy number
00169 
00170 
00171   
00172   itsPhysiOuter = new G4PVPlacement(
00173                                     (G4RotationMatrix*)0,                       // no rotation
00174                                     threeVector1,                      
00175                                     itsOuterBeamPipeLogicalVolume,  // its logical volume
00176                                     itsName+"_inner_bmp_phys",// its name
00177                                     itsMarkerLogicalVolume,   // its mother  volume
00178                                     false,                      // no boolean operation
00179                                     0, BDSGlobalConstants::Instance()->GetCheckOverlaps());                     // copy number
00180 
00181 
00182 
00183 
00184   //Add the physical volumes to a vector which can be used for e.g. geometrical biasing
00185   SetMultiplePhysicalVolumes(itsPhysiInner);
00186   SetMultiplePhysicalVolumes(itsPhysiOuter);
00187 
00188 #ifndef NOUSERLIMITS
00189   itsBeampipeUserLimits =  new G4UserLimits("beampipe cuts");
00190   itsBeampipeUserLimits->SetUserMinEkine(BDSGlobalConstants::Instance()->GetThresholdCutCharged());
00191   itsBeampipeUserLimits->SetUserMaxTime(BDSGlobalConstants::Instance()->GetMaxTime());
00192 
00193   itsInnerBeampipeUserLimits =  new G4UserLimits("inner beampipe cuts");
00194   itsInnerBeampipeUserLimits->SetUserMinEkine(BDSGlobalConstants::Instance()->GetThresholdCutCharged());
00195   itsInnerBeampipeUserLimits->SetUserMaxTime(BDSGlobalConstants::Instance()->GetMaxTime());
00196 
00197   G4double stepfactor = 5;
00198 
00199   itsBeampipeUserLimits->SetMaxAllowedStep(itsLength*stepfactor);
00200   itsInnerBeampipeUserLimits->SetMaxAllowedStep(itsLength*stepfactor);
00201 
00202   itsInnerBeamPipeLogicalVolume->SetUserLimits(itsInnerBeampipeUserLimits);
00203   itsOuterBeamPipeLogicalVolume->SetUserLimits(itsBeampipeUserLimits);
00204 
00205 
00206 #endif
00207 
00208    itsInnerBeamPipeLogicalVolume->SetFieldManager(itsBPFieldMgr,false);
00209    itsOuterBeamPipeLogicalVolume->SetFieldManager(itsBPFieldMgr,false);
00210 
00211   
00212   // now protect the fields inside the marker volume by giving the
00213   // marker a null magnetic field (otherwise G4VPlacement can
00214   // over-ride the already-created fields, by calling 
00215   // G4LogicalVolume::AddDaughter, which calls 
00216   // pDaughterLogical->SetFieldManager(fFieldManager, true) - the
00217   // latter 'true' over-writes all the other fields
00218   
00219   itsMarkerLogicalVolume->
00220     SetFieldManager(BDSGlobalConstants::Instance()->GetZeroFieldManager(),false);
00221   
00222 #ifdef DEBUG
00223   G4cout << "BDSPCLDrift.cc: Finished making beam pipe..." << G4endl;
00224 #endif
00225 }
00226 
00227 void BDSPCLDrift::BuildBpFieldAndStepper(){
00228     // set up the magnetic field and stepper
00229   itsMagField=new BDSMagField(); //Zero magnetic field.
00230   itsEqRhs=new G4Mag_UsualEqRhs(itsMagField);
00231   itsStepper = new BDSDriftStepper(itsEqRhs);
00232 }
00233 
00234 G4VisAttributes* BDSPCLDrift::SetVisAttributes()
00235 {
00236   itsVisAttributes=new G4VisAttributes(G4Colour(0,1,0)); //useless
00237   itsVisAttributes->SetVisibility(false);
00238   itsVisAttributes->SetForceWireframe(true);
00239   return itsVisAttributes;
00240 }
00241 
00242 void BDSPCLDrift::BuildBLMs(){
00243   BDSAcceleratorComponent::BuildBLMs();
00244 }
00245 
00246 BDSPCLDrift::~BDSPCLDrift()
00247 {
00248   delete itsVisAttributes;
00249 //   delete outer_solid;
00250 //   delete inner_solid;
00251 //   delete itsOuterBeamPipeLogicalVolume;
00252 //   delete itsInnerBeamPipeLogicalVolume;
00253 //   delete itsPhysiOuter;
00254 //   delete itsPhysiInner;
00255   delete itsBeampipeVisAtt;
00256   delete itsInnerBeampipeVisAtt;
00257   delete itsMagField;
00258   delete itsEqRhs;
00259   delete itsStepper;
00260 }

Generated on 27 Aug 2013 for BDSIM by  doxygen 1.4.7