src/BDSSextupole.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    Changed StringFromInt to be the BDSGlobal version
00008 */
00009 
00010 #include "BDSGlobalConstants.hh" 
00011 
00012 #include "BDSSextupole.hh"
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 "G4HelixImplicitEuler.hh"
00021 #include "G4CashKarpRKF45.hh"
00022 
00023 #include <map>
00024 
00025 //============================================================
00026 
00027 typedef std::map<G4String,int> LogVolCountMap;
00028 extern LogVolCountMap* LogVolCount;
00029 
00030 typedef std::map<G4String,G4LogicalVolume*> LogVolMap;
00031 extern LogVolMap* LogVol;
00032 
00033 extern BDSMaterials* theMaterials;
00034 //============================================================
00035 
00036 BDSSextupole::BDSSextupole(G4String aName, G4double aLength, 
00037                            G4double bpRad, G4double FeRad,
00038                            G4double BDblPrime, G4double tilt, 
00039                            G4double outR, 
00040                            std::list<G4double> blmLocZ, std::list<G4double> blmLocTheta,
00041                            G4String aTunnelMaterial, G4String aMaterial):
00042   BDSMultipole(aName, aLength, bpRad, FeRad, SetVisAttributes(), blmLocZ, blmLocTheta, aTunnelMaterial, aMaterial),
00043   itsBDblPrime(BDblPrime),
00044   itsStepper(NULL),itsMagField(NULL),itsEqRhs(NULL)
00045 {
00046   SetOuterRadius(outR);
00047   itsTilt=tilt;
00048   itsType="sext";
00049 
00050   if (!(*LogVolCount)[itsName])
00051     {
00052       //
00053       // build external volume
00054       // 
00055       BuildDefaultMarkerLogicalVolume();
00056 
00057       //
00058       //build tunnel
00059       //
00060       if(BDSGlobalConstants::Instance()->GetBuildTunnel()){
00061         BuildTunnel();
00062       }
00063 
00064       //
00065       // build beampipe (geometry + magnetic field)
00066       //
00067       BuildBPFieldAndStepper();
00068       BuildBPFieldMgr(itsStepper,itsMagField);
00069       BuildBeampipe();
00070 
00071       //
00072       // build magnet (geometry + magnetic field)
00073       //
00074       BuildDefaultOuterLogicalVolume(itsLength);
00075       if(BDSGlobalConstants::Instance()->GetIncludeIronMagFields())
00076         {
00077           G4double polePos[4];
00078           G4double Bfield[3];
00079 
00080           //coordinate in GetFieldValue
00081           polePos[0]=-BDSGlobalConstants::Instance()->GetMagnetPoleRadius()*sin(pi/6);
00082           polePos[1]=BDSGlobalConstants::Instance()->GetMagnetPoleRadius()*cos(pi/6);
00083           polePos[2]=0.;
00084           polePos[3]=-999.;//flag to use polePos rather than local track
00085 
00086           itsMagField->GetFieldValue(polePos,Bfield);
00087           G4double BFldIron=
00088             sqrt(Bfield[0]*Bfield[0]+Bfield[1]*Bfield[1])*
00089             BDSGlobalConstants::Instance()->GetMagnetPoleSize()/
00090             (BDSGlobalConstants::Instance()->GetComponentBoxSize()/2-
00091              BDSGlobalConstants::Instance()->GetMagnetPoleRadius());
00092 
00093           // Magnetic flux from a pole is divided in two directions
00094           BFldIron/=2.;
00095 
00096           BuildOuterFieldManager(6, BFldIron,pi/6);
00097         }
00098       //Build the beam loss monitors
00099       BuildBLMs();
00100       //
00101       // define sensitive volumes for hit generation
00102       //
00103       if(BDSGlobalConstants::Instance()->GetSensitiveBeamPipe()){
00104         SetMultipleSensitiveVolumes(itsBeampipeLogicalVolume);
00105       }
00106       if(BDSGlobalConstants::Instance()->GetSensitiveComponents()){
00107         SetMultipleSensitiveVolumes(itsOuterLogicalVolume);
00108       }
00109 
00110       //
00111       // set visualization attributes
00112       //
00113       itsVisAttributes=SetVisAttributes();
00114       itsVisAttributes->SetForceSolid(true);
00115       itsOuterLogicalVolume->SetVisAttributes(itsVisAttributes);
00116 
00117       //
00118       // append marker logical volume to volume map
00119       //
00120       (*LogVolCount)[itsName]=1;
00121       (*LogVol)[itsName]=itsMarkerLogicalVolume;
00122     }
00123   else
00124     {
00125       (*LogVolCount)[itsName]++;
00126       if(BDSGlobalConstants::Instance()->GetSynchRadOn()&& BDSGlobalConstants::Instance()->GetSynchRescale())
00127         {
00128           // with synchrotron radiation, the rescaled magnetic field
00129           // means elements with the same name must have different
00130           //logical volumes, becuase they have different fields
00131           itsName+=BDSGlobalConstants::Instance()->StringFromInt((*LogVolCount)[itsName]);
00132 
00133           //
00134           // build external volume
00135           // 
00136           BuildDefaultMarkerLogicalVolume();
00137 
00138           //
00139           // build beampipe (geometry + magnetic field)
00140           //
00141           BuildBPFieldAndStepper();
00142           BuildBPFieldMgr(itsStepper,itsMagField);
00143           BuildBeampipe();
00144 
00145           //
00146           // build magnet (geometry + magnetic field)
00147           //
00148           BuildDefaultOuterLogicalVolume(itsLength);
00149           if(BDSGlobalConstants::Instance()->GetIncludeIronMagFields())
00150             {
00151               G4double polePos[4];
00152               G4double Bfield[3];
00153               
00154               //coordinate in GetFieldValue
00155               polePos[0]=-BDSGlobalConstants::Instance()->GetMagnetPoleRadius()*sin(pi/6);
00156               polePos[1]=BDSGlobalConstants::Instance()->GetMagnetPoleRadius()*cos(pi/6);
00157               polePos[2]=0.;
00158               polePos[3]=-999.;//flag to use polePos rather than local track
00159 
00160               itsMagField->GetFieldValue(polePos,Bfield);
00161               G4double BFldIron=
00162                 sqrt(Bfield[0]*Bfield[0]+Bfield[1]*Bfield[1])*
00163                 BDSGlobalConstants::Instance()->GetMagnetPoleSize()/
00164                 (BDSGlobalConstants::Instance()->GetComponentBoxSize()/2-
00165                  BDSGlobalConstants::Instance()->GetMagnetPoleRadius());
00166 
00167               // Magnetic flux from a pole is divided in two directions
00168               BFldIron/=2.;
00169               
00170               BuildOuterFieldManager(6, BFldIron,pi/6);
00171             }
00172           //When is SynchRescale(factor) called?
00173 
00174           //
00175           // define sensitive volumes for hit generation
00176           //
00177           if(BDSGlobalConstants::Instance()->GetSensitiveBeamPipe()){
00178             SetMultipleSensitiveVolumes(itsBeampipeLogicalVolume);
00179           }
00180           if(BDSGlobalConstants::Instance()->GetSensitiveComponents()){
00181             SetMultipleSensitiveVolumes(itsOuterLogicalVolume);
00182           }
00183         
00184           //
00185           // set visualization attributes
00186           //
00187           itsVisAttributes=SetVisAttributes();
00188           itsVisAttributes->SetForceSolid(true);
00189           itsOuterLogicalVolume->SetVisAttributes(itsVisAttributes);
00190           
00191           //
00192           // append marker logical volume to volume map
00193           //
00194           (*LogVol)[itsName]=itsMarkerLogicalVolume;
00195         }
00196       else
00197         {
00198           //
00199           // use already defined marker volume
00200           //
00201           itsMarkerLogicalVolume=(*LogVol)[itsName];
00202         }      
00203     }
00204 }
00205 
00206 void BDSSextupole::SynchRescale(G4double factor)
00207 {
00208   itsStepper->SetBDblPrime(factor*itsBDblPrime);
00209   itsMagField->SetBDblPrime(factor*itsBDblPrime);
00210 #ifdef DEBUG 
00211   G4cout << "Sext " << itsName << " has been scaled" << G4endl;
00212 #endif
00213 }
00214 
00215 G4VisAttributes* BDSSextupole::SetVisAttributes()
00216 {
00217   itsVisAttributes=new G4VisAttributes(G4Colour(1,1,0));
00218   return itsVisAttributes;
00219 }
00220 
00221 
00222 void BDSSextupole::BuildBPFieldAndStepper()
00223 {
00224   // set up the magnetic field and stepper
00225   itsMagField=new BDSSextMagField(1*itsBDblPrime); //L Deacon testing field sign 4/7/12
00226   itsEqRhs=new G4Mag_UsualEqRhs(itsMagField);
00227 
00228   itsStepper=new BDSSextStepper(itsEqRhs);
00229   itsStepper->SetBDblPrime(itsBDblPrime);
00230 }
00231 
00232 BDSSextupole::~BDSSextupole()
00233 {
00234   delete itsVisAttributes;
00235   delete itsMagField;
00236   delete itsEqRhs;
00237   delete itsStepper;
00238 }

Generated on 27 Aug 2013 for BDSIM by  doxygen 1.4.7