/home/cern/BDSIM_new/src/BDSDecapole.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" // must be first in include list
00011 
00012 #include "BDSDecapole.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 
00021 #include <map>
00022 
00023 const int DEBUG = 0;
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 BDSDecapole::BDSDecapole(G4String aName, G4double aLength, 
00037                          G4double bpRad, G4double FeRad,
00038                          G4double BQuadPrime, G4double tilt, 
00039                          G4double outR, G4String aMaterial):
00040   BDSMultipole(aName, aLength, bpRad, FeRad, SetVisAttributes(), aMaterial),
00041   itsBQuadPrime(BQuadPrime)
00042 {
00043   SetOuterRadius(outR);
00044   itsTilt=tilt;
00045   itsType="deca";
00046 
00047   if (!(*LogVolCount)[itsName])
00048     {
00049       //
00050       // build external volume
00051       // 
00052       BuildDefaultMarkerLogicalVolume();
00053 
00054       //
00055       // build beampipe (geometry + magnetic field)
00056       //
00057       BuildBPFieldAndStepper();
00058       BuildBPFieldMgr(itsStepper,itsMagField);
00059       BuildBeampipe(itsLength);
00060 
00061       //
00062       // build magnet (geometry + magnetic field)
00063       //
00064       BuildDefaultOuterLogicalVolume(itsLength);
00065       if(BDSGlobals->GetIncludeIronMagFields())
00066         {
00067           G4double polePos[4];
00068           G4double Bfield[3];
00069 
00070           //coordinate in GetFieldValue
00071           polePos[0]=-BDSGlobals->GetMagnetPoleRadius()*sin(pi/10);
00072           polePos[1]=BDSGlobals->GetMagnetPoleRadius()*cos(pi/10);
00073           polePos[2]=0.;
00074           polePos[3]=-999.;//flag to use polePos rather than local track
00075 
00076           itsMagField->GetFieldValue(polePos,Bfield);
00077           G4double BFldIron=
00078             sqrt(Bfield[0]*Bfield[0]+Bfield[1]*Bfield[1])*
00079             BDSGlobals->GetMagnetPoleSize()/
00080             (BDSGlobals->GetComponentBoxSize()/2-
00081              BDSGlobals->GetMagnetPoleRadius());
00082 
00083           // Magnetic flux from a pole is divided in two directions
00084           BFldIron/=2.;
00085 
00086           BuildOuterFieldManager(10, BFldIron,pi/10);
00087         }
00088 
00089       //
00090       // define sensitive volumes for hit generation
00091       //
00092       SetMultipleSensitiveVolumes(itsBeampipeLogicalVolume);
00093       SetMultipleSensitiveVolumes(itsOuterLogicalVolume);
00094 
00095       //
00096       // set visualization attributes
00097       //
00098       itsVisAttributes=SetVisAttributes();
00099       itsVisAttributes->SetForceSolid(true);
00100       itsOuterLogicalVolume->SetVisAttributes(itsVisAttributes);
00101 
00102       //
00103       // append marker logical volume to volume map
00104       //
00105       (*LogVolCount)[itsName]=1;
00106       (*LogVol)[itsName]=itsMarkerLogicalVolume;
00107     }
00108   else
00109     {
00110       (*LogVolCount)[itsName]++;
00111       if(BDSGlobals->GetSynchRadOn()&& BDSGlobals->GetSynchRescale())
00112         {
00113           // with synchrotron radiation, the rescaled magnetic field
00114           // means elements with the same name must have different
00115           //logical volumes, becuase they have different fields
00116           itsName+=BDSGlobals->StringFromInt((*LogVolCount)[itsName]);
00117 
00118           //
00119           // build external volume
00120           // 
00121           BuildDefaultMarkerLogicalVolume();
00122 
00123           //
00124           // build beampipe (geometry + magnetic field)
00125           //
00126           BuildBPFieldAndStepper();
00127           BuildBPFieldMgr(itsStepper,itsMagField);
00128           BuildBeampipe(itsLength);
00129 
00130           //
00131           // build magnet (geometry + magnetic field)
00132           //
00133           BuildDefaultOuterLogicalVolume(itsLength);
00134           if(BDSGlobals->GetIncludeIronMagFields())
00135             {
00136               G4double polePos[4];
00137               G4double Bfield[3];
00138               
00139               //coordinate in GetFieldValue
00140               polePos[0]=-BDSGlobals->GetMagnetPoleRadius()*sin(pi/10);
00141               polePos[1]=BDSGlobals->GetMagnetPoleRadius()*cos(pi/10);
00142               polePos[2]=0.;
00143               polePos[3]=-999.;//flag to use polePos rather than local track
00144 
00145               itsMagField->GetFieldValue(polePos,Bfield);
00146               G4double BFldIron=
00147                 sqrt(Bfield[0]*Bfield[0]+Bfield[1]*Bfield[1])*
00148                 BDSGlobals->GetMagnetPoleSize()/
00149                 (BDSGlobals->GetComponentBoxSize()/2-
00150                  BDSGlobals->GetMagnetPoleRadius());
00151 
00152               // Magnetic flux from a pole is divided in two directions
00153               BFldIron/=2.;
00154               
00155               BuildOuterFieldManager(10, BFldIron,pi/10);
00156             }
00157           //When is SynchRescale(factor) called?
00158 
00159           //
00160           // define sensitive volumes for hit generation
00161           //
00162           SetSensitiveVolume(itsBeampipeLogicalVolume);// for synchrotron
00163           //SetSensitiveVolume(itsOuterLogicalVolume);// for laserwire
00164           
00165           //
00166           // set visualization attributes
00167           //
00168           itsVisAttributes=SetVisAttributes();
00169           itsVisAttributes->SetForceSolid(true);
00170           itsOuterLogicalVolume->SetVisAttributes(itsVisAttributes);
00171           
00172           //
00173           // append marker logical volume to volume map
00174           //
00175           (*LogVol)[itsName]=itsMarkerLogicalVolume;
00176         }
00177       else
00178         {
00179           //
00180           // use already defined marker volume
00181           //
00182           itsMarkerLogicalVolume=(*LogVol)[itsName];
00183         }      
00184     }
00185 }
00186 
00187 void BDSDecapole::SynchRescale(G4double factor)
00188 {
00189   itsStepper->SetBQuadPrime(factor*itsBQuadPrime);
00190   itsMagField->SetBQuadPrime(factor*itsBQuadPrime);
00191   if(DEBUG) G4cout << "Dec " << itsName << " has been scaled" << G4endl;
00192 }
00193 
00194 G4VisAttributes* BDSDecapole::SetVisAttributes()
00195 {
00196   itsVisAttributes=new G4VisAttributes(G4Colour(0,0,1)); //green
00197   return itsVisAttributes;
00198 }
00199 
00200 void BDSDecapole::BuildBPFieldAndStepper()
00201 {
00202   // set up the magnetic field and stepper
00203   itsMagField=new BDSDecMagField(itsBQuadPrime);
00204   itsEqRhs=new G4Mag_UsualEqRhs(itsMagField);
00205   
00206   itsStepper=new BDSDecStepper(itsEqRhs);
00207   itsStepper->SetBQuadPrime(itsBQuadPrime);
00208 }
00209 
00210 BDSDecapole::~BDSDecapole()
00211 {
00212   delete itsVisAttributes;
00213   delete itsMarkerLogicalVolume;
00214   delete itsOuterLogicalVolume;
00215   delete itsPhysiComp;
00216   delete itsMagField;
00217   delete itsEqRhs;
00218   delete itsStepper;
00219 }

Generated on Wed Mar 5 17:25:22 2008 for BDSIM by  doxygen 1.5.3