Class Reference for E1039 Core & Analysis Software
PHG4Reco.cc
Go to the documentation of this file.
1 #include "PHG4Reco.h"
2 
4 #include "PHG4InEvent.h"
5 #include "PHG4PhenixDetector.h"
10 #include "PHG4Subsystem.h"
11 #include "PHG4TrackingAction.h"
12 #include "PHG4UIsession.h"
13 #include "PHG4Utils.h"
14 
15 #include <phgeom/PHGeomUtility.h>
16 
18 
19 #include <phfield/PHFieldUtility.h>
20 #include <phfield/PHFieldConfig_v1.h>
21 #include <phfield/PHFieldConfig_v2.h>
22 #include <phfield/PHFieldConfig_v3.h>
23 
25 
26 #include <phool/PHCompositeNode.h>
27 #include <phool/PHRandomSeed.h>
28 #include <phool/getClass.h>
29 #include <phool/recoConsts.h>
30 
31 #include <phparameter/PHParameters.h>
32 
33 #include <TThread.h>
34 #include <TSystem.h>
35 
36 #include <CLHEP/Random/Random.h>
37 
38 #include <Geant4/G4RunManager.hh>
39 
40 #include <Geant4/G4Material.hh>
41 #include <Geant4/G4NistManager.hh>
42 #include <Geant4/G4OpenGLImmediateX.hh>
43 #include <Geant4/G4Region.hh>
44 #include <Geant4/G4RegionStore.hh>
45 #include <Geant4/G4StepLimiterPhysics.hh>
46 #include <Geant4/G4UIExecutive.hh>
47 #include <Geant4/G4UImanager.hh>
48 #include <Geant4/G4VisExecutive.hh>
49 #include <Geant4/G4Cerenkov.hh>
50 #include <Geant4/G4EmProcessOptions.hh>
51 #include <Geant4/G4EmSaturation.hh>
52 #include <Geant4/G4FastSimulationManager.hh>
53 #include <Geant4/G4HadronicProcessStore.hh>
54 #include <Geant4/G4LossTableManager.hh>
55 #include <Geant4/G4OpAbsorption.hh>
56 #include <Geant4/G4OpBoundaryProcess.hh>
57 #include <Geant4/G4OpMieHG.hh>
58 #include <Geant4/G4OpRayleigh.hh>
59 #include <Geant4/G4OpWLS.hh>
60 #include <Geant4/G4OpticalPhoton.hh>
61 #include <Geant4/G4OpticalPhysics.hh>
62 #include <Geant4/G4PAIModel.hh>
63 #include <Geant4/G4PEEffectFluoModel.hh>
64 #include <Geant4/G4EmProcessOptions.hh>
65 #include <Geant4/G4HadronicProcessStore.hh>
66 #include <Geant4/G4StepLimiterPhysics.hh>
67 #include <Geant4/G4ParticleDefinition.hh>
68 #include <Geant4/G4ParticleTable.hh>
69 #include <Geant4/G4ParticleTypes.hh>
70 #include <Geant4/G4ProcessManager.hh>
71 #include <Geant4/G4Scintillation.hh>
72 #include <Geant4/G4SystemOfUnits.hh>
73 #include <Geant4/G4Version.hh>
74 #include <Geant4/globals.hh>
75 
76 // physics lists
77 #include <Geant4/FTFP_BERT.hh>
78 #include <Geant4/LBE.hh>
79 #include <Geant4/QGSP_BERT.hh>
80 #include <Geant4/QGSP_BIC.hh>
81 #include <Geant4/QGSP_BIC_HP.hh>
82 
83 #if G4VERSION_NUMBER <= 951
84 #define HAVE_LHEP
85 #include <Geant4/LHEP.hh>
86 #endif
87 
88 #if G4VERSION_NUMBER >= 1001
89 #define HAVE_FTFP_BERT_HP
90 #define HAVE_QGSP_BERT_HP
91 #include <Geant4/FTFP_BERT_HP.hh>
92 #include <Geant4/QGSP_BERT_HP.hh>
93 #endif
94 
95 //#include <boost/filesystem.hpp>
96 //#include <boost/foreach.hpp>
97 
98 #include <memory>
99 #include <cassert>
100 #include <cstdlib>
101 #include <cstring>
102 #include <sstream>
103 #include <iterator>
104 
105 using namespace std;
106 
107 static TThread *gui_thread = nullptr;
108 
109 // for the G4 cmd line interface
110 G4UImanager *UImanager = nullptr;
111 
112 // the gui thread
113 void g4guithread(void *ptr);
114 
115 //_________________________________________________________________
116 PHG4Reco::PHG4Reco(const string &name)
117  : SubsysReco(name)
118  , magfield(0.)
119  , magfield_rescale(1.0)
120  , field_(nullptr)
121  , runManager_(nullptr)
122  , uisession_(nullptr)
123  , detector_(nullptr)
124  , eventAction_(nullptr)
125  , steppingAction_(nullptr)
126  , trackingAction_(nullptr)
127  , generatorAction_(nullptr)
128  , visManager(nullptr)
129  , _eta_coverage(1.0)
130  , mapdim(PHFieldConfig::kFieldUniform)
131  , fieldmapfile("")
132  , worldshape("G4Tubs")
133  , worldmaterial("G4_AIR")
134 #if G4VERSION_NUMBER >= 1033
135  , physicslist("FTFP_BERT")
136 #else
137  , physicslist("QGSP_BERT")
138 #endif
139  , save_DST_geometry_(true)
140  , optical_photon_(false)
141  , energy_threshold_(0.0005)
142  , zero_field_line_(99999.)
143  , _timer(PHTimeServer::get()->insert_new(name))
144 {
145  for (int i = 0; i < 3; i++)
146  {
147  WorldSize[i] = 1000.;
148  }
149  return;
150 }
151 
152 //_________________________________________________________________
154 {
155  // one can delete null pointer (it results in a nop), so checking if
156  // they are non zero is not needed
157  delete gui_thread;
158  delete field_;
159  delete runManager_;
160  delete uisession_;
161  delete visManager;
162  while (subsystems_.begin() != subsystems_.end())
163  {
164  delete subsystems_.back();
165  subsystems_.pop_back();
166  }
167 }
168 
169 //_________________________________________________________________
171 {
172  if (verbosity > 0)
173  {
174  cout << "========================= PHG4Reco::Init() ================================" << endl;
175  }
176  unsigned int iseed = PHRandomSeed();
177  G4Seed(iseed); // fixed seed handled in PHRandomSeed()
178 
179  // create GEANT run manager
180  if (verbosity > 1) cout << "PHG4Reco::Init - create run manager" << endl;
181 
182  // redirect GEANT verbosity to nowhere
183 // if (verbosity < 1)
184  if (0)
185  {
186  G4UImanager *uimanager = G4UImanager::GetUIpointer();
187  uisession_ = new PHG4UIsession();
188  uimanager->SetSession(uisession_);
189  uimanager->SetCoutDestination(uisession_);
190  }
191 
192  runManager_ = new G4RunManager();
193 
194  DefineMaterials();
195  // create physics processes
196  G4VModularPhysicsList *myphysicslist = nullptr;
197  if (physicslist == "FTFP_BERT")
198  {
199  myphysicslist = new FTFP_BERT(verbosity);
200  }
201  else if (physicslist == "QGSP_BERT")
202  {
203  myphysicslist = new QGSP_BERT(verbosity);
204  }
205  else if (physicslist == "QGSP_BIC")
206  {
207  myphysicslist = new QGSP_BIC(verbosity);
208  }
209  else if (physicslist == "QGSP_BIC_HP")
210  {
211  setenv("AllowForHeavyElements", "1", 1);
212  myphysicslist = new QGSP_BIC_HP(verbosity);
213  }
214 #ifdef HAVE_LHEP
215  else if (physicslist == "LHEP")
216  {
217  myphysicslist = new LHEP(verbosity);
218  }
219 #endif
220 #ifdef HAVE_FTFP_BERT_HP
221  else if (physicslist == "FTFP_BERT_HP")
222  {
223  setenv("AllowForHeavyElements", "1", 1);
224  myphysicslist = new FTFP_BERT_HP(verbosity);
225  }
226 #endif
227 #ifdef HAVE_QGSP_BERT_HP
228  else if (physicslist == "QGSP_BERT_HP")
229  {
230  setenv("AllowForHeavyElements", "1", 1);
231  myphysicslist = new QGSP_BERT_HP(verbosity);
232  }
233 #endif
234  else
235  {
236  cout << "Physics List " << physicslist << " not implemented" << endl;
237  gSystem->Exit(1);
238  }
239 
240  myphysicslist->RegisterPhysics(new G4StepLimiterPhysics());
241 // initialize cuts so we can ask the world region for it's default
242 // cuts to propagate them to other regions in DefineRegions()
243  myphysicslist->SetCutsWithDefault();
244  runManager_->SetUserInitialization(myphysicslist);
245 
246  DefineRegions();
247 #if G4VERSION_NUMBER >= 1033
248  G4EmSaturation* emSaturation = G4LossTableManager::Instance()->EmSaturation();
249  if (! emSaturation)
250  {
251  cout << PHWHERE << "Could not initialize EmSaturation, Birks constants will fail" << endl;
252  }
253 #endif
254  // initialize registered subsystems
255  //BOOST_FOREACH (SubsysReco *reco, subsystems_)
256  for (SubsysReco *reco : subsystems_)
257  {
258  reco->Init(topNode);
259  }
260 
261  if (verbosity > 0)
262  {
263  cout << "===========================================================================" << endl;
264  }
265  return 0;
266 }
267 
269 {
270  if (verbosity > 1) cout << "PHG4Reco::InitField - create magnetic field setup" << endl;
271 
272  unique_ptr<PHFieldConfig> default_field_cfg(nullptr);
273 
274  std::istringstream ssfieldmapfile(fieldmapfile);
275  std::vector<std::string> fieldmapfiles(std::istream_iterator<std::string>{ssfieldmapfile},
276  std::istream_iterator<std::string>());
277 
278  if (fieldmapfiles.size() == 1)
279  {
280  default_field_cfg.reset(new PHFieldConfig_v1(mapdim, fieldmapfile, magfield_rescale));
281  }
282  else if (fieldmapfiles.size() == 0)
283  {
284  default_field_cfg.reset(new PHFieldConfig_v2(0, magfield * magfield_rescale, 0));
285  }
286  else if (fieldmapfiles.size() == 2)
287  {
288  default_field_cfg.reset(new PHFieldConfig_v3(
289  fieldmapfiles[0], fieldmapfiles[1], 1.0, 1.0, 5.0
290  ));
291  }
292  else if (fieldmapfiles.size() == 5)
293  {
294  default_field_cfg.reset(new PHFieldConfig_v3(
295  fieldmapfiles[0], fieldmapfiles[1], std::stod(fieldmapfiles[2]), std::stod(fieldmapfiles[3]), std::stod(fieldmapfiles[4])
296  ));
297  }
298 
299  if (verbosity > 1) cout << "PHG4Reco::InitField - create magnetic field setup" << endl;
300 
301  PHField * phfield = PHFieldUtility::GetFieldMapNode(default_field_cfg.get(), topNode, Verbosity()+1);
302  assert(phfield);
303 
304  field_ = new G4TBMagneticFieldSetup(phfield);
305 
307 }
308 
310 {
312  ostringstream oss;
313  oss << rc->get_CharFlag("fMagFile") << " "
314  << rc->get_CharFlag("kMagFile") << " "
315  << rc->get_DoubleFlag("FMAGSTR") << " "
316  << rc->get_DoubleFlag("KMAGSTR") << " "
317  << rc->get_DoubleFlag("TMAGSTR");
318  if (verbosity > 0) cout << "PHG4Reco::set_field_map(): " << oss.str() << "." << endl;
320 }
321 
323 {
324  // this is a dumb protection against executing this twice.
325  // we have cases (currently detector display or material scan) where
326  // we need the detector bu have not run any events (who wants to wait
327  // for processing an event if you just want a detector display?).
328  // Then the InitRun is executed from a macro. But if you decide to run events
329  // afterwards, the InitRun is executed by the framework together with all
330  // other modules. This prevents the PHG4Reco::InitRun() to be executed
331  // again in this case
332  static int InitRunDone = 0;
333  if (InitRunDone)
334  {
336  }
337  InitRunDone = 1;
338  if (verbosity > 0)
339  {
340  cout << "========================= PHG4Reco::InitRun() ================================" << endl;
341  }
342 
344 
345  //setup the constant field
346  const int field_ret = InitField(topNode);
347  if (field_ret!=Fun4AllReturnCodes::EVENT_OK)
348  {
349  cout <<"PHG4Reco::InitRun- Error - Failed field init with status = "<<field_ret<<endl;
350  return field_ret;
351  }
352 
353  // initialize registered subsystems
354  //BOOST_FOREACH (SubsysReco *reco, subsystems_)
355  for (SubsysReco *reco : subsystems_)
356  {
357  reco->InitRun(topNode);
358  }
359 
360  // create phenix detector, add subsystems, and register to GEANT
361  if (verbosity > 1) cout << "PHG4Reco::Init - create detector" << endl;
371 
372  rc->set_FloatFlag("WorldSizex", WorldSize[0]);
373  rc->set_FloatFlag("WorldSizey", WorldSize[1]);
374  rc->set_FloatFlag("WorldSizez", WorldSize[2]);
375 
376  double z_offset = 0.;
377  if(fabs(zero_field_line_) < 0.5*WorldSize[2])
378  {
379  z_offset = 0.5*(zero_field_line_ + 0.5*WorldSize[2]); //unit cm
380  cout << "Will update the z-positions of the downstream volumes by zOffset = " << z_offset << "cm." << endl;
381  }
382 
383  //BOOST_FOREACH (PHG4Subsystem *g4sub, subsystems_)
384  for (PHG4Subsystem *g4sub : subsystems_)
385  {
386  //re-calculate the place_z for some of the volumes in the non-field-zone
387  if(fabs(zero_field_line_) < 0.5*WorldSize[2] && g4sub->GetParams() != nullptr) //only true for class inherited from PHG4DetectorSubsystem
388  {
389  double z_in_world = g4sub->GetParams()->get_double_param("place_z");
390  if(z_in_world > zero_field_line_)
391  {
392  g4sub->GetParams()->set_double_param("place_z", z_in_world - z_offset);
393  detector_->AddDetector(g4sub->GetDetector(), 1);
394  }
395  else
396  {
397  detector_->AddDetector(g4sub->GetDetector(), 0);
398  }
399  }
400  else
401  {
402  detector_->AddDetector(g4sub->GetDetector(), 0);
403  }
404  }
405  runManager_->SetUserInitialization(detector_);
406 
407  setupInputEventNodeReader(topNode);
408  // create main event action, add subsystemts and register to GEANT
410 
411  //BOOST_FOREACH (PHG4Subsystem *g4sub, subsystems_)
412  for (PHG4Subsystem *g4sub : subsystems_)
413  {
414  PHG4EventAction *evtact = g4sub->GetEventAction();
415  if (evtact)
416  {
417  eventAction_->AddAction(evtact);
418  }
419  }
420  runManager_->SetUserAction(eventAction_);
421 
422  // create main stepping action, add subsystems and register to GEANT
425  //BOOST_FOREACH (PHG4Subsystem *g4sub, subsystems_)
426  for (PHG4Subsystem *g4sub : subsystems_)
427  {
428  PHG4SteppingAction *action = g4sub->GetSteppingAction();
429  if (action)
430  {
431  if (verbosity > 1)
432  {
433  cout << "Adding steppingaction for " << g4sub->Name() << endl;
434  }
435  steppingAction_->AddAction(g4sub->GetSteppingAction());
436  }
437  }
438  runManager_->SetUserAction(steppingAction_);
439 
440  // create main tracking action, add subsystems and register to GEANT
442  //BOOST_FOREACH (PHG4Subsystem *g4sub, subsystems_)
443  for (PHG4Subsystem *g4sub : subsystems_)
444  {
445  trackingAction_->AddAction(g4sub->GetTrackingAction());
446 
447  // not all subsystems define a user tracking action
448  if (g4sub->GetTrackingAction())
449  {
450  // make tracking manager accessible within user tracking action if defined
451  if (G4TrackingManager *trackingManager = G4EventManager::GetEventManager()->GetTrackingManager())
452  {
453  g4sub->GetTrackingAction()->SetTrackingManagerPointer(trackingManager);
454  }
455  }
456  }
457 
458  runManager_->SetUserAction(trackingAction_);
459 
460  // initialize
461  runManager_->Initialize();
462 
463  if(optical_photon_)
464  {
465  // add cerenkov and optical photon processes
466  // cout << endl << "Ignore the next message - we implemented this correctly" << endl;
467  G4Cerenkov *theCerenkovProcess = new G4Cerenkov("Cerenkov");
468  // cout << "End of bogus warning message" << endl << endl;
469  // G4Scintillation* theScintillationProcess = new G4Scintillation("Scintillation");
470 
471 
472  // if (verbosity > 0)
473  // {
474  // // This segfaults
475  // theCerenkovProcess->DumpPhysicsTable();
476  // }
477 
478  theCerenkovProcess->SetMaxNumPhotonsPerStep(100);
479  theCerenkovProcess->SetMaxBetaChangePerStep(10.0);
480  theCerenkovProcess->SetTrackSecondariesFirst(true);
481 
482  // theScintillationProcess->SetScintillationYieldFactor(1.);
483  // theScintillationProcess->SetTrackSecondariesFirst(true);
484 
485  // Use Birks Correction in the Scintillation process
486 
487  // G4EmSaturation* emSaturation = G4LossTableManager::Instance()->EmSaturation();
488  // theScintillationProcess->AddSaturation(emSaturation);
489 
490  G4ParticleTable *theParticleTable = G4ParticleTable::GetParticleTable();
491  G4ParticleTable::G4PTblDicIterator *_theParticleIterator;
492  _theParticleIterator = theParticleTable->GetIterator();
493  _theParticleIterator->reset();
494  while ((*_theParticleIterator)())
495  {
496  G4ParticleDefinition *particle = _theParticleIterator->value();
497  G4String particleName = particle->GetParticleName();
498  G4ProcessManager *pmanager = particle->GetProcessManager();
499  if (theCerenkovProcess->IsApplicable(*particle))
500  {
501  pmanager->AddProcess(theCerenkovProcess);
502  pmanager->SetProcessOrdering(theCerenkovProcess, idxPostStep);
503  }
504  // if (theScintillationProcess->IsApplicable(*particle))
505  // {
506  // pmanager->AddProcess(theScintillationProcess);
507  // pmanager->SetProcessOrderingToLast(theScintillationProcess, idxAtRest);
508  // pmanager->SetProcessOrderingToLast(theScintillationProcess, idxPostStep);
509  // }
510  }
511  G4ProcessManager *pmanager = G4OpticalPhoton::OpticalPhoton()->GetProcessManager();
512  // G4cout << " AddDiscreteProcess to OpticalPhoton " << G4endl;
513  pmanager->AddDiscreteProcess(new G4OpAbsorption());
514  pmanager->AddDiscreteProcess(new G4OpRayleigh());
515  pmanager->AddDiscreteProcess(new G4OpMieHG());
516  pmanager->AddDiscreteProcess(new G4OpBoundaryProcess());
517  pmanager->AddDiscreteProcess(new G4OpWLS());
518  pmanager->AddDiscreteProcess(new G4PhotoElectricEffect());
519  // pmanager->DumpInfo();
520 
521  // needs large amount of memory which kills central hijing events
522  // store generated trajectories
523  //if( G4TrackingManager* trackingManager = G4EventManager::GetEventManager()->GetTrackingManager() ){
524  // trackingManager->SetStoreTrajectory( true );
525  //}
526  }
527 
528  // quiet some G4 print-outs (EM and Hadronic settings during first event)
529  G4HadronicProcessStore::Instance()->SetVerbose(0);
530  G4LossTableManager::Instance()->SetVerbose(1);
531 
532  if ((verbosity < 1) && (uisession_))
533  {
534  uisession_->Verbosity(1); // let messages after setup come through
535  }
536 
537  // Geometry export to DST
538  if (save_DST_geometry_)
539  {
540  const string filename =
543  cout << "PHG4Reco::InitRun - export geometry to DST via tmp file " << filename << endl;
544 
545  Dump_GDML(filename);
546 
547  PHGeomUtility::ImportGeomFile(topNode, filename);
548 
550  }
551 
552  if (verbosity > 0)
553  {
554  cout << "===========================================================================" << endl;
555  }
556 
557  return 0;
558 }
559 
560 //________________________________________________________________
561 //Dump TGeo File
562 void PHG4Reco::Dump_GDML(const std::string &filename)
563 {
565 }
566 
567 //_________________________________________________________________
568 int PHG4Reco::ApplyCommand(const std::string &cmd)
569 {
570  InitUImanager();
571  int iret = UImanager->ApplyCommand(cmd.c_str());
572  return iret;
573 }
574 //_________________________________________________________________
575 
577 {
578  if (!gui_thread)
579  {
580  InitUImanager();
581  gui_thread = new TThread("G4Gui", g4guithread);
582  gui_thread->Run();
583  return 0;
584  }
585  return 1;
586 }
587 
589 {
590  if (!UImanager)
591  {
592  // Get the pointer to the User Interface manager
593  // Initialize visualization
594  visManager = new G4VisExecutive;
595  visManager->Initialize();
596  UImanager = G4UImanager::GetUIpointer();
597  }
598  return 0;
599 }
600 
601 //_________________________________________________________________
603 {
604  TThread::Lock();
605  // make sure Actions and subsystems have the relevant pointers set
606  PHG4InEvent *ineve = findNode::getClass<PHG4InEvent>(topNode, "PHG4INEVENT");
608 
609  //BOOST_FOREACH (SubsysReco *reco, subsystems_)
610  for (SubsysReco *reco : subsystems_)
611  {
612  if (Verbosity() >= 2)
613  cout << "PHG4Reco::process_event - " << reco->Name() << "->process_event" << endl;
614 
615  try
616  {
617  reco->process_event(topNode);
618  }
619  catch (const exception &e)
620  {
621  cout << PHWHERE << " caught exception thrown during process_event from "
622  << reco->Name() << endl;
623  cout << "error: " << e.what() << endl;
624  TThread::UnLock();
626  }
627  }
628 
629  _timer.get()->restart();
630 
631  // run one event
632  if (Verbosity() >= 2)
633  {
634  cout << " PHG4Reco::process_event - "
635  << "run one event :" << endl;
636  ineve->identify();
637  }
638  runManager_->BeamOn(1);
639  _timer.get()->stop();
640 
641  //BOOST_FOREACH (PHG4Subsystem *g4sub, subsystems_)
642  for (PHG4Subsystem *g4sub : subsystems_)
643  {
644  if (Verbosity() >= 2)
645  cout << " PHG4Reco::process_event - " << g4sub->Name() << "->process_after_geant" << endl;
646  try
647  {
648  g4sub->process_after_geant(topNode);
649  }
650  catch (const exception &e)
651  {
652  cout << PHWHERE << " caught exception thrown during process_after_geant from "
653  << g4sub->Name() << endl;
654  cout << "error: " << e.what() << endl;
655  TThread::UnLock();
657  }
658  }
659  TThread::UnLock();
660  return 0;
661 }
662 
664 {
665  //BOOST_FOREACH (SubsysReco *reco, subsystems_)
666  for (SubsysReco *reco : subsystems_)
667  {
668  reco->ResetEvent(topNode);
669  }
670  return 0;
671 }
672 
673 //_________________________________________________________________
675 {
676  return 0;
677 }
678 
679 void PHG4Reco::Print(const std::string &what) const
680 {
681  //BOOST_FOREACH (SubsysReco *reco, subsystems_)
682  for (SubsysReco *reco : subsystems_)
683  {
684  if (what.empty() || what == "ALL" || (reco->Name()).find(what) != string::npos)
685  {
686  cout << "Printing " << reco->Name() << endl;
687  reco->Print(what);
688  cout << "---------------------------" << endl;
689  }
690  }
691  return;
692 }
693 
695 {
696  PHG4InEvent *ineve = findNode::getClass<PHG4InEvent>(topNode, "PHG4INEVENT");
697  if (!ineve)
698  {
699  PHNodeIterator iter(topNode);
700  PHCompositeNode *dstNode;
701  dstNode = dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "DST"));
702 
703  ineve = new PHG4InEvent();
704  PHIODataNode<PHObject> *newNode = new PHIODataNode<PHObject>(ineve, "PHG4INEVENT", "PHObject");
705  dstNode->addNode(newNode);
706  }
708  runManager_->SetUserAction(generatorAction_);
709  return 0;
710 }
711 
712 void PHG4Reco::setGeneratorAction(G4VUserPrimaryGeneratorAction *action)
713 {
714  if (runManager_)
715  {
716  runManager_->SetUserAction(action);
717  }
718  return;
719 }
720 
721 void PHG4Reco::set_rapidity_coverage(const double eta)
722 {
723  _eta_coverage = eta;
725 }
726 
727 void PHG4Reco::G4Seed(const unsigned int i)
728 {
729  CLHEP::HepRandom::setTheSeed(i);
730  return;
731 }
732 
733 void g4guithread(void *ptr)
734 {
735  TThread::Lock();
736  G4UIExecutive *ui = new G4UIExecutive(0, nullptr, "Xm");
737  //FIXME test file exist
738  //if (ui->IsGUI() && boost::filesystem::exists("gui.mac"))
739  if (ui->IsGUI() && true)
740  {
741  UImanager->ApplyCommand("/control/execute gui.mac");
742  }
743  TThread::UnLock();
744  ui->SessionStart();
745  TThread::Lock();
746  delete ui;
747  TThread::UnLock();
748  gui_thread = nullptr;
749  return;
750 }
751 
752 //____________________________________________________________________________
754 {
755  G4String symbol; //a=mass of a mole;
756  G4double density; //z=mean number of protons;
757  G4double fractionmass;
758 
759  G4int ncomponents, natoms;
760  // this is for FTFP_BERT_HP where the neutron code barfs
761  // if the z difference to the last known element (U) is too large
762  set<G4String> ignoremat;
763  ignoremat.insert("G4_Cf");
764 
765  //load all Materials from the nist database
766  G4NistManager *nist = G4NistManager::Instance();
767  vector<G4String> matnames = nist->GetNistMaterialNames();
768  while (matnames.begin() != matnames.end())
769  {
770  G4String mat = matnames.back();
771  if (ignoremat.find(mat) == ignoremat.end())
772  {
773  nist->FindOrBuildMaterial(mat);
774  }
775  matnames.pop_back();
776  }
777  // home made compounds
778  // making quartz
779  G4Material *quartz = new G4Material("Quartz", density = 2.200 * g / cm3, ncomponents = 2);
780  quartz->AddElement(G4Element::GetElement("Si"), 1);
781  quartz->AddElement(G4Element::GetElement("O"), 2);
782 
783  // gas mixture for the MuID in fsPHENIX. CLS 02-25-14
784  G4Material *IsoButane = new G4Material("Isobutane", 0.00265 * g / cm3, 2);
785  IsoButane->AddElement(G4Element::GetElement("C"), 4);
786  IsoButane->AddElement(G4Element::GetElement("H"), 10);
787 
788  G4Material *MuIDgas = new G4Material("MuIDgas", density = (1.977e-3 * 0.92 + 0.00265 * 0.08) * g / cm3, ncomponents = 2);
789  MuIDgas->AddMaterial(IsoButane, fractionmass = 0.08);
790  MuIDgas->AddMaterial(G4Material::GetMaterial("G4_CARBON_DIOXIDE"), fractionmass = 0.92);
791 
792  // that seems to be the composition of 304 Stainless steel
793  G4Material *StainlessSteel =
794  new G4Material("SS304", density = 7.9 * g / cm3, ncomponents = 8);
795  StainlessSteel->AddElement(G4Element::GetElement("Fe"), 0.70105);
796  StainlessSteel->AddElement(G4Element::GetElement("Cr"), 0.18);
797  StainlessSteel->AddElement(G4Element::GetElement("Ni"), 0.09);
798  StainlessSteel->AddElement(G4Element::GetElement("Mn"), 0.02);
799  StainlessSteel->AddElement(G4Element::GetElement("C"), 0.0007);
800  StainlessSteel->AddElement(G4Element::GetElement("S"), 0.0003);
801  StainlessSteel->AddElement(G4Element::GetElement("Si"), 0.0075);
802  StainlessSteel->AddElement(G4Element::GetElement("P"), 0.00045);
803 
804  G4Material *SS310 =
805  new G4Material("SS310", density = 8.0 * g / cm3, ncomponents = 8);
806  SS310->AddElement(G4Element::GetElement("Fe"), 0.50455);
807  SS310->AddElement(G4Element::GetElement("Cr"), 0.25);
808  SS310->AddElement(G4Element::GetElement("Ni"), 0.20);
809  SS310->AddElement(G4Element::GetElement("Mn"), 0.02);
810  SS310->AddElement(G4Element::GetElement("C"), 0.0025);
811  SS310->AddElement(G4Element::GetElement("S"), 0.015);
812  SS310->AddElement(G4Element::GetElement("Si"), 0.0075);
813  SS310->AddElement(G4Element::GetElement("P"), 0.00045);
814 
815  G4Material *Steel =
816  new G4Material("Steel", density = 7.86 * g / cm3, ncomponents = 5);
817  Steel->AddElement(G4Element::GetElement("Fe"), 0.9834);
818  Steel->AddElement(G4Element::GetElement("Mn"), 0.014);
819  Steel->AddElement(G4Element::GetElement("C"), 0.0017);
820  Steel->AddElement(G4Element::GetElement("S"), 0.00045);
821  Steel->AddElement(G4Element::GetElement("P"), 0.00045);
822 
823  // a36 steel from http://www.matweb.com
824  G4Material *a36 = new G4Material("Steel_A36", density = 7.85 * g / cm3, ncomponents = 5);
825  a36->AddElement(G4Element::GetElement("Fe"), 0.9824);
826  a36->AddElement(G4Element::GetElement("Cu"), 0.002);
827  a36->AddElement(G4Element::GetElement("C"), 0.0025);
828  a36->AddElement(G4Element::GetElement("Mn"), 0.0103);
829  a36->AddElement(G4Element::GetElement("Si"), 0.0028);
830 
831  // 1006 steel from http://www.matweb.com
832  G4Material *steel_1006 = new G4Material("Steel_1006", density = 7.872 * g / cm3, ncomponents = 2);
833  steel_1006->AddElement(G4Element::GetElement("Fe"), 0.996);
834  steel_1006->AddElement(G4Element::GetElement("Mn"), 0.004);
835 
836  // from www.aalco.co.uk
837  G4Material *Al5083 = new G4Material("Al5083", density = 2.65 * g / cm3, ncomponents = 3);
838  Al5083->AddElement(G4Element::GetElement("Mn"), 0.004);
839  Al5083->AddElement(G4Element::GetElement("Mg"), 0.04);
840  Al5083->AddElement(G4Element::GetElement("Al"), 0.956);
841 
842  G4Material *FPC = new G4Material("FPC", 1.542 * g / cm3, 2);
843  FPC->AddMaterial(G4Material::GetMaterial("G4_Cu"), 0.0162);
844  FPC->AddMaterial(G4Material::GetMaterial("G4_KAPTON"), 0.9838);
845 
846  // This is an approximation for the W saturated epoxy of the EMCal.
847  G4Material *W_Epoxy = new G4Material("W_Epoxy", density = 10.2 * g / cm3, ncomponents = 2);
848  W_Epoxy->AddMaterial(G4Material::GetMaterial("G4_W"), fractionmass = 0.5);
849  W_Epoxy->AddMaterial(G4Material::GetMaterial("G4_POLYSTYRENE"), fractionmass = 0.5);
850 
851  //from http://www.physi.uni-heidelberg.de/~adler/TRD/TRDunterlagen/RadiatonLength/tgc2.htm
852  //Epoxy (for FR4 )
853  //density = 1.2*g/cm3;
854  G4Material *Epoxy = new G4Material("Epoxy", 1.2 * g / cm3, ncomponents = 2);
855  Epoxy->AddElement(G4Element::GetElement("H"), natoms = 2);
856  Epoxy->AddElement(G4Element::GetElement("C"), natoms = 2);
857 
858  //FR4 (Glass + Epoxy)
859  density = 1.86 * g / cm3;
860  G4Material *FR4 = new G4Material("FR4", density, ncomponents = 2);
861  FR4->AddMaterial(quartz, fractionmass = 0.528);
862  FR4->AddMaterial(Epoxy, fractionmass = 0.472);
863 // NOMEX (HoneyComb)
864 // density from http://www.fibreglast.com/product/Nomex_Honeycomb_1562/Vacuum_Bagging_Sandwich_Core
865 // 1562: 29 kg/m^3 <-- I guess it is this one
866 // 2562: 48 kg/m^3
867 // chemical composition http://ww2.unime.it/cdlchimind/adm/inviofile/uploads/HP_Pols2.b.pdf
868  density = 29*kg/m3;
869  G4Material *NOMEX = new G4Material("NOMEX", density, ncomponents = 4);
870  NOMEX->AddElement(G4Element::GetElement("C"), natoms = 14);
871  NOMEX->AddElement(G4Element::GetElement("H"), natoms = 10);
872  NOMEX->AddElement(G4Element::GetElement("N"), natoms = 2);
873  NOMEX->AddElement(G4Element::GetElement("O"), natoms = 2);
874  // spacal material. Source : EICROOT/A. Kiselev
875  /*
876  WEpoxyMix 3 12.011 1.008 183.85 6. 1. 74. 12.18 0.029 0.002 0.969
877  1 1 30. .00001
878  0
879  */
880  G4Material *Spacal_W_Epoxy =
881  new G4Material("Spacal_W_Epoxy", density = 12.18 * g / cm3, ncomponents = 3);
882  Spacal_W_Epoxy->AddElement(G4Element::GetElement("C"), 0.029);
883  Spacal_W_Epoxy->AddElement(G4Element::GetElement("H"), 0.002);
884  Spacal_W_Epoxy->AddElement(G4Element::GetElement("W"), 0.969);
885  /*
886 PMMA -3 12.01 1.008 15.99 6. 1. 8. 1.19 3.6 5.7 1.4
887  1 1 20. .00001
888  0
889  */
890  G4Material *PMMA =
891  new G4Material("PMMA", density = 1.19 * g / cm3, ncomponents = 3);
892  PMMA->AddElement(G4Element::GetElement("C"), 3.6 / (3.6 + 5.7 + 1.4));
893  PMMA->AddElement(G4Element::GetElement("H"), 5.7 / (3.6 + 5.7 + 1.4));
894  PMMA->AddElement(G4Element::GetElement("O"), 1.4 / (3.6 + 5.7 + 1.4));
895 
896  G4Material *G10 =
897  new G4Material("G10", density = 1.700 * g / cm3, ncomponents = 4);
898  G10->AddElement(G4Element::GetElement("Si"), natoms = 1);
899  G10->AddElement(G4Element::GetElement("O"), natoms = 2);
900  G10->AddElement(G4Element::GetElement("C"), natoms = 3);
901  G10->AddElement(G4Element::GetElement("H"), natoms = 3);
902 
903  G4Material *CsI =
904  new G4Material("CsI", density = 4.534 * g / cm3, ncomponents = 2);
905  CsI->AddElement(G4Element::GetElement("Cs"), natoms = 1);
906  CsI->AddElement(G4Element::GetElement("I"), natoms = 1);
907  CsI->GetIonisation()->SetMeanExcitationEnergy(553.1 * eV);
908 
909  G4Material *C4F10 =
910  new G4Material("C4F10", density = 0.00973 * g / cm3, ncomponents = 2);
911  C4F10->AddElement(G4Element::GetElement("C"), natoms = 4);
912  C4F10->AddElement(G4Element::GetElement("F"), natoms = 10);
913 
914  G4Material *CF4 = new G4Material("CF4", density = 3.72 * mg / cm3, ncomponents = 2, kStateGas, 288.15 * kelvin, 1 * atmosphere);
915  CF4->AddElement(G4Element::GetElement("C"), natoms = 1);
916  CF4->AddElement(G4Element::GetElement("F"), natoms = 4);
917 
921 
922  const double den_CF4 = CF4->GetDensity() * .1;
923  const double den_G4_Ar = G4Material::GetMaterial("G4_Ar")->GetDensity() * .8;
924  const double den_G4_CARBON_DIOXIDE = G4Material::GetMaterial("G4_CARBON_DIOXIDE")->GetDensity() * .1;
925  const double den = den_CF4 + den_G4_Ar + den_G4_CARBON_DIOXIDE;
926 
927  G4Material *ePHEINX_TPC_Gas = new G4Material("ePHEINX_TPC_Gas", den,
928  ncomponents = 3, kStateGas);
929  ePHEINX_TPC_Gas->AddMaterial(CF4, den_CF4 / den);
930  ePHEINX_TPC_Gas->AddMaterial(G4Material::GetMaterial("G4_Ar"), den_G4_Ar / den);
931  ePHEINX_TPC_Gas->AddMaterial(G4Material::GetMaterial("G4_CARBON_DIOXIDE"),
932  den_G4_CARBON_DIOXIDE / den);
933  // cross checked with original implementation made up of Ne,C,F
934  // this here is very close but makes more sense since it uses Ne and CF4
935  double G4_Ne_frac = 0.9;
936  double CF4_frac = 0.1;
937  const double den_G4_Ne = G4Material::GetMaterial("G4_Ne")->GetDensity();
938  const double den_CF4_2 = CF4->GetDensity();
939  const double den_sphenix_tpc_gas = den_G4_Ne * G4_Ne_frac + den_CF4_2 * CF4_frac;
940  G4Material *sPHENIX_tpc_gas = new G4Material("sPHENIX_TPC_Gas", den_sphenix_tpc_gas, ncomponents = 2, kStateGas);
941  sPHENIX_tpc_gas->AddMaterial(CF4, den_CF4_2 * CF4_frac / den_sphenix_tpc_gas);
942  sPHENIX_tpc_gas->AddMaterial(G4Material::GetMaterial("G4_Ne"), den_G4_Ne * G4_Ne_frac / den_sphenix_tpc_gas);
943  //
944  // CF4
945  //
946  const G4int nEntries_CF4 = 50;
947 
948  G4double PhotonEnergy_CF4[nEntries_CF4] =
949  {5.5 * eV, 5.6 * eV, 5.7 * eV, 5.8 * eV, 5.9 * eV,
950  6.0 * eV, 6.1 * eV, 6.2 * eV, 6.3 * eV, 6.4 * eV,
951  6.5 * eV, 6.6 * eV, 6.7 * eV, 6.8 * eV, 6.9 * eV,
952  7.0 * eV, 7.1 * eV, 7.2 * eV, 7.3 * eV, 7.4 * eV,
953  7.5 * eV, 7.6 * eV, 7.7 * eV, 7.8 * eV, 7.9 * eV,
954  8.0 * eV, 8.1 * eV, 8.2 * eV, 8.4 * eV, 8.6 * eV,
955  8.8 * eV, 9.0 * eV, 9.2 * eV, 9.4 * eV, 9.6 * eV,
956  9.8 * eV, 10.0 * eV, 10.2 * eV, 10.4 * eV, 10.6 * eV,
957  10.8 * eV, 11.0 * eV, 11.2 * eV, 11.3 * eV, 11.4 * eV,
958  11.5 * eV, 11.6 * eV, 11.7 * eV, 11.8 * eV, 11.9 * eV};
959 
960  G4double RefractiveIndex_CF4[nEntries_CF4] =
961  {1.000480, 1.000482, 1.000483, 1.000485, 1.000486,
962  1.000488, 1.000490, 1.000491, 1.000493, 1.000495,
963  1.000497, 1.000498, 1.000500, 1.000502, 1.000504,
964  1.000506, 1.000508, 1.000510, 1.000512, 1.000514,
965  1.000517, 1.000519, 1.000521, 1.000524, 1.000526,
966  1.000529, 1.000531, 1.000534, 1.000539, 1.000545,
967  1.000550, 1.000557, 1.000563, 1.000570, 1.000577,
968  1.000584, 1.000592, 1.000600, 1.000608, 1.000617,
969  1.000626, 1.000636, 1.000646, 1.000652, 1.000657,
970  1.000662, 1.000667, 1.000672, 1.000677, 1.000682};
971 
972  G4double Absorption_CF4[nEntries_CF4] =
973  {1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
974  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
975  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
976  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
977  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
978  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
979  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
980  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
981  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
982  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m};
983 
984  G4MaterialPropertiesTable *MPT_CF4 = new G4MaterialPropertiesTable();
985 
986  MPT_CF4->AddProperty("RINDEX", PhotonEnergy_CF4, RefractiveIndex_CF4, nEntries_CF4)
987  ->SetSpline(true);
988  MPT_CF4->AddProperty("ABSLENGTH", PhotonEnergy_CF4, Absorption_CF4, nEntries_CF4)
989  ->SetSpline(true);
990 
991  CF4->SetMaterialPropertiesTable(MPT_CF4);
992 
993  //
994  // LiF
995  //
996  G4Material *g4_lif = nist->FindOrBuildMaterial("G4_LITHIUM_FLUORIDE");
997  G4Material *LiF = new G4Material("LiF", density = 2.635 * g / cm3, ncomponents = 2);
998  LiF->AddElement(G4Element::GetElement("Li"), natoms = 1);
999  LiF->AddElement(G4Element::GetElement("F"), natoms = 1);
1000 
1001  if (verbosity > 1)
1002  {
1003  cout << "g4_lif: " << g4_lif << endl;
1004  cout << "LiF: " << LiF << endl;
1005  }
1006 
1007  const G4int nEntries_LiF = 50;
1008 
1009  G4double PhotonEnergy_LiF[nEntries_LiF] =
1010  {5.5 * eV, 5.6 * eV, 5.7 * eV, 5.8 * eV, 5.9 * eV,
1011  6.0 * eV, 6.1 * eV, 6.2 * eV, 6.3 * eV, 6.4 * eV,
1012  6.5 * eV, 6.6 * eV, 6.7 * eV, 6.8 * eV, 6.9 * eV,
1013  7.0 * eV, 7.1 * eV, 7.2 * eV, 7.3 * eV, 7.4 * eV,
1014  7.5 * eV, 7.6 * eV, 7.7 * eV, 7.8 * eV, 7.9 * eV,
1015  8.0 * eV, 8.1 * eV, 8.2 * eV, 8.4 * eV, 8.6 * eV,
1016  8.8 * eV, 9.0 * eV, 9.2 * eV, 9.4 * eV, 9.6 * eV,
1017  9.8 * eV, 10.0 * eV, 10.2 * eV, 10.4 * eV, 10.6 * eV,
1018  10.8 * eV, 11.0 * eV, 11.2 * eV, 11.3 * eV, 11.4 * eV,
1019  11.5 * eV, 11.6 * eV, 11.7 * eV, 11.8 * eV, 11.9 * eV};
1020 
1021  G4double RefractiveIndex_LiF[nEntries_LiF] =
1022  {1.42709, 1.42870, 1.42998, 1.43177, 1.43368,
1023  1.43520, 1.43736, 1.43907, 1.44088, 1.44279,
1024  1.44481, 1.44694, 1.44920, 1.45161, 1.45329,
1025  1.45595, 1.45781, 1.46077, 1.46285, 1.46503,
1026  1.46849, 1.47093, 1.47349, 1.47618, 1.47901,
1027  1.48198, 1.48511, 1.48841, 1.49372, 1.50152,
1028  1.50799, 1.51509, 1.52290, 1.53152, 1.54108,
1029  1.54805, 1.55954, 1.56799, 1.58202, 1.59243,
1030  1.60382, 1.61632, 1.63010, 1.63753, 1.64536,
1031  1.65363, 1.66236, 1.67159, 1.68139, 1.69178};
1032 
1033  G4double Absorption_LiF[nEntries_LiF] =
1034  {1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
1035  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
1036  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
1037  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
1038  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
1039  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
1040  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
1041  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
1042  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m,
1043  1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m, 1.e4 * m};
1044 
1045  G4MaterialPropertiesTable *MPT_LiF = new G4MaterialPropertiesTable();
1046 
1047  MPT_LiF->AddProperty("RINDEX", PhotonEnergy_LiF, RefractiveIndex_LiF, nEntries_LiF)
1048  ->SetSpline(true);
1049  MPT_LiF->AddProperty("ABSLENGTH", PhotonEnergy_LiF, Absorption_LiF, nEntries_LiF)
1050  ->SetSpline(true);
1051 
1052  LiF->SetMaterialPropertiesTable(MPT_LiF);
1053 
1054  //
1055  // CsI
1056  //
1057  const G4int nEntries_CsI = 50;
1058 
1059  G4double PhotonEnergy_CsI[nEntries_CsI] =
1060  {5.5 * eV, 5.6 * eV, 5.7 * eV, 5.8 * eV, 5.9 * eV,
1061  6.0 * eV, 6.1 * eV, 6.2 * eV, 6.3 * eV, 6.4 * eV,
1062  6.5 * eV, 6.6 * eV, 6.7 * eV, 6.8 * eV, 6.9 * eV,
1063  7.0 * eV, 7.1 * eV, 7.2 * eV, 7.3 * eV, 7.4 * eV,
1064  7.5 * eV, 7.6 * eV, 7.7 * eV, 7.8 * eV, 7.9 * eV,
1065  8.0 * eV, 8.1 * eV, 8.2 * eV, 8.4 * eV, 8.6 * eV,
1066  8.8 * eV, 9.0 * eV, 9.2 * eV, 9.4 * eV, 9.6 * eV,
1067  9.8 * eV, 10.0 * eV, 10.2 * eV, 10.4 * eV, 10.6 * eV,
1068  10.8 * eV, 11.0 * eV, 11.2 * eV, 11.3 * eV, 11.4 * eV,
1069  11.5 * eV, 11.6 * eV, 11.7 * eV, 11.8 * eV, 11.9 * eV};
1070 
1071  G4double RefractiveIndex_CsI[nEntries_CsI] =
1072  {1., 1., 1., 1., 1.,
1073  1., 1., 1., 1., 1.,
1074  1., 1., 1., 1., 1.,
1075  1., 1., 1., 1., 1.,
1076  1., 1., 1., 1., 1.,
1077  1., 1., 1., 1., 1.,
1078  1., 1., 1., 1., 1.,
1079  1., 1., 1., 1., 1.,
1080  1., 1., 1., 1., 1.,
1081  1., 1., 1., 1., 1.};
1082 
1083  G4double Absorption_CsI[nEntries_CsI] =
1084  {0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m,
1085  0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m,
1086  0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m,
1087  0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m,
1088  0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m,
1089  0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m,
1090  0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m,
1091  0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m,
1092  0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m,
1093  0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m, 0.0000001 * m};
1094 
1095  G4MaterialPropertiesTable *MPT_CsI = new G4MaterialPropertiesTable();
1096 
1097  MPT_CsI->AddProperty("RINDEX", PhotonEnergy_CsI, RefractiveIndex_CsI, nEntries_CsI)
1098  ->SetSpline(true);
1099  MPT_CsI->AddProperty("ABSLENGTH", PhotonEnergy_CsI, Absorption_CsI, nEntries_CsI)
1100  ->SetSpline(true);
1101 
1102  CsI->SetMaterialPropertiesTable(MPT_CsI);
1103 
1104  // define P10 Gas which will be used for TPC Benchmarking
1105  G4Material *P10 =
1106  new G4Material("P10", density = 1.74 * mg / cm3, ncomponents = 3); // @ 0K, 1atm
1107  P10->AddElement(G4Element::GetElement("Ar"), fractionmass = 0.9222);
1108  P10->AddElement(G4Element::GetElement("C"), fractionmass = 0.0623);
1109  P10->AddElement(G4Element::GetElement("H"), fractionmass = 0.0155);
1110 
1111  //------------------------------
1112  // for Modular RICH (mRICH)
1113  //------------------------------
1114 
1115  // mRICH_Air_Opt ---------------
1116  const int mRICH_nEntries_Air_Opt=2;
1117 
1118  G4double mRICH_PhotonEnergy_Air_Opt[mRICH_nEntries_Air_Opt] = {2.034*eV, 4.136*eV};
1119  G4double mRICH_RefractiveIndex_Air_Opt[mRICH_nEntries_Air_Opt] = {1.00, 1.00};
1120 
1121  G4MaterialPropertiesTable* mRICH_Air_Opt_MPT = new G4MaterialPropertiesTable();
1122  mRICH_Air_Opt_MPT->AddProperty("RINDEX", mRICH_PhotonEnergy_Air_Opt, mRICH_RefractiveIndex_Air_Opt, mRICH_nEntries_Air_Opt);
1123 
1124  G4Material* mRICH_Air_Opt = new G4Material("mRICH_Air_Opt", density=1.29*mg/cm3, ncomponents=2);
1125  mRICH_Air_Opt->AddElement(G4Element::GetElement("N"), fractionmass=70.*perCent);
1126  mRICH_Air_Opt->AddElement(G4Element::GetElement("O"), fractionmass=30.*perCent);
1127  mRICH_Air_Opt->SetMaterialPropertiesTable(mRICH_Air_Opt_MPT);
1128 
1129  // mRICH_Acrylic ---------------
1130  const int mRICH_nEntries1=32;
1131 
1132  //same photon energy array as aerogel 1
1133  G4double mRICH_PhotonEnergy[mRICH_nEntries1] =
1134  { 2.034*eV, 2.068*eV, 2.103*eV, 2.139*eV, // 610, 600, 590, 580, (nm)
1135  2.177*eV, 2.216*eV, 2.256*eV, 2.298*eV, // 570, 560, 550, 540,
1136  2.341*eV, 2.386*eV, 2.433*eV, 2.481*eV, // 530, 520, 510, 500,
1137  2.532*eV, 2.585*eV, 2.640*eV, 2.697*eV, // 490, 480, 470, 460,
1138  2.757*eV, 2.820*eV, 2.885*eV, 2.954*eV, // 450, 440, 430, 420,
1139  3.026*eV, 3.102*eV, 3.181*eV, 3.265*eV, // 410, 400, 390, 380,
1140  3.353*eV, 3.446*eV, 3.545*eV, 3.649*eV, // 370, 360, 350, 340,
1141  3.760*eV, 3.877*eV, 4.002*eV, 4.136*eV }; // 330, 320, 310, 300.
1142 
1143  G4double mRICH_AcRefractiveIndex[mRICH_nEntries1] =
1144  { 1.4902, 1.4907, 1.4913, 1.4918, 1.4924, // 610, 600, 590, 580, 570,
1145  1.4930, 1.4936, 1.4942, 1.4948, 1.4954, // 560, 550, 540, 530, 520, (this line is interpolated)
1146  1.4960, 1.4965, 1.4971, 1.4977, 1.4983, // 510, 500, 490, 480, 470,
1147  1.4991, 1.5002, 1.5017, 1.5017, 1.5017, // 460, 450, 440, 430, 420,
1148  1.5017, 1.5017, 1.5017, 1.5017, 1.5017, // 410,
1149  1.5017, 1.5017, 1.5017, 1.5017, 1.5017, // 360, look up values below 435
1150  1.5017, 1.5017};
1151 
1152  G4double mRICH_AcAbsorption[mRICH_nEntries1] =
1153  {25.25*cm, 25.25*cm, 25.25*cm, 25.25*cm,
1154  25.25*cm, 25.25*cm, 25.25*cm, 25.25*cm,
1155  25.25*cm, 25.25*cm, 25.25*cm, 25.25*cm,
1156  25.25*cm, 25.25*cm, 25.25*cm, 25.25*cm,
1157  25.25*cm, 25.25*cm, 25.25*cm, 25.25*cm,
1158  25.25*cm, 00.667*cm, 00.037*cm, 00.333*cm,
1159  00.001*cm, 00.001*cm, 00.001*cm, 00.001*cm,
1160  00.001*cm, 00.001*cm, 00.001*cm, 00.001*cm};
1161 
1162  G4MaterialPropertiesTable* mRICH_Ac_myMPT = new G4MaterialPropertiesTable();
1163  mRICH_Ac_myMPT->AddProperty("RINDEX" ,mRICH_PhotonEnergy, mRICH_AcRefractiveIndex,mRICH_nEntries1);
1164  mRICH_Ac_myMPT->AddProperty("ABSLENGTH", mRICH_PhotonEnergy, mRICH_AcAbsorption,mRICH_nEntries1);
1165 
1166  G4Material* mRICH_Acrylic=new G4Material("mRICH_Acrylic", density=1.19*g/cm3, ncomponents=3);
1167  mRICH_Acrylic->AddElement(G4Element::GetElement("C"), natoms=5);
1168  mRICH_Acrylic->AddElement(G4Element::GetElement("H"), natoms=8); // molecular ratios
1169  mRICH_Acrylic->AddElement(G4Element::GetElement("O"), natoms=2);
1170  mRICH_Acrylic->SetMaterialPropertiesTable(mRICH_Ac_myMPT);
1171 
1172  // mRICH_Agel1 -----------------
1173  // using same photon energy array as mRICH_Acrylic
1174 
1175  G4double mRICH_Agel1RefractiveIndex[mRICH_nEntries1] =
1176  { 1.02435, 1.0244, 1.02445, 1.0245, 1.02455,
1177  1.0246, 1.02465, 1.0247, 1.02475, 1.0248,
1178  1.02485, 1.02492, 1.025, 1.02505, 1.0251,
1179  1.02518, 1.02522, 1.02530, 1.02535, 1.0254,
1180  1.02545, 1.0255, 1.02555, 1.0256, 1.02568,
1181  1.02572, 1.0258, 1.02585, 1.0259, 1.02595,
1182  1.026, 1.02608};
1183 
1184  G4double mRICH_Agel1Absorption[mRICH_nEntries1] = //from Hubert
1185  { 3.448*m, 4.082*m, 6.329*m, 9.174*m, 12.346*m, 13.889*m,
1186  15.152*m, 17.241*m, 18.868*m, 20.000*m, 26.316*m, 35.714*m,
1187  45.455*m, 47.619*m, 52.632*m, 52.632*m, 55.556*m, 52.632*m,
1188  52.632*m, 47.619*m, 45.455*m, 41.667*m, 37.037*m, 33.333*m,
1189  30.000*m, 28.500*m, 27.000*m, 24.500*m, 22.000*m, 19.500*m,
1190  17.500*m, 14.500*m };
1191 
1192  G4double mRICH_Agel1Rayleigh[mRICH_nEntries1];
1193  //const G4double AerogelTypeAClarity = 0.00719*micrometer*micrometer*micrometer*micrometer/cm;
1194  const G4double AerogelTypeAClarity = 0.0020*micrometer*micrometer*micrometer*micrometer/cm;
1195  G4double Cparam = AerogelTypeAClarity*cm/(micrometer*micrometer*micrometer*micrometer);
1196  G4double PhotMomWaveConv = 1239*eV*nm;
1197 
1198  if(Cparam != 0.0 ) {
1199  for(int i=0; i<mRICH_nEntries1; i++ ){
1200  G4double ephoton = mRICH_PhotonEnergy[i];
1201  //In the following the 1000 is to convert form nm to micrometer
1202  G4double wphoton=(PhotMomWaveConv/ephoton)/(1000.0*nm);
1203  mRICH_Agel1Rayleigh[i]=(std::pow(wphoton,4))/Cparam;
1204  }
1205  }
1206 
1207  G4MaterialPropertiesTable* mRICH_Agel1_myMPT = new G4MaterialPropertiesTable();
1208  mRICH_Agel1_myMPT->AddProperty("RINDEX", mRICH_PhotonEnergy, mRICH_Agel1RefractiveIndex, mRICH_nEntries1);
1209  mRICH_Agel1_myMPT->AddProperty("ABSLENGTH", mRICH_PhotonEnergy, mRICH_Agel1Absorption,mRICH_nEntries1);
1210  mRICH_Agel1_myMPT->AddProperty("RAYLEIGH", mRICH_PhotonEnergy, mRICH_Agel1Rayleigh, mRICH_nEntries1); //Need table of rayleigh Scattering!!!
1211  mRICH_Agel1_myMPT->AddConstProperty("SCINTILLATIONYIELD",0./MeV);
1212  mRICH_Agel1_myMPT->AddConstProperty("RESOLUTIONSCALE",1.0);
1213 
1214  G4Material* mRICH_Aerogel1 = new G4Material("mRICH_Aerogel1", density=0.02*g/cm3, ncomponents=2);
1215  mRICH_Aerogel1->AddElement(G4Element::GetElement("Si"), natoms=1);
1216  mRICH_Aerogel1->AddElement(G4Element::GetElement("O"), natoms=2);
1217  mRICH_Aerogel1->SetMaterialPropertiesTable(mRICH_Agel1_myMPT);
1218 
1219  // mRICH_Agel2 -----------------
1220  const int mRICH_nEntries2=50;
1221 
1222  G4double mRICH_Agel2PhotonEnergy[mRICH_nEntries2]=
1223  {1.87855*eV,1.96673*eV,2.05490*eV,2.14308*eV,2.23126*eV,
1224  2.31943*eV,2.40761*eV,2.49579*eV,2.58396*eV,2.67214*eV,
1225  2.76032*eV,2.84849*eV,2.93667*eV,3.02485*eV,3.11302*eV,
1226  3.20120*eV,3.28938*eV,3.37755*eV,3.46573*eV,3.55391*eV,
1227  3.64208*eV,3.73026*eV,3.81844*eV,3.90661*eV,3.99479*eV,
1228  4.08297*eV,4.17114*eV,4.25932*eV,4.34750*eV,4.43567*eV,
1229  4.52385*eV,4.61203*eV,4.70020*eV,4.78838*eV,4.87656*eV,
1230  4.96473*eV,5.05291*eV,5.14109*eV,5.22927*eV,5.31744*eV,
1231  5.40562*eV,5.49380*eV,5.58197*eV,5.67015*eV,5.75833*eV,
1232  5.84650*eV,5.93468*eV,6.02286*eV,6.11103*eV,6.19921*eV };
1233 
1234  G4double mRICH_Agel2RefractiveIndex[mRICH_nEntries2] =
1235  {1.02825,1.02829,1.02834,1.02839,1.02844,
1236  1.02849,1.02854,1.02860,1.02866,1.02872,
1237  1.02878,1.02885,1.02892,1.02899,1.02906,
1238  1.02914,1.02921,1.02929,1.02938,1.02946,
1239  1.02955,1.02964,1.02974,1.02983,1.02993,
1240  1.03003,1.03014,1.03025,1.03036,1.03047,
1241  1.03059,1.03071,1.03084,1.03096,1.03109,
1242  1.03123,1.03137,1.03151,1.03166,1.03181,
1243  1.03196,1.03212,1.03228,1.03244,1.03261,
1244  1.03279,1.03297,1.03315,1.03334,1.03354};
1245 
1246  G4double mRICH_Agel2Absorption[mRICH_nEntries2] = //from Marco
1247  {17.5000*cm,17.7466*cm,17.9720*cm,18.1789*cm,18.3694*cm,
1248  18.5455*cm,18.7086*cm,18.8602*cm,19.0015*cm,19.1334*cm,
1249  19.2569*cm,19.3728*cm,19.4817*cm,19.5843*cm,19.6810*cm,
1250  19.7725*cm,19.8590*cm,19.9410*cm,20.0188*cm,20.0928*cm,
1251  18.4895*cm,16.0174*cm,13.9223*cm,12.1401*cm,10.6185*cm,
1252  9.3147*cm,8.1940*cm,7.2274*cm,6.3913*cm,5.6659*cm,
1253  5.0347*cm,4.4841*cm,4.0024*cm,3.5801*cm,3.2088*cm,
1254  2.8817*cm,2.5928*cm,2.3372*cm,2.1105*cm,1.9090*cm,
1255  1.7296*cm,1.5696*cm,1.4266*cm,1.2986*cm,1.1837*cm,
1256  1.0806*cm,0.9877*cm,0.9041*cm,0.8286*cm,0.7603*cm };
1257 
1258  G4double mRICH_Agel2Rayleigh[mRICH_nEntries2] = //from Marco
1259  {35.1384*cm, 29.24805*cm, 24.5418*cm, 20.7453*cm, 17.6553*cm,
1260  15.1197*cm, 13.02345*cm, 11.2782*cm, 9.81585*cm, 8.58285*cm,
1261  7.53765*cm, 6.6468*cm, 5.88375*cm, 5.22705*cm, 4.6596*cm,
1262  4.167*cm, 3.73785*cm, 3.36255*cm, 3.03315*cm, 2.7432*cm,
1263  2.487*cm, 2.26005*cm, 2.05845*cm, 1.87875*cm, 1.71825*cm,
1264  1.57455*cm, 1.44555*cm, 1.3296*cm, 1.2249*cm, 1.1304*cm,
1265  1.04475*cm, 0.9672*cm, 0.89655*cm, 0.83235*cm, 0.77385*cm,
1266  0.7203*cm, 0.67125*cm, 0.6264*cm, 0.58515*cm, 0.54735*cm,
1267  0.51255*cm, 0.48045*cm, 0.45075*cm, 0.4233*cm, 0.39795*cm,
1268  0.37455*cm, 0.3528*cm, 0.33255*cm, 0.3138*cm, 0.29625*cm};
1269 
1270  G4MaterialPropertiesTable* mRICH_Agel2MPT = new G4MaterialPropertiesTable();
1271  mRICH_Agel2MPT->AddProperty("RINDEX", mRICH_Agel2PhotonEnergy, mRICH_Agel2RefractiveIndex, mRICH_nEntries2);
1272  mRICH_Agel2MPT->AddProperty("ABSLENGTH", mRICH_Agel2PhotonEnergy, mRICH_Agel2Absorption,mRICH_nEntries2);
1273  mRICH_Agel2MPT->AddProperty("RAYLEIGH", mRICH_Agel2PhotonEnergy, mRICH_Agel2Rayleigh, mRICH_nEntries2);
1274  mRICH_Agel2MPT->AddConstProperty("SCINTILLATIONYIELD",0./MeV);
1275  mRICH_Agel2MPT->AddConstProperty("RESOLUTIONSCALE",1.0);
1276 
1277  G4Material* mRICH_Aerogel2 = new G4Material("mRICH_Aerogel2", density=0.02*g/cm3, ncomponents=2);
1278  mRICH_Aerogel2->AddElement(G4Element::GetElement("Si"), natoms=1);
1279  mRICH_Aerogel2->AddElement(G4Element::GetElement("O"), natoms=2);
1280  mRICH_Aerogel2->SetMaterialPropertiesTable(mRICH_Agel2MPT);
1281 
1282  // mRICH_Borosilicate ----------
1283  //using the same photon energy array as mRICH_Acrylic
1284 
1285  G4double mRICH_glassRefractiveIndex[mRICH_nEntries1] =
1286  { 1.47, 1.47, 1.47, 1.47, 1.47,
1287  1.47, 1.47, 1.47, 1.47, 1.47,
1288  1.47, 1.47, 1.47, 1.47, 1.47,
1289  1.47, 1.47, 1.47, 1.47, 1.47,
1290  1.47, 1.47, 1.47, 1.47, 1.47,
1291  1.47, 1.47, 1.47, 1.47, 1.47,
1292  1.47, 1.47};
1293 
1294  G4double mRICH_glassAbsorption[mRICH_nEntries1] =
1295  {4.25*cm, 4.25*cm, 4.25*cm, 4.25*cm,
1296  4.25*cm, 4.25*cm, 4.25*cm, 4.25*cm,
1297  4.25*cm, 4.25*cm, 4.25*cm, 4.25*cm,
1298  4.25*cm, 4.25*cm, 4.25*cm, 4.25*cm,
1299  4.25*cm, 4.25*cm, 4.25*cm, 4.25*cm,
1300  4.25*cm, 00.667*cm, 00.037*cm, 00.333*cm,
1301  00.001*cm, 00.001*cm, 00.001*cm, 00.001*cm,
1302  00.001*cm, 00.001*cm, 00.001*cm, 00.001*cm};
1303 
1304  G4MaterialPropertiesTable* mRICH_glass_myMPT = new G4MaterialPropertiesTable();
1305  //same photon energy array as aerogel 1
1306  mRICH_glass_myMPT->AddProperty("RINDEX", mRICH_PhotonEnergy, mRICH_glassRefractiveIndex, mRICH_nEntries1);
1307  mRICH_glass_myMPT->AddProperty("ABSLENGTH", mRICH_PhotonEnergy, mRICH_glassAbsorption, mRICH_nEntries1);
1308 
1309  const G4Material *G4_Pyrex_Glass = nist->FindOrBuildMaterial("G4_Pyrex_Glass");
1310  G4Material *mRICH_Borosilicate= new G4Material("mRICH_Borosilicate",G4_Pyrex_Glass->GetDensity(),G4_Pyrex_Glass,G4_Pyrex_Glass->GetState(),G4_Pyrex_Glass->GetTemperature(),G4_Pyrex_Glass->GetPressure());
1311  mRICH_Borosilicate->SetMaterialPropertiesTable(mRICH_glass_myMPT);
1312 
1313  // mRICH_Air ------------------
1314  //using same photon energy array as mRICH_Acrylic
1315 
1316  G4double mRICH_AirRefractiveIndex[mRICH_nEntries1] =
1317  { 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00,
1318  1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00,
1319  1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00,
1320  1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00,
1321  1.00, 1.00, 1.00, 1.00 };
1322 
1323  G4MaterialPropertiesTable* mRICH_Air_myMPT = new G4MaterialPropertiesTable();
1324  mRICH_Air_myMPT->AddProperty("RINDEX", mRICH_PhotonEnergy, mRICH_AirRefractiveIndex , mRICH_nEntries1);
1325  const G4Material *G4_AIR = G4Material::GetMaterial("G4_AIR");
1326  G4Material *mRICH_Air= new G4Material("mRICH_Air",G4_AIR->GetDensity(),G4_AIR,G4_AIR->GetState(),G4_AIR->GetTemperature(),G4_AIR->GetPressure());
1327  mRICH_Air->SetMaterialPropertiesTable(mRICH_Air_myMPT);
1328 
1329  //E1039 materials
1330  G4Element* eleD = new G4Element("eleD", "eleD", 1.0, 2.014);
1331 
1332  G4Material* SQ_LH2 = new G4Material("SQ_LH2", density = 0.07066*g/cm3, ncomponents = 1);
1333  SQ_LH2->AddElement(G4Element::GetElement("H"), natoms = 2);
1334 
1335  G4Material* SQ_LD2 = new G4Material("SQ_LD2", density = 0.15707*g/cm3, ncomponents = 1);
1336  SQ_LD2->AddElement(eleD, natoms = 2);
1337 
1338  G4Material* SQ_NH3 = new G4Material("SQ_NH3", density = 0.917*g/cm3, ncomponents = 2);
1339  SQ_NH3->AddElement(G4Element::GetElement("N"), natoms = 1);
1340  SQ_NH3->AddElement(G4Element::GetElement("H"), natoms = 3);
1341 
1342  G4Material* SQ_ND3 = new G4Material("SQ_ND3", density = 0.987*g/cm3, ncomponents = 2);
1343  SQ_ND3->AddElement(G4Element::GetElement("N"), natoms = 1);
1344  SQ_ND3->AddElement(eleD, natoms = 3);
1345 
1346  G4Material* SQ_boratedPoly = new G4Material("SQ_boratedPoly", density = 0.95000146*g/cm3, ncomponents = 3);
1347  SQ_boratedPoly->AddElement(G4Element::GetElement("B"), fractionmass = 0.7032045573308*perCent);
1348  SQ_boratedPoly->AddElement(G4Element::GetElement("C"), fractionmass = 4.0128564211902*perCent);
1349  SQ_boratedPoly->AddElement(G4Element::GetElement("H"), fractionmass = 95.283939021478*perCent);
1350 
1351  G4Material* SQ_Scintillator = new G4Material("SQ_Scintillator", density = 1.02300158*g/cm3, ncomponents = 2);
1352  SQ_Scintillator->AddElement(G4Element::GetElement("C"), fractionmass = 47.368421052631*perCent);
1353  SQ_Scintillator->AddElement(G4Element::GetElement("H"), fractionmass = 52.631578947368*perCent);
1354 
1355  G4Material* SQ_ArCO2 = new G4Material("SQ_ArCO2", density = 0.001822, ncomponents = 3);
1356  SQ_ArCO2->AddElement(G4Element::GetElement("Ar"), fractionmass = 77.651981*perCent);
1357  SQ_ArCO2->AddElement(G4Element::GetElement("C"), fractionmass = 5.83366*perCent);
1358  SQ_ArCO2->AddElement(G4Element::GetElement("O"), fractionmass = 16.514359*perCent);
1359 
1360  G4Material* SQ_P08CF4 = new G4Material("SQ_P08CF4", density = 0.00177504, ncomponents = 4);
1361  SQ_P08CF4->AddElement(G4Element::GetElement("Ar"), fractionmass = 89.121162862311*perCent);
1362  SQ_P08CF4->AddElement(G4Element::GetElement("C"), fractionmass = 3.2270886143930*perCent);
1363  SQ_P08CF4->AddElement(G4Element::GetElement("H"), fractionmass = 0.7168203310375*perCent);
1364  SQ_P08CF4->AddElement(G4Element::GetElement("F"), fractionmass = 6.9349281922577*perCent);
1365 }
1366 
1367 void
1369 {
1370  const G4RegionStore* theRegionStore = G4RegionStore::GetInstance();
1371  G4ProductionCuts *gcuts = new G4ProductionCuts(*(theRegionStore->GetRegion("DefaultRegionForTheWorld")->GetProductionCuts()));
1372  G4Region *tpcregion = new G4Region("REGION_TPCGAS");
1373  tpcregion->SetProductionCuts(gcuts);
1374 #if G4VERSION_NUMBER >= 1033
1375 // Use this from the new G4 version 10.03 on
1376 // add the PAI model to the TPCGAS region
1377 // undocumented, painfully digged out with debugger by tracing what
1378 // is done for command "/process/em/AddPAIRegion all TPCGAS PAI"
1379  G4EmParameters *g4emparams = G4EmParameters::Instance();
1380  g4emparams->AddPAIModel("all", "REGION_TPCGAS", "PAI");
1381 #endif
1382  return;
1383 }
1384 
1385 PHG4Subsystem *
1386 PHG4Reco::getSubsystem(const string &name)
1387 {
1388  //BOOST_FOREACH (PHG4Subsystem *subsys, subsystems_)
1389  for (PHG4Subsystem *subsys : subsystems_)
1390  {
1391  if (subsys->Name() == name)
1392  {
1393  if (verbosity > 0)
1394  {
1395  cout << "Found Subsystem " << name << endl;
1396  }
1397  return subsys;
1398  }
1399  }
1400  cout << "Could not find Subsystem " << name << endl;
1401  return nullptr;
1402 }
1403 
1404 void
1406 {
1407  if (runManager_)
1408  {
1409  runManager_->SetVerboseLevel(i);
1410  }
1411 }
void g4guithread(void *ptr)
Definition: PHG4Reco.cc:733
static TThread * gui_thread
Definition: PHG4Reco.cc:107
G4UImanager * UImanager
Definition: PHG4Reco.cc:110
int verbosity
The verbosity level. 0 means not verbose at all.
Definition: Fun4AllBase.h:75
virtual int Verbosity() const
Gets the verbosity of this module.
Definition: Fun4AllBase.h:64
G4FieldManager * GetDummyFieldManager()
PHBoolean addNode(PHNode *)
PHFieldConfig_v1 implements field configuration information for input a field map file.
PHFieldConfig_v2 implements field configuration information for uniform field model.
PHFieldConfig_v3 implements field configuration information for input a field map file.
PHFieldConfig store field configuration information.
Definition: PHFieldConfig.h:20
static PHField * GetFieldMapNode(const PHFieldConfig *default_config=nullptr, PHCompositeNode *topNode=nullptr, const int verbosity=0)
Get transient PHField from DST nodes. If not found, make a new one based on default_config.
transient DST object for field storage and access
Definition: PHField.h:14
virtual double get_DoubleFlag(const std::string &name) const
Definition: PHFlag.cc:49
virtual void set_FloatFlag(const std::string &name, const float flag)
Definition: PHFlag.cc:111
virtual const std::string get_CharFlag(const std::string &flag) const
Definition: PHFlag.cc:13
static void Dump_GDML(const std::string &filename, G4VPhysicalVolume *vol, PHCompositeNode *topNode=nullptr)
save the current Geant4 geometry to GDML file. Reading PHG4GDMLConfig from topNode
virtual void identify(std::ostream &os=std::cout) const
Definition: PHG4InEvent.cc:134
this is the main detector construction class, passed to geant to construct the entire phenix detector
G4VPhysicalVolume * GetPhysicalVolume(void)
void AddDetector(PHG4Detector *detector, int zero_field=0)
register a detector. This is called in PHG4Reco::Init based on which detectors are found on the tree
void SetWorldSizeX(const G4double sx)
void Verbosity(int verb)
void SetWorldShape(const std::string &s)
void SetWorldSizeY(const G4double sy)
void SetWorldSizeZ(const G4double sz)
void SetZeroFieldManager(G4FieldManager *man)
void SetWorldMaterial(const std::string &s)
void SetZeroFieldStartZ(const G4double z)
void AddAction(PHG4EventAction *action)
register an action. This is called in PHG4Reco::Init based on which actions are found on the tree
void AddAction(PHG4SteppingAction *action)
register an action. This is called in PHG4Reco::Init based on which actions are found on the tree
void AddAction(PHG4TrackingAction *action)
register an action. This is called in PHG4Reco::Init based on which actions are found on the tree
void SetInEvent(PHG4InEvent *const inevt)
set top node (from where particle list is retrieved for passing to geant
void set_field_map()
Definition: PHG4Reco.cc:309
void DefineRegions()
Definition: PHG4Reco.cc:1368
double zero_field_line_
Definition: PHG4Reco.h:188
PHG4PhenixTrackingAction * trackingAction_
pointer to main tracking action
Definition: PHG4Reco.h:161
PHG4PhenixEventAction * eventAction_
pointer to main event action
Definition: PHG4Reco.h:155
double _eta_coverage
Definition: PHG4Reco.h:173
int ResetEvent(PHCompositeNode *)
Clean up after each event.
Definition: PHG4Reco.cc:663
int StartGui()
start the gui
Definition: PHG4Reco.cc:576
SubsystemList subsystems_
Definition: PHG4Reco.h:168
PHG4Reco(const std::string &name="PHG4RECO")
constructor
Definition: PHG4Reco.cc:116
std::string worldshape
Definition: PHG4Reco.h:176
double energy_threshold_
Definition: PHG4Reco.h:187
void Dump_GDML(const std::string &filename)
Definition: PHG4Reco.cc:562
int Init(PHCompositeNode *)
full initialization
Definition: PHG4Reco.cc:170
void setGeneratorAction(G4VUserPrimaryGeneratorAction *action)
Definition: PHG4Reco.cc:712
PHTimeServer::timer _timer
module timer.
Definition: PHG4Reco.h:183
PHG4PrimaryGeneratorAction * generatorAction_
event generator (read from PHG4INEVENT node)
Definition: PHG4Reco.h:164
int InitField(PHCompositeNode *topNode)
Definition: PHG4Reco.cc:268
float magfield
Definition: PHG4Reco.h:138
bool optical_photon_
speed up options
Definition: PHG4Reco.h:186
bool save_DST_geometry_
Definition: PHG4Reco.h:180
std::string physicslist
Definition: PHG4Reco.h:178
void set_rapidity_coverage(const double eta)
Definition: PHG4Reco.cc:721
PHG4PhenixDetector * detector_
pointer to detector
Definition: PHG4Reco.h:152
virtual ~PHG4Reco()
destructor
Definition: PHG4Reco.cc:153
int ApplyCommand(const std::string &cmd)
interface to G4 cmd interpreter
Definition: PHG4Reco.cc:568
G4VisManager * visManager
Definition: PHG4Reco.h:171
std::string worldmaterial
Definition: PHG4Reco.h:177
void DefineMaterials()
Definition: PHG4Reco.cc:753
double WorldSize[3]
Definition: PHG4Reco.h:140
PHFieldConfig::FieldConfigTypes mapdim
Definition: PHG4Reco.h:174
int InitRun(PHCompositeNode *topNode)
Definition: PHG4Reco.cc:322
void Print(const std::string &what=std::string()) const
print info
Definition: PHG4Reco.cc:679
G4TBMagneticFieldSetup * field_
magnetic field
Definition: PHG4Reco.h:143
int End(PHCompositeNode *)
end of run method
Definition: PHG4Reco.cc:674
int process_event(PHCompositeNode *)
event processing method
Definition: PHG4Reco.cc:602
int InitUImanager()
Definition: PHG4Reco.cc:588
static void G4Seed(const unsigned int i)
Definition: PHG4Reco.cc:727
PHG4PhenixSteppingAction * steppingAction_
pointer to main stepping action
Definition: PHG4Reco.h:158
float magfield_rescale
Definition: PHG4Reco.h:139
void G4Verbosity(const int i)
Definition: PHG4Reco.cc:1405
G4RunManager * runManager_
pointer to geant run manager
Definition: PHG4Reco.h:146
std::string fieldmapfile
Definition: PHG4Reco.h:175
int setupInputEventNodeReader(PHCompositeNode *)
Definition: PHG4Reco.cc:694
PHG4UIsession * uisession_
pointer to geant ui session
Definition: PHG4Reco.h:149
PHG4Subsystem * getSubsystem(const std::string &name)
Definition: PHG4Reco.cc:1386
void Verbosity(int verb)
Definition: PHG4UIsession.h:12
static void SetPseudoRapidityCoverage(const double eta)
Definition: PHG4Utils.cc:27
static std::string GenerateGeometryFileName(const std::string &filename_extension="gdml")
static int ImportGeomFile(PHCompositeNode *topNode, const std::string &geometry_file)
TGeo ROOT/GDML/Macro file -> DST node with automatic file type discrimination based on file names.
static bool RemoveGeometryFile(const std::string &file_name)
delete the geometry file after use
PHNode * findFirst(const std::string &, const std::string &)
PHTimer * get(void)
Definition: PHTimeServer.h:43
PHTimer server for accessing external information.
Definition: PHTimeServer.h:29
void restart()
Restart timer.
Definition: PHTimer.h:70
void stop()
stops the counter
Definition: PHTimer.h:61
static recoConsts * instance()
Definition: recoConsts.cc:7
string cmd
Definition: submit_bnl.py:63
#define PHWHERE
Definition: phool.h:23