Public Member Functions | Public Attributes | Protected Attributes | List of all members
LBNERunManager Class Reference

#include <LBNERunManager.hh>

Inheritance diagram for LBNERunManager:

Public Member Functions

 LBNERunManager ()
 
virtual ~LBNERunManager ()
 
virtual void InitializeGeometry ()
 
virtual void InitializePhysics ()
 
virtual void BeamOn (G4int n_event, const char *macroFile=0, G4int n_select=-1)
 
G4int GetNumberOfEvents ()
 
void SetNptInputFileName (G4String &aName)
 
LBNEPrimaryGeneratorActionGetLBNEPrimaryGeneratorAction ()
 
LBNESteppingActionGetLBNESteppingManager ()
 
void SetUseFlukaInput (bool t)
 
void SetUseMarsInput (bool t)
 
bool GetUseFlukaInput () const
 
bool GetUseMarsInput () const
 
G4String GetNptInputFileName () const
 
void SetOutputNtpFileName (G4String &aName)
 
G4String GetOutputNtpFileName () const
 
void SetOutputDk2NuFileName (G4String &aName)
 
G4String GetOutputDk2NuFileName () const
 
void SetOutputASCIIFileName (G4String &aName)
 
G4String GetOutputASCIIFileName () const
 
bool GetCreateOutput () const
 
void SetCreateOutput (bool t)
 
bool GetCreateDk2NuOutput () const
 
void SetCreateDk2NuOutput (bool t)
 
bool GetCreateASCIIOutput () const
 
void SetCreateASCIIOutput (bool t)
 
bool GetCreateTrkPlaneOutput () const
 
void SetCreateTrkPlaneOutput (bool t)
 
bool GetUseRealisticNearDetectorVolume () const
 
void SetUseRealisticNearDetectorVolume (bool t)
 
bool GetCreateTrackingTargetOutASCII () const
 
void SetCreateTrackingTargetOutASCII (bool t)
 
G4String GetPhysicsListName () const
 
void SetPhysicsListName (G4String &t)
 
G4String GetDetectorLocationFileName () const
 
void SetDetectorLocationFileName (G4String &sDetLocName)
 
bool GetCreateTargetOutput () const
 
void SetCreateTargetOutput (bool t)
 
bool GetCreateTrkPlaneH1Output () const
 
bool GetCreateTrkPlaneH2Output () const
 
bool GetCreateTrkPlaneDPOutput () const
 
void SetCreateTrkPlaneH1Output (bool t)
 
void SetCreateTrkPlaneH2Output (bool t)
 
void SetCreateTrkPlaneDPOutput (bool t)
 
bool GetCreateAlcoveTrackingOutput () const
 
void SetCreateAlcoveTrackingOutput (bool t)
 
void SetNumberOfEventsLBNE (int n)
 
int GetNumberOfEventsLBNE ()
 
int GetVerboseLevel () const
 

Public Attributes

G4int nEvents
 

Protected Attributes

bool fGeometryIntializedHere
 
bool fPhysicsInitializedHere
 
G4String fPhysicsListName
 
LBNEEventActionfLBNEEventAction
 
LBNESteppingActionfLBNESteppingAction
 
LBNEStackingActionfLBNEStackingAction
 
LBNETrackingActionfLBNETrackingAction
 
LBNERunActionfLBNERunAction
 
G4String fMarsOrFlukaInputFileName
 
G4String fAsciiOutputFileName
 
G4String fNptOutputFileName
 
G4String fDk2NuOutputFileName
 
G4String fDetectorLocationFileName
 
bool fUseFluka
 
bool fUseMars
 
bool fCreateOutput
 
bool fCreateDk2NuOutput
 
bool fCreateAsciiOutput
 
bool fCreateTrkPlaneOutput
 
bool fCreateTrkPlaneH1Output
 
bool fCreateTrkPlaneH2Output
 
bool fCreateTrkPlaneDPOutput
 
bool fCreateTrackingOutput
 
bool fCreateTrackingTargetOutASCII
 
bool fCreateTargetOutput
 
bool fCreateAlcoveTrackingOutput
 
bool fUseRealisticNearDetectorVolume
 
int fVerboseLevel
 

Detailed Description

Definition at line 17 of file LBNERunManager.hh.

Constructor & Destructor Documentation

LBNERunManager::LBNERunManager ( )

Definition at line 30 of file LBNERunManager.cc.

30  :
31 nEvents(50000),
34 fPhysicsListName("QGSP_BERT"),
45 fUseFluka(false),
46 fUseMars(false),
47 fCreateOutput(true),
48 fCreateDk2NuOutput(false),
49 fCreateAsciiOutput(false),
51 fCreateTrkPlaneH1Output(false), //Amit Bashyal To turn on and off the Horn 1 output
52 fCreateTrkPlaneH2Output(false), //Amit Bashyal To turn on and off the Horn 2 output
53 fCreateTrkPlaneDPOutput(false), //Amit Bashyal To turn on and off the Decay Pipe output
55 fCreateTargetOutput(false), //Amit Bashyal To turn on and off the Target Output
58 
59 //:primaryGeneratorAction(0)
60 {
62 
63  if(fVerboseLevel > 0)
64  {
65  std::cout << "LBNERunManager Constructor Called." << std::endl;
66  }
67 }
bool fPhysicsInitializedHere
bool fCreateAlcoveTrackingOutput
bool fCreateTrkPlaneOutput
LBNEStackingAction * fLBNEStackingAction
G4String fNptOutputFileName
bool fCreateTrackingTargetOutASCII
bool fCreateTrkPlaneDPOutput
LBNERunAction * fLBNERunAction
bool fCreateTrkPlaneH2Output
LBNETrackingAction * fLBNETrackingAction
G4String fMarsOrFlukaInputFileName
bool fGeometryIntializedHere
bool fCreateTrkPlaneH1Output
G4String fDetectorLocationFileName
bool fUseRealisticNearDetectorVolume
LBNEEventAction * fLBNEEventAction
LBNESteppingAction * fLBNESteppingAction
G4String fPhysicsListName
G4String fAsciiOutputFileName
G4String fDk2NuOutputFileName
QTextStream & endl(QTextStream &s)
LBNERunManager::~LBNERunManager ( )
virtual

Definition at line 69 of file LBNERunManager.cc.

70 {
71  if(fVerboseLevel > 0)
72  {
73  std::cout << "LBNERunManager Destructor Called." << std::endl;
74  }
75 
76 }
QTextStream & endl(QTextStream &s)

Member Function Documentation

void LBNERunManager::BeamOn ( G4int  n_event,
const char *  macroFile = 0,
G4int  n_select = -1 
)
virtual

Definition at line 332 of file LBNERunManager.cc.

333 {
334  G4bool runOn(true);
335  G4bool cond = ConfirmBeamOnCondition();
336  if(cond)
337  {
338  if(runOn)
339  {
340  TStopwatch *sWatch=new TStopwatch();
341  sWatch->Start(true);
342  // Over write the number of events with our own
343  numberOfEventToBeProcessed = nEvents;
344 
345  //
346  //RunInitializtion() calls a lot of functions one of which is
347  //LBNERunAction::BeginOfRunAction(). Put all run Initiliation code here
348  //like opening the Output Ntuple and checking to make sure that
349  //the simulation is setup correctly and it is OK to Run.
350  //
351  G4cout << G4endl;
352  std::cout << "********************************************************" << std::endl;
353  std::cout << "********************************************************" << std::endl;
354  std::cout << "*****LBNERunManager::BeamOn - With "<<nEvents << " Events " << std::endl;
355  std::cout << "***** (argument n_event = " << n_event <<" ignored, O.K.) " << std::endl;
356  std::cout << "********************************************************" << std::endl;
357  G4RunManager::BeamOn(nEvents, macroFile, n_select);
358  std::cout << "********************************************************" << std::endl;
359  std::cout << "***LBNERunManager::BeamOn - ...Done ***" << std::endl;
360  std::cout << "********************************************************" << std::endl;
361  std::cout << "********************************************************" << std::endl;
362  G4cout << G4endl;
363  sWatch->Stop();
364  G4cout << "LBNERunManager::BeamOn - Run Summary..." << G4endl;
365  G4cout << " Processed "<< nEvents << " beam particles in ";
366  G4cout.precision(3);
367  G4cout << sWatch->RealTime() << " s = "
368  << sWatch->RealTime()/60. << " min "
369  << sWatch->RealTime()/3600. << " hr " << G4endl;
370  delete sWatch;
371 
372  }
373  }
374 }
QTextStream & endl(QTextStream &s)
bool LBNERunManager::GetCreateAlcoveTrackingOutput ( ) const
inline

Definition at line 128 of file LBNERunManager.hh.

bool fCreateAlcoveTrackingOutput
bool LBNERunManager::GetCreateASCIIOutput ( ) const
inline

Definition at line 95 of file LBNERunManager.hh.

95 { return fCreateAsciiOutput; }
bool LBNERunManager::GetCreateDk2NuOutput ( ) const
inline

Definition at line 93 of file LBNERunManager.hh.

93 { return fCreateDk2NuOutput; }
bool LBNERunManager::GetCreateOutput ( ) const
inline

Definition at line 91 of file LBNERunManager.hh.

91 { return fCreateOutput; }
bool LBNERunManager::GetCreateTargetOutput ( ) const
inline

Definition at line 120 of file LBNERunManager.hh.

120 {return fCreateTargetOutput;}
bool LBNERunManager::GetCreateTrackingTargetOutASCII ( ) const
inline

Definition at line 101 of file LBNERunManager.hh.

bool fCreateTrackingTargetOutASCII
bool LBNERunManager::GetCreateTrkPlaneDPOutput ( ) const
inline

Definition at line 124 of file LBNERunManager.hh.

124 {return fCreateTrkPlaneDPOutput;}
bool fCreateTrkPlaneDPOutput
bool LBNERunManager::GetCreateTrkPlaneH1Output ( ) const
inline

Definition at line 122 of file LBNERunManager.hh.

122 {return fCreateTrkPlaneH1Output;}
bool fCreateTrkPlaneH1Output
bool LBNERunManager::GetCreateTrkPlaneH2Output ( ) const
inline

Definition at line 123 of file LBNERunManager.hh.

123 {return fCreateTrkPlaneH2Output;}
bool fCreateTrkPlaneH2Output
bool LBNERunManager::GetCreateTrkPlaneOutput ( ) const
inline

Definition at line 97 of file LBNERunManager.hh.

97 { return fCreateTrkPlaneOutput; }
bool fCreateTrkPlaneOutput
G4String LBNERunManager::GetDetectorLocationFileName ( ) const
inline

Definition at line 118 of file LBNERunManager.hh.

118 { return fDetectorLocationFileName; }
G4String fDetectorLocationFileName
LBNEPrimaryGeneratorAction* LBNERunManager::GetLBNEPrimaryGeneratorAction ( )
inline

Definition at line 71 of file LBNERunManager.hh.

72  { return dynamic_cast<LBNEPrimaryGeneratorAction*>(userPrimaryGeneratorAction); }
LBNESteppingAction* LBNERunManager::GetLBNESteppingManager ( )
inline

Definition at line 74 of file LBNERunManager.hh.

75  { return dynamic_cast<LBNESteppingAction*> (userSteppingAction); }
G4String LBNERunManager::GetNptInputFileName ( ) const
inline

Definition at line 82 of file LBNERunManager.hh.

G4String fMarsOrFlukaInputFileName
G4int LBNERunManager::GetNumberOfEvents ( )
inline

Definition at line 26 of file LBNERunManager.hh.

26  {
27  return numberOfEventToBeProcessed;
28  }
int LBNERunManager::GetNumberOfEventsLBNE ( )
inline

Definition at line 133 of file LBNERunManager.hh.

133 { return nEvents;}
G4String LBNERunManager::GetOutputASCIIFileName ( ) const
inline

Definition at line 89 of file LBNERunManager.hh.

89 {return fAsciiOutputFileName;}
G4String fAsciiOutputFileName
G4String LBNERunManager::GetOutputDk2NuFileName ( ) const
inline

Definition at line 87 of file LBNERunManager.hh.

87 {return fDk2NuOutputFileName;}
G4String fDk2NuOutputFileName
G4String LBNERunManager::GetOutputNtpFileName ( ) const
inline

Definition at line 85 of file LBNERunManager.hh.

85 {return fNptOutputFileName;}
G4String fNptOutputFileName
G4String LBNERunManager::GetPhysicsListName ( ) const
inline

Definition at line 116 of file LBNERunManager.hh.

116 {return fPhysicsListName; }
G4String fPhysicsListName
bool LBNERunManager::GetUseFlukaInput ( ) const
inline

Definition at line 80 of file LBNERunManager.hh.

80 {return fUseFluka;}
bool LBNERunManager::GetUseMarsInput ( ) const
inline

Definition at line 81 of file LBNERunManager.hh.

81 {return fUseMars;}
bool LBNERunManager::GetUseRealisticNearDetectorVolume ( ) const
inline

Definition at line 99 of file LBNERunManager.hh.

bool fUseRealisticNearDetectorVolume
int LBNERunManager::GetVerboseLevel ( ) const
inline

Definition at line 135 of file LBNERunManager.hh.

135 {return fVerboseLevel; }
void LBNERunManager::InitializeGeometry ( )
virtual

Definition at line 79 of file LBNERunManager.cc.

79  {
80  if (fGeometryIntializedHere) return;
81  std::cout <<
82  "LBNERunManager InitializeGeometry: LBNE Beam line should be already constructed... " << std::endl;
83 // G4GeometryManager::GetInstance()->CloseGeometry(); //
84  G4RunManager::InitializeGeometry();
85  std::cout << " ... Geometry has been closed. keep going... " << std::endl;
86  std::cout << " Geant4 is at stage " << G4StateManager::GetStateManager()->GetCurrentState() << std::endl;
88 }
bool fGeometryIntializedHere
QTextStream & endl(QTextStream &s)
void LBNERunManager::InitializePhysics ( )
virtual

Definition at line 90 of file LBNERunManager.cc.

90  {
91  std::cout<<G4State_PreInit<<"\t"<<G4State_Init<<"\t"<<G4State_Idle<<G4endl;
92  std::cout << " InitPhys::Geant4 is at stage " << G4StateManager::GetStateManager()->GetCurrentState() << std::endl;
93  if (fPhysicsInitializedHere) return;
94 
95  //in g4.10.3 new state was introduced (see release notes: http://geant4-data.web.cern.ch/geant4-data/ReleaseNotes/ReleaseNotes4.10.3.html
96  //code crashes here with g4.10.3
97  //with previous version of geant4 stage was PreInit (0) here, while now it is Init
98  std::cout << " Before::Geant4 is at stage " << G4StateManager::GetStateManager()->GetCurrentState() << std::endl;
99  G4StateManager* stateManager = G4StateManager::GetStateManager();
100  stateManager->SetNewState(G4State_PreInit);
101 
102  G4PhysListFactory factory;
103  G4VModularPhysicsList *phys = 0;
104  if(factory.IsReferencePhysList(fPhysicsListName)){
105  phys = factory.GetReferencePhysList(fPhysicsListName);
106  } else {
107  G4cout << "Couldn't find physics list name " << fPhysicsListName << G4endl;
108  G4cout << "Exiting" << G4endl;
109  exit(0);
110  }
111  this->SetUserInitialization(phys);
112  this->SetUserAction(new LBNEPrimaryGeneratorAction);
113 
114  G4RunManager::InitializePhysics();
115 
116  std::cout << " After::Geant4 is at stage " << G4StateManager::GetStateManager()->GetCurrentState() << std::endl;
117  G4UImanager* UI = G4UImanager::GetUIpointer();
118 // ---> April 2006. D.Jaffe modifications to set K0L,K+,K- branching fractions
119 
120  // K0L branching fractions from Table 3 of hep-ex/0512039
121  // "Status of the Cabibbo Angle", E. Blucher et al.,
122  // Proceedings of the CKM 2005 Workshop (WG1), UC San Diego, 15-18 March 2005
123  // Average takes into account KTeV, KLOE, NA48 results
124  // KTeV PRL 93, 181802. PR D70 092006
125  // KLOE hep-ex/0408027, PL B632 (2006) 43
126  // NA48 PL B602 (2004) 41
127 
128  // K+ branching fractions from PDG2005 update and section 3.2.1 of hep-ex/0512039
129  // for average Ke3 from BNL E865 (PRL 91, 261802) and NA48.
130  // B(Ke3) = 5.14+-0.06 % hep-ex/0512039
131  // B(Ke3) = 4.93+-0.07 % PDG2005 (takes into account BNL and some much older results)
132 
133  G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
134  if (!particleTable->contains("kaon0L")) {
135  std::ostringstream mStrStr;
136  mStrStr << " The Particle table does not know about kaon Long.. Something not initialized ";
137  G4String mStr(mStrStr.str());
138  G4Exception("LBNERunManager::InitializePhysics", " ", FatalException, mStr.c_str());
139  }
140  G4cout<<G4endl;
141  G4cout << "Setting correct branching fractions and form factors for "
142  << "K0L, K+ and K- decays (MINOS-doc-1786)"<<G4endl;
143 
144  UI->ApplyCommand("/particle/select kaon0L"); // ----------- K0L
145  UI->ApplyCommand("/particle/property/decay/select 0"); //3pi0
146  UI->ApplyCommand("/particle/property/decay/br 19.72e-2");
147  UI->ApplyCommand("/particle/property/decay/select 1"); //pi-,e+,nue
148  UI->ApplyCommand("/particle/property/decay/br 20.20e-2");
149  UI->ApplyCommand("/particle/property/decay/select 2"); //pi+,e-,nue
150  UI->ApplyCommand("/particle/property/decay/br 20.20e-2");
151  UI->ApplyCommand("/particle/property/decay/select 3"); //pi-,mu+,nue
152  UI->ApplyCommand("/particle/property/decay/br 13.495e-2");
153  UI->ApplyCommand("/particle/property/decay/select 4"); //pi+,mu-,nue
154  UI->ApplyCommand("/particle/property/decay/br 13.495e-2");
155  UI->ApplyCommand("/particle/property/decay/select 5"); //pi+,pi-,pi0
156  UI->ApplyCommand("/particle/property/decay/br 12.53e-2");
157 
158  UI->ApplyCommand("/particle/select kaon+"); // ------------- K+
159  UI->ApplyCommand("/particle/property/decay/select 0"); //mu,nu
160  UI->ApplyCommand("/particle/property/decay/br 63.39e-2");
161  UI->ApplyCommand("/particle/property/decay/select 1"); //pi,pi0
162  UI->ApplyCommand("/particle/property/decay/br 21.03e-2");
163  UI->ApplyCommand("/particle/property/decay/select 2"); //3pi
164  UI->ApplyCommand("/particle/property/decay/br 5.59e-2");
165  UI->ApplyCommand("/particle/property/decay/select 3"); //pi0,e,nu
166  UI->ApplyCommand("/particle/property/decay/br 5.14e-2");
167  UI->ApplyCommand("/particle/property/decay/select 4"); //pi0,mu,nu
168  UI->ApplyCommand("/particle/property/decay/br 3.30e-2");
169  UI->ApplyCommand("/particle/property/decay/select 5"); //pi0,pi0,pi
170  UI->ApplyCommand("/particle/property/decay/br 1.757e-2");
171 
172  UI->ApplyCommand("/particle/select kaon-"); // ------------- K-
173  UI->ApplyCommand("/particle/property/decay/select 0"); //mu,nu
174  UI->ApplyCommand("/particle/property/decay/br 63.39e-2");
175  UI->ApplyCommand("/particle/property/decay/select 1"); //pi,pi0
176  UI->ApplyCommand("/particle/property/decay/br 21.03e-2");
177  UI->ApplyCommand("/particle/property/decay/select 2"); //3pi
178  UI->ApplyCommand("/particle/property/decay/br 5.59e-2");
179  UI->ApplyCommand("/particle/property/decay/select 3"); //pi0,e,nu
180  UI->ApplyCommand("/particle/property/decay/br 5.14e-2");
181  UI->ApplyCommand("/particle/property/decay/select 4"); //pi0,mu,nu
182  UI->ApplyCommand("/particle/property/decay/br 3.30e-2");
183  UI->ApplyCommand("/particle/property/decay/select 5"); //pi0,pi0,pi
184  UI->ApplyCommand("/particle/property/decay/br 1.757e-2");
185  // <--- end of djaffe april06 additions
186 
187  // ---> April 2006. Adding Davids modifications here instead of changing
188  // G4KL3DecayChannel.cc
189  G4ParticleTable::G4PTblDicIterator *ptiter=particleTable->GetIterator();
190 
191  G4double lamp = 2.81e-2;
192  G4double lam0 = 1.64e-2;
193  G4double xi_kp= 12.377 * (lam0 - lamp);
194  G4double xi_k0= 11.713 * (lam0 - lamp);
195  G4double pLambda=0.;
196  G4double pXi0=0.;
197 
198  ptiter->reset();
199  // G4cout <<G4endl;
200  // G4cout <<" Decay mode \t \t \t \t Lambda \t \t \t Xi0 "<<G4endl;
201  while ((*ptiter)()){
202  G4ParticleDefinition *particle = ptiter->value();
203  if (particle==G4KaonPlus::KaonPlusDefinition()||
204  particle==G4KaonMinus::KaonMinusDefinition()){
205  G4DecayTable *decayTable=particle->GetDecayTable();
206  for (G4int ii=0; ii<decayTable->entries(); ++ii){
207  if (decayTable->GetDecayChannel(ii)->GetKinematicsName()=="KL3 Decay"){
208  G4KL3DecayChannel *decayChannel=(G4KL3DecayChannel *)decayTable->GetDecayChannel(ii);
209 
210  //G4cout<<particle->GetParticleName();
211  for (G4int jj=0; jj<decayChannel->GetNumberOfDaughters();++jj){
212  // if (jj==0) G4cout << " -> ";
213  // else G4cout << " + ";
214  // G4cout<<decayChannel->GetDaughter(jj)->GetParticleName();
215  if (decayChannel->GetDaughter(jj)==G4Positron::PositronDefinition()){
216  pLambda=lamp;
217  pXi0=xi_kp;
218  }
219  if (decayChannel->GetDaughter(jj)==G4Electron::ElectronDefinition()){
220  pLambda=lamp;
221  pXi0=xi_kp;
222  }
223  if (decayChannel->GetDaughter(jj)==G4MuonPlus::MuonPlusDefinition()){
224  pLambda=lamp;
225  pXi0=xi_kp;
226  }
227  if (decayChannel->GetDaughter(jj)==G4MuonMinus::MuonMinusDefinition()){
228  pLambda=lamp;
229  pXi0=xi_kp;
230  }
231  }//loop over daughters
232  // G4cout<<"\t \t" << decayChannel->GetDalitzParameterLambda() << " -> " << pLambda<< " "
233  // <<"\t"<< decayChannel->GetDalitzParameterXi() << " -> " << pXi0 <<G4endl ;
234  decayChannel->SetDalitzParameter(pLambda,pXi0);
235  }//check if its KL3 Decay
236  }//loop over decay channels
237  }
238  if (particle==G4KaonZeroLong::KaonZeroLongDefinition()){
239  G4DecayTable *decayTable=particle->GetDecayTable();
240  for (G4int ii=0; ii<decayTable->entries(); ++ii){
241  if (decayTable->GetDecayChannel(ii)->GetKinematicsName()=="KL3 Decay"){
242  G4KL3DecayChannel *decayChannel=(G4KL3DecayChannel *)decayTable->GetDecayChannel(ii);
243  // G4cout<<particle->GetParticleName();
244  for (G4int jj=0; jj<decayChannel->GetNumberOfDaughters();++jj){
245  // if (jj==0) G4cout << " -> ";
246  // else G4cout << " + ";
247  // G4cout<<decayChannel->GetDaughter(jj)->GetParticleName();
248  if (decayChannel->GetDaughter(jj)==G4Positron::PositronDefinition()||
249  decayChannel->GetDaughter(jj)==G4Electron::ElectronDefinition()){
250  pLambda=lamp;
251  pXi0=xi_k0;
252  }
253  if (decayChannel->GetDaughter(jj)==G4MuonPlus::MuonPlusDefinition()||
254  decayChannel->GetDaughter(jj)==G4MuonMinus::MuonMinusDefinition()){
255  pLambda=lamp;
256  pXi0=xi_k0;
257  }
258  }//loop over daughters
259  // G4cout<<"\t \t" << decayChannel->GetDalitzParameterLambda() << " -> " << pLambda<< " "
260  // <<"\t"<< decayChannel->GetDalitzParameterXi() << " -> " << pXi0 <<G4endl ;
261 
262  decayChannel->SetDalitzParameter(pLambda,pXi0);
263 
264  }//check if its KL3 Decay
265  }
266  }
267  }//loop over particles
268  // done with kaon form factors
269  // Setting the maximum step size to 1 cm
270  G4ParticleTable* ptbl = G4ParticleTable::GetParticleTable();
271  G4ParticleTable::G4PTblDicIterator* piter = ptbl->GetIterator();
272  G4StepLimiter* slim = new G4StepLimiter("StepLimiter");
273 
274  piter->reset();
275  while ( (*piter)() ) {
276  G4ParticleDefinition* pdef = piter->value();
277  G4ProcessManager* pmgr = pdef->GetProcessManager();
278 
279  // add user limit processes for steps and special cuts
280  if ( pmgr ) {
281  pmgr->AddProcess( slim, -1, -1, 3);
282  }
283  }
285  G4RunManager::Initialize();
286 
287  std::cerr << " G4RunManager::Initialize invoked, state is now "
288  << G4StateManager::GetStateManager()->GetCurrentState() << std::endl;
289 
291  G4RunManager::SetUserAction(fLBNEEventAction);
292 
294  G4RunManager::SetUserAction(fLBNESteppingAction);
295 
297  G4RunManager::SetUserAction(fLBNETrackingAction);
298 
300  G4RunManager::SetUserAction(fLBNERunAction);
301 
303  G4RunManager::SetUserAction(fLBNEStackingAction);
304 
305  if (fVerboseLevel > 0) std::cerr << " All User Actions created and declared " << std::endl;
306 
307  if (!G4RunManager::ConfirmBeamOnCondition()) {
308  std::ostringstream mStrStr;
309  mStrStr << " Expecting to be able to turn on the beam at this stage. Not the case ! Fatal ";
310  G4String mStr(mStrStr.str());
311  G4Exception("LBNERunManager::InitializePhysics", " ", FatalException, mStr.c_str());
312  } else {
313  std::cout << " Expecting to be able to turn on the beam at this stage " << std::endl;
314  }
315  // Implementation of Perfect Focusing: as a genuiune G4Process.
316 //
317  G4ParticleDefinition *piPlus = particleTable->FindParticle(211);
319  G4ProcessManager *piPlusProcM = piPlus->GetProcessManager();
320  piPlusProcM->AddDiscreteProcess(pfProc);
321  G4ParticleDefinition *piMinus = particleTable->FindAntiParticle(211);
322  G4ProcessManager *piMinusProcM = piMinus->GetProcessManager();
323  piMinusProcM->AddDiscreteProcess(pfProc);
324  G4ParticleDefinition *kPlus = particleTable->FindParticle(321);
325  G4ProcessManager *kPlusProcM = kPlus->GetProcessManager();
326  kPlusProcM->AddDiscreteProcess(pfProc);
327  G4ParticleDefinition *kMinus = particleTable->FindAntiParticle(321);
328  G4ProcessManager *kMinusProcM = kMinus->GetProcessManager();
329  kMinusProcM->AddDiscreteProcess(pfProc);
330 
331  }
bool fPhysicsInitializedHere
LBNEStackingAction * fLBNEStackingAction
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:45
LBNERunAction * fLBNERunAction
LBNETrackingAction * fLBNETrackingAction
LBNEEventAction * fLBNEEventAction
LBNESteppingAction * fLBNESteppingAction
G4String fPhysicsListName
QTextStream & endl(QTextStream &s)
void LBNERunManager::SetCreateAlcoveTrackingOutput ( bool  t)
inline

Definition at line 129 of file LBNERunManager.hh.

void LBNERunManager::SetCreateASCIIOutput ( bool  t)
inline

Definition at line 96 of file LBNERunManager.hh.

96 { fCreateAsciiOutput = t;}
void LBNERunManager::SetCreateDk2NuOutput ( bool  t)
inline

Definition at line 94 of file LBNERunManager.hh.

94 { fCreateDk2NuOutput = t;}
void LBNERunManager::SetCreateOutput ( bool  t)
inline

Definition at line 92 of file LBNERunManager.hh.

92 { fCreateOutput = t;}
void LBNERunManager::SetCreateTargetOutput ( bool  t)
inline

Definition at line 121 of file LBNERunManager.hh.

121 {fCreateTargetOutput = t;}
void LBNERunManager::SetCreateTrackingTargetOutASCII ( bool  t)
inline

Definition at line 102 of file LBNERunManager.hh.

102  {
104  if (t) {
105  const LBNESteppingAction *stepAct
106  = reinterpret_cast<const LBNESteppingAction *>(G4RunManager::GetUserSteppingAction());
108  const LBNETrackingAction *trackAct
109  = reinterpret_cast<const LBNETrackingAction *>(G4RunManager::GetUserTrackingAction());
110  trackAct->OpenHadronAtVertex();
111  const LBNEStackingAction *stackAct
112  = reinterpret_cast<const LBNEStackingAction *>(G4RunManager::GetUserStackingAction());
113  stackAct->OpenHadronAtVertex();
114  }
115  }
void OpenHadronAtVertex() const
bool fCreateTrackingTargetOutASCII
void InitiateHadronFluxFromTargetASCII() const
void OpenHadronAtVertex() const
void LBNERunManager::SetCreateTrkPlaneDPOutput ( bool  t)
inline

Definition at line 127 of file LBNERunManager.hh.

void LBNERunManager::SetCreateTrkPlaneH1Output ( bool  t)
inline

Definition at line 125 of file LBNERunManager.hh.

bool fCreateTrkPlaneH1Output
void LBNERunManager::SetCreateTrkPlaneH2Output ( bool  t)
inline

Definition at line 126 of file LBNERunManager.hh.

void LBNERunManager::SetCreateTrkPlaneOutput ( bool  t)
inline

Definition at line 98 of file LBNERunManager.hh.

bool fCreateTrkPlaneOutput
void LBNERunManager::SetDetectorLocationFileName ( G4String &  sDetLocName)
inline

Definition at line 119 of file LBNERunManager.hh.

119 { fDetectorLocationFileName = sDetLocName;}
G4String fDetectorLocationFileName
void LBNERunManager::SetNptInputFileName ( G4String &  aName)
inline

Definition at line 69 of file LBNERunManager.hh.

69 { fMarsOrFlukaInputFileName = aName; }
G4String fMarsOrFlukaInputFileName
void LBNERunManager::SetNumberOfEventsLBNE ( int  n)
inline

Definition at line 132 of file LBNERunManager.hh.

132 {nEvents = n;}
void LBNERunManager::SetOutputASCIIFileName ( G4String &  aName)
inline

Definition at line 88 of file LBNERunManager.hh.

88 {fAsciiOutputFileName = aName;}
G4String fAsciiOutputFileName
void LBNERunManager::SetOutputDk2NuFileName ( G4String &  aName)
inline

Definition at line 86 of file LBNERunManager.hh.

86 {fDk2NuOutputFileName = aName;}
G4String fDk2NuOutputFileName
void LBNERunManager::SetOutputNtpFileName ( G4String &  aName)
inline

Definition at line 84 of file LBNERunManager.hh.

84 {fNptOutputFileName = aName;}
G4String fNptOutputFileName
void LBNERunManager::SetPhysicsListName ( G4String &  t)
inline

Definition at line 117 of file LBNERunManager.hh.

117 {fPhysicsListName = t;}
G4String fPhysicsListName
void LBNERunManager::SetUseFlukaInput ( bool  t)
inline

Definition at line 77 of file LBNERunManager.hh.

77 { fUseFluka = t; }
void LBNERunManager::SetUseMarsInput ( bool  t)
inline

Definition at line 78 of file LBNERunManager.hh.

78 { fUseMars = t; }
void LBNERunManager::SetUseRealisticNearDetectorVolume ( bool  t)
inline

Definition at line 100 of file LBNERunManager.hh.

Member Data Documentation

G4String LBNERunManager::fAsciiOutputFileName
protected

Definition at line 47 of file LBNERunManager.hh.

bool LBNERunManager::fCreateAlcoveTrackingOutput
protected

Definition at line 63 of file LBNERunManager.hh.

bool LBNERunManager::fCreateAsciiOutput
protected

Definition at line 55 of file LBNERunManager.hh.

bool LBNERunManager::fCreateDk2NuOutput
protected

Definition at line 54 of file LBNERunManager.hh.

bool LBNERunManager::fCreateOutput
protected

Definition at line 53 of file LBNERunManager.hh.

bool LBNERunManager::fCreateTargetOutput
protected

Definition at line 62 of file LBNERunManager.hh.

bool LBNERunManager::fCreateTrackingOutput
protected

Definition at line 60 of file LBNERunManager.hh.

bool LBNERunManager::fCreateTrackingTargetOutASCII
protected

Definition at line 61 of file LBNERunManager.hh.

bool LBNERunManager::fCreateTrkPlaneDPOutput
protected

Definition at line 59 of file LBNERunManager.hh.

bool LBNERunManager::fCreateTrkPlaneH1Output
protected

Definition at line 57 of file LBNERunManager.hh.

bool LBNERunManager::fCreateTrkPlaneH2Output
protected

Definition at line 58 of file LBNERunManager.hh.

bool LBNERunManager::fCreateTrkPlaneOutput
protected

Definition at line 56 of file LBNERunManager.hh.

G4String LBNERunManager::fDetectorLocationFileName
protected

Definition at line 50 of file LBNERunManager.hh.

G4String LBNERunManager::fDk2NuOutputFileName
protected

Definition at line 49 of file LBNERunManager.hh.

bool LBNERunManager::fGeometryIntializedHere
protected

Definition at line 34 of file LBNERunManager.hh.

LBNEEventAction* LBNERunManager::fLBNEEventAction
protected

Definition at line 38 of file LBNERunManager.hh.

LBNERunAction* LBNERunManager::fLBNERunAction
protected

Definition at line 42 of file LBNERunManager.hh.

LBNEStackingAction* LBNERunManager::fLBNEStackingAction
protected

Definition at line 40 of file LBNERunManager.hh.

LBNESteppingAction* LBNERunManager::fLBNESteppingAction
protected

Definition at line 39 of file LBNERunManager.hh.

LBNETrackingAction* LBNERunManager::fLBNETrackingAction
protected

Definition at line 41 of file LBNERunManager.hh.

G4String LBNERunManager::fMarsOrFlukaInputFileName
protected

Definition at line 46 of file LBNERunManager.hh.

G4String LBNERunManager::fNptOutputFileName
protected

Definition at line 48 of file LBNERunManager.hh.

bool LBNERunManager::fPhysicsInitializedHere
protected

Definition at line 35 of file LBNERunManager.hh.

G4String LBNERunManager::fPhysicsListName
protected

Definition at line 36 of file LBNERunManager.hh.

bool LBNERunManager::fUseFluka
protected

Definition at line 51 of file LBNERunManager.hh.

bool LBNERunManager::fUseMars
protected

Definition at line 52 of file LBNERunManager.hh.

bool LBNERunManager::fUseRealisticNearDetectorVolume
protected

Definition at line 64 of file LBNERunManager.hh.

int LBNERunManager::fVerboseLevel
protected

Definition at line 65 of file LBNERunManager.hh.

G4int LBNERunManager::nEvents

Definition at line 30 of file LBNERunManager.hh.


The documentation for this class was generated from the following files: