/home/cern/BDSIM_new/src/BDSKicker.cc

00001 #include "BDSGlobalConstants.hh" // must be first in include list
00002 
00003 #include "BDSKicker.hh"
00004 #include "G4VisAttributes.hh"
00005 #include "G4LogicalVolume.hh"
00006 #include "G4VPhysicalVolume.hh"
00007 #include "G4UserLimits.hh"
00008 #include "G4TransportationManager.hh"
00009 
00010 #include <map>
00011 
00012 const int DEBUG = 0;
00013 
00014 //============================================================
00015 
00016 typedef std::map<G4String,int> LogVolCountMap;
00017 extern LogVolCountMap* LogVolCount;
00018 
00019 typedef std::map<G4String,G4LogicalVolume*> LogVolMap;
00020 extern LogVolMap* LogVol;
00021 
00022 extern BDSMaterials* theMaterials;
00023 
00024 //============================================================
00025 
00026 BDSKicker::BDSKicker(G4String aName, G4double aLength, 
00027                      G4double bpRad, G4double FeRad,
00028                      G4double bField, G4double angle, G4double outR,
00029                      G4double tilt, G4double bGrad, 
00030                      G4String aMaterial, G4int nSegments):
00031   BDSMultipole(aName, aLength, bpRad, FeRad, SetVisAttributes(), aMaterial,
00032                0, 0, angle)
00033 {
00034   SetOuterRadius(outR);
00035   itsTilt=tilt;
00036   itsBField=bField;
00037   itsBGrad=bGrad;
00038   if (tilt==0)
00039     itsType="hkick";
00040   else if (tilt==pi/2)
00041     itsType="vkick";
00042   else
00043     itsType="kick";
00044 
00045   if (!(*LogVolCount)[itsName])
00046     {
00047       //
00048       // build external volume
00049       // 
00050       BuildDefaultMarkerLogicalVolume();
00051 
00052       //
00053       // build beampipe (geometry + magnetic field)
00054       //
00055       BuildBPFieldAndStepper();
00056       BuildBPFieldMgr(itsStepper,itsMagField);
00057       BuildBeampipe(itsLength);
00058 
00059       //
00060       // build magnet (geometry + magnetic field)
00061       //
00062       BuildDefaultOuterLogicalVolume(itsLength);
00063       if(BDSGlobals->GetIncludeIronMagFields())
00064         {
00065           G4double polePos[4];
00066           G4double Bfield[3];
00067 
00068           //coordinate in GetFieldValue
00069           polePos[0]=0.;
00070           polePos[1]=BDSGlobals->GetMagnetPoleRadius();
00071           polePos[2]=0.;
00072           polePos[0]=-999.;//flag to use polePos rather than local track
00073           itsMagField->GetFieldValue(polePos,Bfield);
00074           G4double BFldIron=
00075             sqrt(Bfield[0]*Bfield[0]+Bfield[1]*Bfield[1])*
00076             BDSGlobals->GetMagnetPoleSize()/
00077             (BDSGlobals->GetComponentBoxSize()/2-
00078              BDSGlobals->GetMagnetPoleRadius());
00079 
00080           // Magnetic flux from a pole is divided in two directions
00081           BFldIron/=2.;
00082 
00083           BuildOuterFieldManager(2, BFldIron,0);
00084         }
00085 
00086       //
00087       // define sensitive volumes for hit generation
00088       //
00089       SetMultipleSensitiveVolumes(itsBeampipeLogicalVolume);
00090       SetMultipleSensitiveVolumes(itsOuterLogicalVolume);
00091 
00092       //
00093       // set visualization attributes
00094       //
00095       itsVisAttributes->SetForceSolid(true);
00096       itsOuterLogicalVolume->SetVisAttributes(itsVisAttributes);
00097 
00098       //
00099       // append marker logical volume to volume map
00100       //
00101       (*LogVolCount)[itsName]=1;
00102       (*LogVol)[itsName]=itsMarkerLogicalVolume;
00103     }
00104   else
00105     {
00106       (*LogVolCount)[itsName]++;
00107 
00108       //
00109       // no rescaling ???
00110       //
00111  
00112       //
00113       // use already defined marker volume
00114       //
00115       itsMarkerLogicalVolume=(*LogVol)[itsName];
00116     }
00117 }
00118 
00119 
00120 G4VisAttributes* BDSKicker::SetVisAttributes()
00121 {
00122   itsVisAttributes=new G4VisAttributes(G4Colour(0,0,1));
00123   return itsVisAttributes;
00124 }
00125 
00126 
00127 void BDSKicker::BuildBPFieldAndStepper()
00128 {
00129   // set up the magnetic field and stepper
00130   G4ThreeVector Bfield(0.,-itsBField,0.); // note the - sign...
00131   itsMagField=new BDSSbendMagField(Bfield,itsLength,itsAngle);
00132   
00133   itsEqRhs=new G4Mag_UsualEqRhs(itsMagField);  
00134   
00135   itsStepper = new myQuadStepper(itsEqRhs);
00136   itsStepper->SetBField(-itsBField); // note the - sign...
00137   itsStepper->SetBGrad(itsBGrad);
00138 }
00139 
00140 BDSKicker::~BDSKicker()
00141 {
00142   
00143   delete itsVisAttributes;
00144   delete itsMarkerLogicalVolume;
00145   delete itsOuterLogicalVolume;
00146   delete itsPhysiComp;
00147   delete itsMagField;
00148   delete itsEqRhs;
00149   delete itsStepper;
00150 }

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