src/BDSRegions.cc

00001 #include "G4Electron.hh"
00002 #include "G4Positron.hh"
00003 
00004 #include "BDSGlobalConstants.hh"
00005 #include "G4VUserDetectorConstruction.hh"
00006 #include "globals.hh"
00007 
00008 #include "BDSWorld.hh"
00009 #include "BDSMaterials.hh"
00010 #include "BDSBeamline.hh"
00011 
00012 #include "G4Region.hh"
00013 
00014 #include "G4IStore.hh"
00015 #include "G4GeometrySampler.hh"
00016 
00017 //GFlash parameterisation
00018 #include "GFlashHomoShowerParameterisation.hh"
00019 #include "G4FastSimulationManager.hh"
00020 #include "BDSShowerModel.hh"
00021 #include "GFlashHitMaker.hh"
00022 #include "GFlashParticleBounds.hh"
00023 #include "BDSRegions.hh"
00024 
00025 #include "G4UserLimits.hh"
00026 #include "G4Region.hh"
00027 #include "G4ProductionCuts.hh"
00028 
00029 #include "G4Tubs.hh"
00030 #include "G4Box.hh"
00031 #include "G4LogicalVolume.hh"
00032 #include "G4VPhysicalVolume.hh"
00033 #include "G4PVPlacement.hh"
00034 #include "G4UniformMagField.hh"
00035 #include "G4TransportationManager.hh"
00036 #include "G4PropagatorInField.hh"
00037 #include "G4SDManager.hh"
00038 #include "G4RunManager.hh"
00039 #include "G4ScoringBox.hh"
00040 #include "G4ScoringManager.hh"
00041 #include "G4PSCellFlux3D.hh"
00042 #include "BDSScoreWriter.hh"
00043 
00044 #include "G4VisAttributes.hh"
00045 #include "G4Colour.hh"
00046 #include "globals.hh"
00047 #include "G4ios.hh"
00048 #include <iostream>
00049 #include <list>
00050 #include <map>
00051 #include "BDSAcceleratorComponent.hh"
00052 
00053 #include "G4Navigator.hh"
00054 #include "G4UniformMagField.hh"
00055 
00056 #include "G4Material.hh"
00057 #include "BDSEnergyCounterSD.hh"
00058 
00059 extern G4int gflash;
00060 extern G4double gflashemax;
00061 extern G4double gflashemin;
00062 
00063 BDSRegions::BDSRegions(){
00064   buildRegions();
00065 }
00066 
00067 BDSRegions::~BDSRegions(){
00068   delete _precisionRegion;
00069   delete _worldRegion;
00070 }
00071 
00072 void BDSRegions::buildRegions(){
00073   buildPrecisionRegion();
00074   buildGFlashRegion();
00075 }
00076 
00077 void BDSRegions::buildPrecisionRegion(){
00078   _precisionRegion = new G4Region("precisionRegion");
00079   _precisionProductionCuts = new G4ProductionCuts();
00080   
00081   if(BDSGlobalConstants::Instance()->GetProdCutPhotonsP()>0)
00082     _precisionProductionCuts->SetProductionCut(BDSGlobalConstants::Instance()->GetProdCutPhotonsP(),G4ProductionCuts::GetIndex("gamma"));
00083   
00084   if(BDSGlobalConstants::Instance()->GetProdCutElectronsP()>0)
00085     _precisionProductionCuts->SetProductionCut(BDSGlobalConstants::Instance()->GetProdCutElectronsP(),G4ProductionCuts::GetIndex("e-"));
00086   
00087   if(BDSGlobalConstants::Instance()->GetProdCutPositronsP()>0)
00088     _precisionProductionCuts->SetProductionCut(BDSGlobalConstants::Instance()->GetProdCutPositronsP(),G4ProductionCuts::GetIndex("e+"));
00089   
00090   _precisionRegion->SetProductionCuts(_precisionProductionCuts);
00091 #ifndef NOUSERLIMITS
00092   //  _precisionRegion->SetUserLimits(BDSWorld::Instance()->userLimits());
00093 #endif
00094 }
00095 
00096 void BDSRegions::buildGFlashRegion(){
00097   _gFlashParticleBounds  = new GFlashParticleBounds();              // Energy Cuts to kill particles                                                                
00098   _gFlashParticleBounds->SetMaxEneToParametrise(*G4Electron::ElectronDefinition(),gflashemax*GeV);
00099   _gFlashParticleBounds->SetMinEneToParametrise(*G4Electron::ElectronDefinition(),gflashemin*GeV);
00100   _gFlashParticleBounds->SetEneToKill(*G4Electron::ElectronDefinition(),BDSGlobalConstants::Instance()->GetThresholdCutCharged());
00101   
00102   _gFlashParticleBounds->SetMaxEneToParametrise(*G4Positron::PositronDefinition(),gflashemax*GeV);
00103   _gFlashParticleBounds->SetMinEneToParametrise(*G4Positron::PositronDefinition(),gflashemin*GeV);
00104   _gFlashParticleBounds->SetEneToKill(*G4Positron::PositronDefinition(),BDSGlobalConstants::Instance()->GetThresholdCutCharged());
00105 
00106   _gFlashParticleBoundsVac  = new GFlashParticleBounds();              // Energy Cuts to kill particles                                                                
00107   _gFlashParticleBoundsVac->SetMaxEneToParametrise(*G4Electron::ElectronDefinition(),0*GeV);
00108   _gFlashParticleBoundsVac->SetMaxEneToParametrise(*G4Positron::PositronDefinition(),0*GeV);
00109 
00110   G4cout << "BDSRegions:buildGFlashRegion() -  _gFlashParticleBounds - min E - electron: " << _gFlashParticleBounds->GetMinEneToParametrise(*G4Electron::ElectronDefinition())/GeV<< " GeV" << G4endl;
00111   G4cout << "BDSRegions:buildGFlashRegion() -  _gFlashParticleBounds - max E - electron: " << _gFlashParticleBounds->GetMaxEneToParametrise(*G4Electron::ElectronDefinition())/GeV<< G4endl;
00112   G4cout << "BDSRegions:buildGFlashRegion() -  _gFlashParticleBounds - kill E - electron: " << _gFlashParticleBounds->GetEneToKill(*G4Electron::ElectronDefinition())/GeV<< G4endl;
00113 
00114   G4cout << "BDSRegions:buildGFlashRegion() -  _gFlashParticleBounds - min E - positron: " << _gFlashParticleBounds->GetMinEneToParametrise(*G4Positron::PositronDefinition())/GeV<< G4endl;
00115   G4cout << "BDSRegions:buildGFlashRegion() -  _gFlashParticleBounds - max E - positron: " << _gFlashParticleBounds->GetMaxEneToParametrise(*G4Positron::PositronDefinition())/GeV<< G4endl;
00116   G4cout << "BDSRegions:buildGFlashRegion() -  _gFlashParticleBounds - kill E - positron: " << _gFlashParticleBounds->GetEneToKill(*G4Positron::PositronDefinition())/GeV<< G4endl;
00117 
00118  _gFlashHitMaker = new GFlashHitMaker();                    // Makes the EnergieSpots 
00119 }
00120 
00121 void BDSRegions::buildGFlashRegion(BDSAcceleratorComponent* /*var*/){
00122 
00123   /*
00124   vector<G4LogicalVolume*> MultipleSensVols = var->GetMultipleSensitiveVolumes();
00125   if( ( var->GetType()!="sampler" && var->GetType()!="csampler" )
00126       && MultipleSensVols.size()>0)
00127     {
00128       for(G4int i=0; i<(G4int)MultipleSensVols.size(); i++)
00129         {
00130           if((MultipleSensVols.at(i)->GetRegion() != precisionRegion) && (var->GetType()==_ELEMENT)){//If not in the precision region....
00131             //      G4cout << "...adding " << MultipleSensVols[i]->GetName() << " to gFlashRegion" << G4endl;
00132             G4String rname = "gFlashRegion_" + MultipleSensVols[i]->GetName();
00133             gFlashRegion.push_back(new G4Region(rname.c_str()));
00134             G4String mname = "fastShowerModel" + rname;
00135             G4cout << "...making parameterisation..." << G4endl;
00136                   _gFlashFastShowerModel.push_back(new BDSShowerModel(mname.c_str(),gFlashRegion.back()));
00137                   _gFlashParameterisation.push_back(new GFlashHomoShowerParameterisation(BDSMaterials::Instance()->GetMaterial(MultipleSensVols[i]->GetMaterial()->GetName().c_str()))); 
00138                   _gFlashFastShowerModel.back()->SetParameterisation(*_gFlashParameterisation.back());
00139                   _gFlashFastShowerModel.back()->SetParticleBounds(*_gFlashParticleBounds) ;
00140                   _gFlashFastShowerModel.back()->SetHitMaker(*_gFlashHitMaker);
00141                   if(MultipleSensVols[i]->GetMaterial()->GetState()!=kStateGas){ //If the region material state is not gas, associate with a parameterisation
00142                     _gFlashFastShowerModel.back()->SetFlagParamType(1);//Turn on the parameterisation for e-m showers starting in sensitive material and fitting in the current stack.
00143                     _gFlashFastShowerModel.back()->SetFlagParticleContainment(1);//Turn on containment
00144                   } else {
00145                     _gFlashFastShowerModel.back()->SetFlagParamType(0);
00146                     _gFlashFastShowerModel.back()->SetFlagParticleContainment(0);
00147                     
00148                   }
00149                   MultipleSensVols[i]->SetRegion(gFlashRegion.back());
00150                   gFlashRegion.back()->AddRootLogicalVolume(MultipleSensVols[i]);
00151           }               
00152         }
00153     }
00154 */
00155 }
00156 
00157 
00158 void BDSRegions::buildGasRegion(){
00159   _gasRegion = new G4Region("gasRegion");
00160   _gasProductionCuts = new G4ProductionCuts();
00161   _gasProductionCuts->SetProductionCut(1*m,G4ProductionCuts::GetIndex("gamma"));
00162   _gasProductionCuts->SetProductionCut(1*m,G4ProductionCuts::GetIndex("e-"));
00163   _gasProductionCuts->SetProductionCut(1*m,G4ProductionCuts::GetIndex("e+"));
00164   _gasRegion->SetProductionCuts(_gasProductionCuts);
00165 }

Generated on 27 Aug 2013 for BDSIM by  doxygen 1.4.7