src/BDSRBend.cc

00001 #include "BDSGlobalConstants.hh" 
00002 
00003 #include "BDSRBend.hh"
00004 #include "G4Tubs.hh"
00005 #include "G4Box.hh"
00006 #include "G4Trd.hh"
00007 #include "G4IntersectionSolid.hh"
00008 #include "G4SubtractionSolid.hh"
00009 #include "G4UnionSolid.hh"
00010 #include "G4VisAttributes.hh"
00011 #include "G4LogicalVolume.hh"
00012 #include "G4VPhysicalVolume.hh"
00013 #include "G4UserLimits.hh"
00014 #include "G4TransportationManager.hh"
00015 
00016 #include <map>
00017 
00018 //============================================================
00019 
00020 typedef std::map<G4String,int> LogVolCountMap;
00021 extern LogVolCountMap* LogVolCount;
00022 
00023 typedef std::map<G4String,G4LogicalVolume*> LogVolMap;
00024 extern LogVolMap* LogVol;
00025 
00026 extern G4RotationMatrix* RotY90;
00027 extern G4RotationMatrix* RotYM90;
00028 
00029 //============================================================
00030 
00031 BDSRBend::BDSRBend(G4String aName, G4double aLength, 
00032                    G4double bpRad, G4double FeRad,
00033                    G4double bField, G4double angle, G4double outR,
00034                    std::list<G4double> blmLocZ, std::list<G4double> blmLocTheta,
00035                    G4double tilt, G4double bGrad, 
00036                    G4String aTunnelMaterial, G4String aMaterial):
00037   BDSMultipole(aName, aLength, bpRad, FeRad, SetVisAttributes(), blmLocZ, blmLocTheta, aTunnelMaterial, aMaterial,
00038                0, 0, angle),
00039   markerSolidVolume(NULL),rbendRectangleSolidVolume(NULL),rbendRectangleLogicalVolume(NULL),
00040   middleBeampipeLogicalVolume(NULL),middleInnerBPLogicalVolume(NULL),endsBeampipeLogicalVolume(NULL),
00041   endsInnerBPLogicalVolume(NULL),endsBeampipeUserLimits(NULL),endsInnerBeampipeUserLimits(NULL),
00042   innerBeampipeVisAtt(NULL),beampipeVisAtt(NULL),itsStepper(NULL),itsMagField(NULL),itsEqRhs(NULL)
00043 {
00044   SetOuterRadius(outR);
00045   itsTilt=tilt;
00046   itsBField=bField;
00047   itsBGrad=bGrad;
00048   //  itsBpRadius*=2; //bending magnet so double beam pipe radius
00049   itsType="rbend";
00050   
00051   //magnetic field length for rbend
00052   itsMagFieldLength = std::min (
00053                                 ((itsLength/itsAngle)*sin(itsAngle/2)
00054                                  - fabs(cos(itsAngle/2))*outR*tan(itsAngle/2)/2), 
00055                                 ((itsLength/itsAngle)*sin(itsAngle/2)
00056                                  + fabs(cos(itsAngle/2))*outR*tan(itsAngle/2)/2)
00057                                 );
00058   
00059   G4cout << "BDSRBend>> rbend itsMagFieldLength = " << itsMagFieldLength << G4endl;
00060 
00061   if (!(*LogVolCount)[itsName])
00062     {
00063       //
00064       // build external volume
00065       // 
00066       BuildRBMarkerLogicalVolume();
00067 
00068       //
00069       //build tunnel 
00070       if(BDSGlobalConstants::Instance()->GetBuildTunnel()){
00071         BuildTunnel(); //Geometry problem, do not build tunnel
00072       }
00073       
00074       //
00075       // build beampipe (geometry + magnetic field)
00076       //
00077       BuildBPFieldAndStepper();
00078       BuildBPFieldMgr(itsStepper,itsMagField);
00079       BuildRBBeampipe();
00080 
00081       //
00082       // build magnet (geometry + magnetic field)
00083       //
00084       BuildRBOuterLogicalVolume();
00085       if(BDSGlobalConstants::Instance()->GetIncludeIronMagFields())
00086         {
00087           G4double polePos[4];
00088           G4double Bfield[3];
00089 
00090           //coordinate in GetFieldValue
00091           polePos[0]=0.;
00092           polePos[1]=BDSGlobalConstants::Instance()->GetMagnetPoleRadius();
00093           polePos[2]=0.;
00094           polePos[3]=-999.;//flag to use polePos rather than local track
00095 
00096           itsMagField->GetFieldValue(polePos,Bfield);
00097           G4double BFldIron=
00098             sqrt(Bfield[0]*Bfield[0]+Bfield[1]*Bfield[1])*
00099             BDSGlobalConstants::Instance()->GetMagnetPoleSize()/
00100             (BDSGlobalConstants::Instance()->GetComponentBoxSize()/2-
00101              BDSGlobalConstants::Instance()->GetMagnetPoleRadius());
00102 
00103           // Magnetic flux from a pole is divided in two directions
00104           BFldIron/=2.;
00105           
00106           BuildOuterFieldManager(2, BFldIron,pi/2);
00107         }
00108       BuildBLMs();
00109       //
00110       // define sensitive volumes for hit generation
00111       //
00112       if(BDSGlobalConstants::Instance()->GetSensitiveBeamPipe()){
00113         SetMultipleSensitiveVolumes(itsBeampipeLogicalVolume);
00114       }
00115       if(BDSGlobalConstants::Instance()->GetSensitiveComponents()){
00116         SetMultipleSensitiveVolumes(itsOuterLogicalVolume);
00117       }
00118 
00119       //
00120       // set visualization attributes
00121       //
00122       itsVisAttributes=SetVisAttributes();
00123       itsVisAttributes->SetForceSolid(true);
00124       itsOuterLogicalVolume->SetVisAttributes(itsVisAttributes);
00125 
00126       //
00127       // append marker logical volume to volume map
00128       //
00129       (*LogVolCount)[itsName]=1;
00130       (*LogVol)[itsName]=itsMarkerLogicalVolume;
00131     }
00132   else
00133     {
00134       (*LogVolCount)[itsName]++;
00135       if(BDSGlobalConstants::Instance()->GetSynchRadOn()&& BDSGlobalConstants::Instance()->GetSynchRescale())
00136         {
00137           // with synchrotron radiation, the rescaled magnetic field
00138           // means elements with the same name must have different
00139           // logical volumes, because they have different fields
00140           itsName+=BDSGlobalConstants::Instance()->StringFromInt((*LogVolCount)[itsName]);
00141 
00142           //
00143           // build external volume
00144           // 
00145           BuildRBMarkerLogicalVolume();
00146           
00147           //
00148           // build beampipe (geometry + magnetic field)
00149           //
00150           BuildBPFieldAndStepper();
00151           BuildBPFieldMgr(itsStepper,itsMagField);
00152           BuildRBBeampipe();
00153 
00154           //
00155           // build magnet (geometry + magnetic field)
00156           //
00157           BuildRBOuterLogicalVolume();
00158           if(BDSGlobalConstants::Instance()->GetIncludeIronMagFields())
00159             {
00160               G4double polePos[4];
00161               G4double Bfield[3];
00162               
00163               //coordinate in GetFieldValue
00164               polePos[0]=0.;
00165               polePos[1]=BDSGlobalConstants::Instance()->GetMagnetPoleRadius();
00166               polePos[2]=0.;
00167               polePos[3]=-999.;//flag to use polePos rather than local track
00168               
00169               itsMagField->GetFieldValue(polePos,Bfield);
00170               G4double BFldIron=
00171                 sqrt(Bfield[0]*Bfield[0]+Bfield[1]*Bfield[1])*
00172                 BDSGlobalConstants::Instance()->GetMagnetPoleSize()/
00173                 (BDSGlobalConstants::Instance()->GetComponentBoxSize()/2-
00174                  BDSGlobalConstants::Instance()->GetMagnetPoleRadius());
00175               
00176               // Magnetic flux from a pole is divided in two directions
00177               BFldIron/=2.;
00178 
00179               BuildOuterFieldManager(2, BFldIron,pi/2);
00180             }
00181           //When is SynchRescale(factor) called?
00182           
00183           //
00184           // define sensitive volumes for hit generation
00185           //
00186           SetSensitiveVolume(middleBeampipeLogicalVolume);// for synchrotron
00187           SetSensitiveVolume(endsBeampipeLogicalVolume);// for synchrotron
00188           SetSensitiveVolume(itsOuterLogicalVolume);// for laserwire
00189 
00190           //
00191           // set visualization attributes
00192           //
00193           itsVisAttributes=SetVisAttributes();
00194           itsVisAttributes->SetForceSolid(true);
00195           itsOuterLogicalVolume->SetVisAttributes(itsVisAttributes);
00196 
00197           //
00198           // append marker logical volume to volume map
00199           //
00200           (*LogVol)[itsName]=itsMarkerLogicalVolume;
00201         }
00202       else
00203         {
00204           //
00205           // use already defined marker volume
00206           //
00207           itsMarkerLogicalVolume=(*LogVol)[itsName];
00208         }      
00209     }
00210 }
00211 
00212 void BDSRBend::SynchRescale(G4double factor)
00213 {
00214   // rescale B field and gradient by same factor
00215   itsStepper->SetBGrad(itsBGrad*factor);
00216   itsStepper->SetBField(-itsBField*factor);
00217   // note that there are no methods to set the BDSRBendMagField as this
00218   // class does not do anything with the BFields.
00219   // not true when I will use Geant4 propagation
00220 #ifdef DEBUG
00221   G4cout << "Sbend " << itsName << " has been scaled" << G4endl;
00222 #endif
00223 }
00224 
00225 G4VisAttributes* BDSRBend::SetVisAttributes()
00226 {
00227   itsVisAttributes = new G4VisAttributes(G4Colour(0,0,1)); //blue
00228   return itsVisAttributes;
00229 }
00230 
00231 void BDSRBend::BuildBPFieldAndStepper()
00232 {
00233   // set up the magnetic field and stepper
00234   G4ThreeVector Bfield(0.,-itsBField,0.);
00235   itsMagField=new BDSSbendMagField(Bfield,itsMagFieldLength,itsAngle);
00236 
00237   itsEqRhs=new G4Mag_UsualEqRhs(itsMagField);  
00238   
00239   itsStepper = new myQuadStepper(itsEqRhs);
00240   itsStepper->SetBField(-itsBField); // note the - sign...
00241   itsStepper->SetBGrad(itsBGrad);
00242 }
00243 
00244 void BDSRBend::BuildRBMarkerLogicalVolume()
00245 {
00246   if (markerSolidVolume==0) {
00247 
00248     G4double boxSize=BDSGlobalConstants::Instance()->GetComponentBoxSize()+BDSGlobalConstants::Instance()->GetTunnelRadius();
00249 
00250     G4double xHalfLengthMinus = (itsLength/itsAngle)*sin(itsAngle/2)
00251       - fabs(cos(itsAngle/2))*boxSize*tan(itsAngle/2)/2;
00252     
00253     G4double xHalfLengthPlus = (itsLength/itsAngle)*sin(itsAngle/2)
00254       + fabs(cos(itsAngle/2))*boxSize*tan(itsAngle/2)/2;
00255 
00256     markerSolidVolume = new G4Trd(itsName+"_marker",
00257                                   xHalfLengthPlus,     // x hlf lgth at +z
00258                                   xHalfLengthMinus,    // x hlf lgth at -z
00259                                   boxSize/2,           // y hlf lgth at +z
00260                                   boxSize/2,           // y hlf lgth at -z
00261                                   fabs(cos(itsAngle/2))*boxSize/2);// z hlf lgth
00262     
00263 
00264     rbendRectangleSolidVolume = new G4Trd(itsName+"_rbend_rectangle",
00265                                           itsMagFieldLength,     
00266                                           itsMagFieldLength,    
00267                                           boxSize/2,          
00268                                           boxSize/2,           
00269                                           fabs(cos(itsAngle/2))*boxSize/2);
00270     
00271   }
00272 
00273   G4String LocalLogicalName=itsName;
00274   
00275   itsMarkerLogicalVolume=    
00276     new G4LogicalVolume(markerSolidVolume,
00277                         BDSMaterials::Instance()->GetMaterial(BDSGlobalConstants::Instance()->GetVacuumMaterial()),
00278                         LocalLogicalName+"_marker");
00279 
00280   rbendRectangleLogicalVolume=    
00281     new G4LogicalVolume(rbendRectangleSolidVolume,
00282                         BDSMaterials::Instance()->GetMaterial(BDSGlobalConstants::Instance()->GetVacuumMaterial()),
00283                         LocalLogicalName+"_rbend_rectangle");
00284 #ifndef NOUSERLIMITS
00285   itsMarkerUserLimits = new G4UserLimits(DBL_MAX,DBL_MAX,DBL_MAX);
00286   itsMarkerUserLimits->SetMaxAllowedStep(itsMagFieldLength*5);
00287   itsMarkerLogicalVolume->SetUserLimits(itsMarkerUserLimits);
00288   rbendRectangleLogicalVolume->SetUserLimits(itsMarkerUserLimits);
00289 #endif
00290   //
00291   // zero field in the marker volume
00292   //
00293   itsMarkerLogicalVolume->
00294     SetFieldManager(BDSGlobalConstants::Instance()->GetZeroFieldManager(),false);
00295 }
00296 
00297 
00298 // construct a beampipe for r bend
00299 void BDSRBend::BuildRBBeampipe()
00300 {
00301   //
00302   // use default beampipe material
00303   //
00304   G4Material *material =  BDSMaterials::Instance()->GetMaterial( BDSGlobalConstants::Instance()->GetPipeMaterialName());
00305   
00306   //
00307   // compute some geometrical parameters
00308   //
00309   G4double bpThickness = BDSGlobalConstants::Instance()->GetBeampipeThickness();
00310   G4double boxSize = BDSGlobalConstants::Instance()->GetComponentBoxSize();
00311 
00312   G4double xHalfLengthMinus =
00313     (itsLength/itsAngle)*sin(itsAngle/2)
00314     - fabs(cos(itsAngle/2)) * boxSize * tan(itsAngle/2)/2;
00315 
00316 
00317   G4double xHalfLengthPlus =
00318     (itsLength/itsAngle)*sin(itsAngle/2)
00319     + fabs(cos(itsAngle/2)) * boxSize * tan(itsAngle/2)/2;
00320 
00321 
00322   G4double tubLen = std::max(xHalfLengthPlus,xHalfLengthMinus);
00323 
00324   //
00325   // build beampipe
00326   //
00327   G4Tubs *pipeTubsEnv = new G4Tubs(itsName+"_pipe_outer_env",
00328                                    itsBpRadius+BDSGlobalConstants::Instance()->GetLengthSafety()/2.0, // inner R
00329                                    itsBpRadius+bpThickness,             // outer R
00330                                    tubLen,                  // length
00331                                    0,                       // starting phi
00332                                    twopi * rad );           // delta phi
00333   
00334   G4Tubs *pipeInnerEnv = new G4Tubs(itsName+"_pipe_inner_env",
00335                                     0,                       // inner R
00336                                     itsBpRadius, // outer R
00337                                     tubLen,                  // length
00338                                     0,                       // starting phi
00339                                     twopi * rad );           // delta phi
00340 
00341   G4IntersectionSolid *pipeTubs =
00342     new G4IntersectionSolid(itsName+"_pipe_outer",
00343                             pipeTubsEnv,
00344                             markerSolidVolume,
00345                             RotYM90,
00346                             (G4ThreeVector)0);
00347   
00348   G4IntersectionSolid *pipeInner =
00349     new G4IntersectionSolid(itsName+"_pipe_inner",
00350                             pipeInnerEnv, 
00351                             markerSolidVolume,
00352                             RotYM90,
00353                             (G4ThreeVector)0);
00354 
00355   G4IntersectionSolid *pipeMiddleTubs =
00356     new G4IntersectionSolid(itsName+"_pipe_middle_outer",
00357                             pipeTubsEnv,
00358                             rbendRectangleSolidVolume,
00359                             RotYM90,
00360                             (G4ThreeVector)0);
00361   
00362   G4IntersectionSolid *pipeMiddleInner =
00363     new G4IntersectionSolid(itsName+"_pipe_middle_inner",
00364                             pipeInnerEnv, 
00365                             rbendRectangleSolidVolume,
00366                             RotYM90,
00367                             (G4ThreeVector)0);
00368 
00369   G4SubtractionSolid *pipeEndsTubsTmp =
00370     new G4SubtractionSolid(itsName+"_pipe_ends_outer_tmp",
00371                             pipeTubsEnv,
00372                             rbendRectangleSolidVolume,
00373                             RotYM90,
00374                             (G4ThreeVector)0);
00375   
00376   G4SubtractionSolid *pipeEndsInnerTmp =
00377     new G4SubtractionSolid(itsName+"_pipe_ends_inner_tmp",
00378                             pipeInnerEnv, 
00379                             rbendRectangleSolidVolume,
00380                             RotYM90,
00381                             (G4ThreeVector)0);
00382 
00383   G4IntersectionSolid *pipeEndsTubs =
00384     new G4IntersectionSolid(itsName+"_pipe_ends_outer",
00385                             pipeEndsTubsTmp,
00386                             markerSolidVolume,
00387                             RotYM90,
00388                             (G4ThreeVector)0);
00389   
00390   G4IntersectionSolid *pipeEndsInner =
00391     new G4IntersectionSolid(itsName+"_pipe_ends_inner",
00392                             pipeEndsInnerTmp, 
00393                             markerSolidVolume,
00394                             RotYM90,
00395                             (G4ThreeVector)0);
00396 
00397   itsBeampipeLogicalVolume=     
00398     new G4LogicalVolume(pipeTubs,
00399                         material,
00400                         itsName+"_bmp_logical");
00401   
00402   itsInnerBPLogicalVolume=      
00403     new G4LogicalVolume(pipeInner,
00404                         BDSMaterials::Instance()->GetMaterial(BDSGlobalConstants::Instance()->GetVacuumMaterial()),
00405                         itsName+"_bmp_Inner_log");
00406 
00407   middleBeampipeLogicalVolume=  
00408     new G4LogicalVolume(pipeMiddleTubs,
00409                         material,
00410                         itsName+"_bmp_logical_middle");
00411   
00412  middleInnerBPLogicalVolume=    
00413     new G4LogicalVolume(pipeMiddleInner,
00414                         BDSMaterials::Instance()->GetMaterial(BDSGlobalConstants::Instance()->GetVacuumMaterial()),
00415                         itsName+"_bmp_Inner_log_middle");
00416 
00417  endsBeampipeLogicalVolume=     
00418     new G4LogicalVolume(pipeEndsTubs,
00419                         material,
00420                         itsName+"_bmp_logical_ends");
00421   
00422  endsInnerBPLogicalVolume=      
00423     new G4LogicalVolume(pipeEndsInner,
00424                         BDSMaterials::Instance()->GetMaterial(BDSGlobalConstants::Instance()->GetVacuumMaterial()),
00425                         itsName+"_bmp_Inner_log_ends");
00426 
00427 
00428 
00429    //
00430   // set magnetic field inside beampipe
00431   //
00432  middleBeampipeLogicalVolume->SetFieldManager(BDSGlobalConstants::Instance()->GetZeroFieldManager(),false);
00433  middleInnerBPLogicalVolume->SetFieldManager(itsBPFieldMgr,false);
00434 
00435 
00436   G4VPhysicalVolume* PhysiInner;
00437   PhysiInner = 
00438     new G4PVPlacement(
00439                       RotY90,                  // rotation
00440                       (G4ThreeVector)0,                // at (0,0,0)
00441                       middleInnerBPLogicalVolume, // its logical volume
00442                       itsName+"_InnerBmp",     // its name
00443                       itsMarkerLogicalVolume,  // its mother volume
00444                       false,                   // no booleanm operation
00445                       0, BDSGlobalConstants::Instance()->GetCheckOverlaps());                  // copy number
00446 
00447 
00448   G4VPhysicalVolume* PhysiComp;
00449   PhysiComp =
00450     new G4PVPlacement(
00451                       RotY90,                   // rotation
00452                       (G4ThreeVector)0,                 // at (0,0,0)
00453                       middleBeampipeLogicalVolume, // its logical volume
00454                       itsName+"_bmp",           // its name
00455                       itsMarkerLogicalVolume,   // its mother volume
00456                       false,                    // no boolean operation
00457                       0, BDSGlobalConstants::Instance()->GetCheckOverlaps()
00458 );                      // copy number
00459 
00460   G4VPhysicalVolume* PhysiInnerEnds;
00461   PhysiInnerEnds = 
00462     new G4PVPlacement(
00463                       RotY90,                  // rotation
00464                       (G4ThreeVector)0,               // at (0,0,0)
00465                       endsInnerBPLogicalVolume, // its logical volume
00466                       itsName+"_InnerBmp",     // its name
00467                       itsMarkerLogicalVolume,  // its mother volume
00468                       false,                   // no booleanm operation
00469                       0, BDSGlobalConstants::Instance()->GetCheckOverlaps() );                 // copy number
00470   
00471   G4VPhysicalVolume* PhysiCompEnds;
00472   PhysiCompEnds =
00473     new G4PVPlacement(
00474                       RotY90,                   // rotation
00475                       (G4ThreeVector)0,                 // at (0,0,0)
00476                       endsBeampipeLogicalVolume, // its logical volume
00477                       itsName+"_bmp",           // its name
00478                       itsMarkerLogicalVolume,   // its mother volume
00479                       false,                    // no boolean operation
00480                       0, BDSGlobalConstants::Instance()->GetCheckOverlaps());                   // copy number
00481   
00482   SetMultiplePhysicalVolumes(PhysiInner);
00483   SetMultiplePhysicalVolumes(PhysiComp);
00484   SetMultiplePhysicalVolumes(PhysiInnerEnds);
00485   SetMultiplePhysicalVolumes(PhysiCompEnds);
00486   
00487 #ifndef NOUSERLIMITS
00488   //
00489   // set user limits for stepping, tracking and propagation in B field
00490   //
00491   itsBeampipeUserLimits =
00492     new G4UserLimits("beampipe cuts",DBL_MAX,DBL_MAX,DBL_MAX,
00493                      BDSGlobalConstants::Instance()->GetThresholdCutCharged());
00494   itsBeampipeUserLimits->SetMaxAllowedStep(itsMagFieldLength);
00495   itsBeampipeUserLimits->SetUserMaxTime(BDSGlobalConstants::Instance()->GetMaxTime());
00496   middleBeampipeLogicalVolume->SetUserLimits(itsBeampipeUserLimits);
00497 
00498   G4double endsMaxAllowedStep = itsMagFieldLength;
00499 
00500   endsBeampipeUserLimits =
00501     new G4UserLimits("beampipe ends cuts",DBL_MAX,DBL_MAX,DBL_MAX,
00502                      BDSGlobalConstants::Instance()->GetThresholdCutCharged());
00503   endsBeampipeUserLimits->SetMaxAllowedStep(endsMaxAllowedStep);
00504   endsBeampipeUserLimits->SetUserMaxTime(BDSGlobalConstants::Instance()->GetMaxTime());
00505   endsBeampipeLogicalVolume->SetUserLimits(endsBeampipeUserLimits);
00506   
00507   itsInnerBeampipeUserLimits =
00508     new G4UserLimits("inner beampipe cuts",DBL_MAX,DBL_MAX,DBL_MAX,
00509                      BDSGlobalConstants::Instance()->GetThresholdCutCharged());
00510   itsInnerBeampipeUserLimits->SetMaxAllowedStep(itsMagFieldLength);
00511   itsInnerBeampipeUserLimits->SetUserMaxTime(BDSGlobalConstants::Instance()->GetMaxTime());
00512   middleInnerBPLogicalVolume->SetUserLimits(itsInnerBeampipeUserLimits);
00513 
00514   endsInnerBeampipeUserLimits =
00515     new G4UserLimits("inner beampipe ends cuts",DBL_MAX,DBL_MAX,DBL_MAX,
00516                      BDSGlobalConstants::Instance()->GetThresholdCutCharged());
00517   endsInnerBeampipeUserLimits->SetMaxAllowedStep(endsMaxAllowedStep);
00518   endsInnerBeampipeUserLimits->SetUserMaxTime(BDSGlobalConstants::Instance()->GetMaxTime());
00519   endsInnerBPLogicalVolume->SetUserLimits(endsInnerBeampipeUserLimits);
00520 #endif
00521 
00522   //
00523   // set visualization attributes
00524   //
00525   innerBeampipeVisAtt =  new G4VisAttributes(G4Colour(0., 0., 0));
00526   innerBeampipeVisAtt->SetForceSolid(true);
00527   middleInnerBPLogicalVolume->SetVisAttributes(innerBeampipeVisAtt);
00528   endsInnerBPLogicalVolume->SetVisAttributes(innerBeampipeVisAtt);
00529 
00530   beampipeVisAtt = new G4VisAttributes(G4Colour(0.4, 0.4, 0.4));
00531   beampipeVisAtt->SetForceSolid(true);
00532   middleBeampipeLogicalVolume->SetVisAttributes(beampipeVisAtt);
00533   endsBeampipeLogicalVolume->SetVisAttributes(beampipeVisAtt);
00534 }
00535 
00536 void BDSRBend::BuildRBOuterLogicalVolume(G4bool OuterMaterialIsVacuum){
00537 
00538   G4Material* material;
00539   if(itsMaterial != "")
00540     material = BDSMaterials::Instance()->GetMaterial(itsMaterial);
00541   else
00542     material = BDSMaterials::Instance()->GetMaterial("Iron");
00543 
00544   G4double boxSize = BDSGlobalConstants::Instance()->GetComponentBoxSize();
00545  
00546   G4double xHalfLengthMinus = (itsLength/itsAngle)*sin(itsAngle/2)
00547     - fabs(cos(itsAngle/2)) * boxSize * tan(itsAngle/2)/2;
00548 
00549   G4double xHalfLengthPlus = (itsLength/itsAngle)*sin(itsAngle/2)
00550     + fabs(cos(itsAngle/2)) * boxSize * tan(itsAngle/2)/2;
00551 
00552   G4double tubLen = std::max(xHalfLengthPlus,xHalfLengthMinus);
00553   
00554   G4Tubs *magTubsEnv =
00555     new G4Tubs(itsName+"_solid_env",
00556                itsInnerIronRadius+BDSGlobalConstants::Instance()->GetLengthSafety()/2.0, // inner R + overlap safety
00557                itsOuterR,          // outer R
00558                tubLen,             // length
00559                0,                  // starting phi
00560                twopi * rad );      // delta phi
00561   
00562   G4IntersectionSolid *magTubs =
00563     new G4IntersectionSolid(itsName+"_solid",
00564                             magTubsEnv,
00565                             rbendRectangleSolidVolume,
00566                             RotYM90,
00567                             (G4ThreeVector)0); 
00568 
00569   if(OuterMaterialIsVacuum)
00570     {
00571       itsOuterLogicalVolume = 
00572         new G4LogicalVolume(magTubs,
00573                             BDSMaterials::Instance()->GetMaterial(BDSGlobalConstants::Instance()->GetVacuumMaterial()),
00574                             itsName+"_outer");
00575     }
00576   else
00577     {
00578       itsOuterLogicalVolume=
00579         new G4LogicalVolume(magTubs,
00580                             material,
00581                             itsName+"_outer");
00582     }
00583 
00584   itsPhysiComp =
00585     new G4PVPlacement(
00586                       RotY90,                 // rotation
00587                       (G4ThreeVector)0,                      // at (0,0,0)
00588                       itsOuterLogicalVolume,  // its logical volume
00589                       itsName+"_solid",       // its name
00590                       itsMarkerLogicalVolume, // its mother  volume
00591                       false,                  // no boolean operation
00592                       0, BDSGlobalConstants::Instance()->GetCheckOverlaps());                     // copy number
00593 
00594   SetMultiplePhysicalVolumes(itsPhysiComp);
00595 
00596 #ifndef NOUSERLIMITS
00597   itsOuterUserLimits =
00598     new G4UserLimits("multipole cut",DBL_MAX,DBL_MAX,BDSGlobalConstants::Instance()->GetMaxTime(),
00599                      BDSGlobalConstants::Instance()->GetThresholdCutCharged());
00600   itsOuterUserLimits->SetMaxAllowedStep(itsMagFieldLength);
00601   itsOuterLogicalVolume->SetUserLimits(itsOuterUserLimits);
00602 #endif
00603 }
00604 
00605 BDSRBend::~BDSRBend()
00606 {
00607   delete itsVisAttributes;
00608   delete innerBeampipeVisAtt;
00609   delete beampipeVisAtt;
00610 //   delete markerSolidVolume;
00611 //   delete rbendRectangleSolidVolume;
00612 //   delete rbendRectangleLogicalVolume;
00613 //   delete middleBeampipeLogicalVolume;
00614 //   delete middleInnerBPLogicalVolume;
00615 //   delete endsBeampipeLogicalVolume;
00616 //   delete endsInnerBPLogicalVolume;
00617 //   delete itsBeampipeUserLimits;
00618 //   delete endsBeampipeUserLimits;
00619 //   delete endsInnerBeampipeUserLimits;
00620   delete itsMagField;
00621   delete itsEqRhs;
00622   delete itsStepper;
00623 }

Generated on 27 Aug 2013 for BDSIM by  doxygen 1.4.7