/home/cern/BDSIM_new/src/BDSRBend.cc

00001 #include "BDSGlobalConstants.hh" // must be first in include list
00002 
00003 #include "BDSRBend.hh"
00004 #include "G4Tubs.hh"
00005 #include "G4IntersectionSolid.hh"
00006 #include "G4VisAttributes.hh"
00007 #include "G4LogicalVolume.hh"
00008 #include "G4VPhysicalVolume.hh"
00009 #include "G4UserLimits.hh"
00010 #include "G4TransportationManager.hh"
00011 
00012 #include <map>
00013 
00014 const int DEBUG = 0;
00015 
00016 //============================================================
00017 
00018 typedef std::map<G4String,int> LogVolCountMap;
00019 extern LogVolCountMap* LogVolCount;
00020 
00021 typedef std::map<G4String,G4LogicalVolume*> LogVolMap;
00022 extern LogVolMap* LogVol;
00023 
00024 extern BDSMaterials* theMaterials;
00025 
00026 extern G4RotationMatrix* RotY90;
00027 extern G4RotationMatrix* RotYM90;
00028 
00029 //============================================================
00030 
00031 // NOTE: BDSRBend is constructed by passing the geometric length,
00032 // but internally uses the arc length
00033 
00034 BDSRBend::BDSRBend(G4String aName, G4double aLength, 
00035                    G4double bpRad, G4double FeRad,
00036                    G4double bField, G4double angle, G4double outR,
00037                    G4double tilt, G4double bGrad, 
00038                    G4String aMaterial, G4int nSegments):
00039   BDSMultipole(aName, aLength, bpRad, FeRad, SetVisAttributes(), aMaterial,
00040                0, 0, angle)
00041 {
00042   SetOuterRadius(outR);
00043   itsTilt=tilt;
00044   itsBField=bField;
00045   itsBGrad=bGrad;
00046   itsType="rbend";
00047 
00048   if (!(*LogVolCount)[itsName])
00049     {
00050       //
00051       // build external volume
00052       // 
00053       BuildRBMarkerLogicalVolume();
00054 
00055       //
00056       // build beampipe (geometry + magnetic field)
00057       //
00058       BuildBPFieldAndStepper();
00059       BuildBPFieldMgr(itsStepper,itsMagField);
00060       BuildRBBeampipe();
00061 
00062       //
00063       // build magnet (geometry + magnetic field)
00064       //
00065       BuildRBOuterLogicalVolume();
00066       if(BDSGlobals->GetIncludeIronMagFields())
00067         {
00068           G4double polePos[4];
00069           G4double Bfield[3];
00070 
00071           //coordinate in GetFieldValue
00072           polePos[0]=0.;
00073           polePos[1]=BDSGlobals->GetMagnetPoleRadius();
00074           polePos[2]=0.;
00075           polePos[3]=-999.;//flag to use polePos rather than local track
00076 
00077           itsMagField->GetFieldValue(polePos,Bfield);
00078           G4double BFldIron=
00079           sqrt(Bfield[0]*Bfield[0]+Bfield[1]*Bfield[1])*
00080           BDSGlobals->GetMagnetPoleSize()/
00081           (BDSGlobals->GetComponentBoxSize()/2-
00082           BDSGlobals->GetMagnetPoleRadius());
00083 
00084           // Magnetic flux from a pole is divided in two directions
00085           BFldIron/=2.;
00086           
00087           BuildOuterFieldManager(2, BFldIron,pi/2);
00088         }
00089 
00090       //
00091       // define sensitive volumes for hit generation
00092       //
00093       SetMultipleSensitiveVolumes(itsBeampipeLogicalVolume);
00094       SetMultipleSensitiveVolumes(itsOuterLogicalVolume);
00095 
00096       //
00097       // set visualization attributes
00098       //
00099       itsVisAttributes=SetVisAttributes();
00100       itsVisAttributes->SetForceSolid(true);
00101       itsOuterLogicalVolume->SetVisAttributes(itsVisAttributes);
00102 
00103       //
00104       // append marker logical volume to volume map
00105       //
00106       (*LogVolCount)[itsName]=1;
00107       (*LogVol)[itsName]=itsMarkerLogicalVolume;
00108     }
00109   else
00110     {
00111       (*LogVolCount)[itsName]++;
00112       if(BDSGlobals->GetSynchRadOn()&& BDSGlobals->GetSynchRescale())
00113         {
00114           // with synchrotron radiation, the rescaled magnetic field
00115           // means elements with the same name must have different
00116           // logical volumes, because they have different fields
00117           itsName+=BDSGlobals->StringFromInt((*LogVolCount)[itsName]);
00118 
00119           //
00120           // build external volume
00121           // 
00122           BuildRBMarkerLogicalVolume();
00123           
00124           //
00125           // build beampipe (geometry + magnetic field)
00126           //
00127           BuildBPFieldAndStepper();
00128           BuildBPFieldMgr(itsStepper,itsMagField);
00129           BuildRBBeampipe();
00130 
00131           //
00132           // build magnet (geometry + magnetic field)
00133           //
00134           BuildRBOuterLogicalVolume();
00135           if(BDSGlobals->GetIncludeIronMagFields())
00136             {
00137               G4double polePos[4];
00138               G4double Bfield[3];
00139               
00140               //coordinate in GetFieldValue
00141               polePos[0]=0.;
00142               polePos[1]=BDSGlobals->GetMagnetPoleRadius();
00143               polePos[2]=0.;
00144               polePos[3]=-999.;//flag to use polePos rather than local track
00145               
00146               itsMagField->GetFieldValue(polePos,Bfield);
00147               G4double BFldIron=
00148                 sqrt(Bfield[0]*Bfield[0]+Bfield[1]*Bfield[1])*
00149                 BDSGlobals->GetMagnetPoleSize()/
00150                 (BDSGlobals->GetComponentBoxSize()/2-
00151                  BDSGlobals->GetMagnetPoleRadius());
00152               
00153               // Magnetic flux from a pole is divided in two directions
00154               BFldIron/=2.;
00155 
00156               BuildOuterFieldManager(2, BFldIron,pi/2);
00157             }
00158           //When is SynchRescale(factor) called?
00159           
00160           //
00161           // define sensitive volumes for hit generation
00162           //
00163           SetSensitiveVolume(itsBeampipeLogicalVolume);// for synchrotron
00164           //SetSensitiveVolume(itsOuterLogicalVolume);// for laserwire
00165 
00166           //
00167           // set visualization attributes
00168           //
00169           itsVisAttributes=SetVisAttributes();
00170           itsVisAttributes->SetForceSolid(true);
00171           itsOuterLogicalVolume->SetVisAttributes(itsVisAttributes);
00172 
00173           //
00174           // append marker logical volume to volume map
00175           //
00176           (*LogVol)[itsName]=itsMarkerLogicalVolume;
00177         }
00178       else
00179         {
00180           //
00181           // use already defined marker volume
00182           //
00183           itsMarkerLogicalVolume=(*LogVol)[itsName];
00184         }      
00185     }
00186 }
00187 
00188 void BDSRBend::SynchRescale(G4double factor)
00189 {
00190   // rescale B field and gradient by same factor
00191   itsStepper->SetBGrad(itsBGrad*factor);
00192   itsStepper->SetBField(-itsBField*factor);
00193   if(DEBUG) G4cout << "Rbend " << itsName << " has been scaled" << G4endl;
00194 }
00195 
00196 G4VisAttributes* BDSRBend::SetVisAttributes()
00197 {
00198   itsVisAttributes = new G4VisAttributes(G4Colour(0,0,1)); //blue
00199   return itsVisAttributes;
00200 }
00201 
00202 void BDSRBend::BuildBPFieldAndStepper()
00203 {
00204   // set up the magnetic field and stepper
00205   G4ThreeVector Bfield(0.,-itsBField,0.);
00206   itsMagField=new BDSSbendMagField(Bfield,GetArcLength(),itsAngle);
00207 
00208   itsEqRhs=new G4Mag_UsualEqRhs(itsMagField);  
00209   
00210   itsStepper = new myQuadStepper(itsEqRhs); // note the - sign...
00211   itsStepper->SetBField(-itsBField);
00212   itsStepper->SetBGrad(itsBGrad);
00213 }
00214 
00215 void BDSRBend::BuildRBMarkerLogicalVolume()
00216 {
00217   //  if (markerSolidVolume==0) {
00218 
00219     G4double boxSize=BDSGlobals->GetComponentBoxSize();
00220 
00221     // itsLength = geometrical length
00222     G4double xHalfLengthMinus = itsLength/2
00223       - fabs(cos(itsAngle/2))*boxSize*tan(itsAngle/2)/2
00224       + BDSGlobals->GetLengthSafety()/2;
00225     
00226     G4double xHalfLengthPlus = itsLength/2
00227       + fabs(cos(itsAngle/2))*boxSize*tan(itsAngle/2)/2
00228       + BDSGlobals->GetLengthSafety()/2;
00229 
00230     markerSolidVolume = new G4Trd(itsName+"_marker",
00231                                   xHalfLengthPlus,     // x hlf lgth at +z
00232                                   xHalfLengthMinus,    // x hlf lgth at -z
00233                                   boxSize/2,           // y hlf lgth at +z
00234                                   boxSize/2,           // y hlf lgth at -z
00235                                   fabs(cos(itsAngle/2))*boxSize/2);// z hlf len
00236     //  }
00237 
00238   G4String LocalLogicalName=itsName;
00239   
00240   itsMarkerLogicalVolume=    
00241     new G4LogicalVolume(markerSolidVolume,
00242                         theMaterials->GetMaterial("Vacuum"),
00243                         LocalLogicalName+"_marker");
00244 
00245   itsMarkerUserLimits = new G4UserLimits(DBL_MAX,DBL_MAX,DBL_MAX);
00246   itsMarkerUserLimits->SetMaxAllowedStep(GetArcLength());
00247   itsMarkerLogicalVolume->SetUserLimits(itsMarkerUserLimits);
00248 
00249   //
00250   // zero field in the marker volume
00251   //
00252   itsMarkerLogicalVolume->
00253     SetFieldManager(BDSGlobals->GetZeroFieldManager(),false);
00254 }
00255 
00256 
00257 // construct a beampipe for sector bend
00258 void BDSRBend::BuildRBBeampipe()
00259 {
00260   //
00261   // use default beampipe material
00262   //
00263   G4Material *material =  theMaterials->GetMaterial( BDSGlobals->GetPipeMaterialName());
00264   
00265   //
00266   // compute some geometrical parameters
00267   //
00268   G4double bpThickness = BDSGlobals->GetBeampipeThickness();
00269   G4double boxSize = BDSGlobals->GetComponentBoxSize();
00270 
00271   // itsLength = geometrical length
00272   G4double xHalfLengthMinus = itsLength/2
00273     - fabs(cos(itsAngle/2)) * boxSize * tan(itsAngle/2)/2
00274     + BDSGlobals->GetLengthSafety()/2;
00275 
00276   G4double xHalfLengthPlus = itsLength/2
00277     + fabs(cos(itsAngle/2)) * boxSize * tan(itsAngle/2)/2
00278     + BDSGlobals->GetLengthSafety()/2;
00279 
00280   G4double tubLen = std::max(xHalfLengthPlus,xHalfLengthMinus);
00281 
00282   //
00283   // build beampipe
00284   //
00285   G4Tubs *pipeTubsEnv = new G4Tubs(itsName+"_pipe_outer_env",
00286                                    itsBpRadius-bpThickness, // inner R
00287                                    itsBpRadius,             // outer R
00288                                    tubLen,                  // length
00289                                    0,                       // starting phi
00290                                    twopi * rad );           // delta phi
00291   
00292   G4Tubs *pipeInnerEnv = new G4Tubs(itsName+"_pipe_inner_env",
00293                                     0,                       // inner R
00294                                     itsBpRadius-bpThickness, // outer R
00295                                     tubLen,                  // length
00296                                     0,                       // starting phi
00297                                     twopi * rad );           // delta phi
00298 
00299   G4IntersectionSolid *pipeTubs =
00300     new G4IntersectionSolid(itsName+"_pipe_outer",
00301                             pipeTubsEnv,
00302                             markerSolidVolume,
00303                             RotYM90,
00304                             0);
00305   
00306   G4IntersectionSolid *pipeInner =
00307     new G4IntersectionSolid(itsName+"_pipe_inner",
00308                             pipeInnerEnv, 
00309                             markerSolidVolume,
00310                             RotYM90,
00311                             0);
00312 
00313   itsBeampipeLogicalVolume=     
00314     new G4LogicalVolume(pipeTubs,
00315                         material,
00316                         itsName+"_bmp_logical");
00317   
00318   itsInnerBPLogicalVolume=      
00319     new G4LogicalVolume(pipeInner,
00320                         theMaterials->GetMaterial("Vacuum"),
00321                         itsName+"_bmp_Inner_log");
00322 
00323   G4VPhysicalVolume* PhysiInner;
00324   PhysiInner = 
00325     new G4PVPlacement(
00326                       RotY90,                  // rotation
00327                       0,                       // at (0,0,0)
00328                       itsInnerBPLogicalVolume, // its logical volume
00329                       itsName+"_InnerBmp",     // its name
00330                       itsMarkerLogicalVolume,  // its mother volume
00331                       false,                   // no booleanm operation
00332                       0);                      // copy number
00333   
00334   G4VPhysicalVolume* PhysiComp;
00335   PhysiComp =
00336     new G4PVPlacement(
00337                       RotY90,                   // rotation
00338                       0,                        // at (0,0,0)
00339                       itsBeampipeLogicalVolume, // its logical volume
00340                       itsName+"_bmp",           // its name
00341                       itsMarkerLogicalVolume,   // its mother volume
00342                       false,                    // no boolean operation
00343                       0);                       // copy number
00344   
00345   //
00346   // set user limits for stepping, tracking and propagation in B field
00347   //
00348   itsBeampipeUserLimits =
00349     new G4UserLimits("beampipe cuts",DBL_MAX,DBL_MAX,DBL_MAX,
00350                      BDSGlobals->GetThresholdCutCharged());
00351   itsBeampipeUserLimits->SetMaxAllowedStep(GetArcLength());
00352   itsBeampipeLogicalVolume->SetUserLimits(itsBeampipeUserLimits);
00353   
00354   itsInnerBeampipeUserLimits =
00355     new G4UserLimits("inner beampipe cuts",DBL_MAX,DBL_MAX,DBL_MAX,
00356                      BDSGlobals->GetThresholdCutCharged());
00357   itsInnerBeampipeUserLimits->SetMaxAllowedStep(GetArcLength());
00358   itsInnerBPLogicalVolume->SetUserLimits(itsInnerBeampipeUserLimits);
00359 
00360   //
00361   // set magnetic field inside beampipe
00362   //
00363   itsBeampipeLogicalVolume->SetFieldManager(BDSGlobals->GetZeroFieldManager(),false);
00364   itsInnerBPLogicalVolume->SetFieldManager(itsBPFieldMgr,false);
00365 
00366   //
00367   // set visualization attributes
00368   //
00369   G4VisAttributes* VisAtt =  new G4VisAttributes(G4Colour(0., 0., 0));
00370   VisAtt->SetForceSolid(true);
00371   itsInnerBPLogicalVolume->SetVisAttributes(VisAtt);
00372 
00373   G4VisAttributes* VisAtt1 = new G4VisAttributes(G4Colour(0.4, 0.4, 0.4));
00374   VisAtt1->SetForceSolid(true);
00375   itsBeampipeLogicalVolume->SetVisAttributes(VisAtt1);
00376 }
00377 
00378 void BDSRBend::BuildRBOuterLogicalVolume(G4bool OuterMaterialIsVacuum){
00379 
00380   G4Material* material;
00381   if(itsMaterial != "")
00382     material = theMaterials->GetMaterial(itsMaterial);
00383   else
00384     material = theMaterials->GetMaterial("Iron");
00385 
00386   G4double boxSize = BDSGlobals->GetComponentBoxSize();
00387 
00388   // itsLength = geometrical length
00389   G4double xHalfLengthMinus = itsLength/2
00390     - fabs(cos(itsAngle/2)) * boxSize * tan(itsAngle/2)/2
00391     + BDSGlobals->GetLengthSafety()/2;
00392 
00393   G4double xHalfLengthPlus = itsLength/2
00394     + fabs(cos(itsAngle/2)) * boxSize * tan(itsAngle/2)/2
00395     + BDSGlobals->GetLengthSafety()/2;
00396 
00397   G4double tubLen = std::max(xHalfLengthPlus,xHalfLengthMinus);
00398   
00399   G4Tubs *magTubsEnv =
00400     new G4Tubs(itsName+"_solid_env",
00401                itsInnerIronRadius, // inner R
00402                itsOuterR,          // outer R
00403                tubLen,             // length
00404                0,                  // starting phi
00405                twopi * rad );      // delta phi
00406   
00407   G4IntersectionSolid *magTubs =
00408     new G4IntersectionSolid(itsName+"_solid",
00409                             magTubsEnv,
00410                             markerSolidVolume,
00411                             RotYM90,
00412                             0); 
00413 
00414   if(OuterMaterialIsVacuum)
00415     {
00416       itsOuterLogicalVolume = 
00417         new G4LogicalVolume(magTubs,
00418                             theMaterials->GetMaterial("Vacuum"),
00419                             itsName+"_outer");
00420     }
00421   else
00422     {
00423       itsOuterLogicalVolume=
00424         new G4LogicalVolume(magTubs,
00425                             material,
00426                             itsName+"_outer");
00427     }
00428 
00429   itsPhysiComp =
00430     new G4PVPlacement(
00431                       RotY90,                 // rotation
00432                       0,                      // at (0,0,0)
00433                       itsOuterLogicalVolume,  // its logical volume
00434                       itsName+"_solid",       // its name
00435                       itsMarkerLogicalVolume, // its mother  volume
00436                       false,                  // no boolean operation
00437                       0);                     // copy number
00438 
00439   itsOuterUserLimits =
00440     new G4UserLimits("multipole cut",DBL_MAX,DBL_MAX,DBL_MAX,
00441                      BDSGlobals->GetThresholdCutCharged());
00442   itsOuterUserLimits->SetMaxAllowedStep(GetArcLength());
00443   itsOuterLogicalVolume->SetUserLimits(itsOuterUserLimits);
00444 }
00445 
00446 BDSRBend::~BDSRBend()
00447 {
00448   delete itsVisAttributes;
00449   delete itsMarkerLogicalVolume;
00450   delete itsOuterLogicalVolume;
00451   delete itsPhysiComp;
00452   delete itsMagField;
00453   delete itsEqRhs;
00454   delete itsStepper;
00455   if (markerSolidVolume) delete markerSolidVolume;
00456 }

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