FairRoot/PandaRoot
gem_material_sim1.C
Go to the documentation of this file.
1 //root macro to run the simulation with RadLen to measure Radiation Length for GEM materials
2 // After this , run gem_material_ana1.C to draw the histograms , ( x vs y vs ELoss )
3 //
4 Int_t gem_material_sim1( Double_t momentum = 15., Int_t nEvents = 1000000, int verboseLevel = 0)
5  {
6  //FileNames
7  TString OutputFile = "gem_mat_sim_GemPion.root";
8  //TString OutputFile = "gem_mat_sim_GemPipe.root";
9  //TString OutputFile = "gem_mat_sim_withoutpipe.root";
10  //TString OutputFile = "gem_mat_sim_GPSM.root";
11  TString ParOutputFile = "gem_mat_simparams_GemPion.root";
12  //TString ParOutputFile = "gem_mat_simparams_GemPipe.root";
13  // TString ParOutputFile = "gem_mat_simparams_withoutpipe.root";
14  //TString ParOutputFile = "gem_mat_simparams_GPSM.root";
15 
16  TString SimEngine = "TGeant4";
17  Double_t BeamMomentum = momentum;
18  TString MediaFile = "media_pnd.geo";
19  gDebug = 0;
20 
21  //------------------------------------------------------------------
22  TStopwatch timer;
23  timer.Start();
24 
25  // Create the Simulation run manager--------------------------------
26  FairRunSim *fRun = new FairRunSim();
27  fRun->SetName(SimEngine.Data() );
28  fRun->SetOutputFile(OutputFile.Data());
29  fRun->SetBeamMom(BeamMomentum);
30  fRun->SetMaterials(MediaFile.Data());
31  fRun->SetUseFairLinks(kTRUE);
32  FairRuntimeDb *rtdb=fRun->GetRuntimeDb();
33 
34  //---------------------Set Parameter output-------------------------
36  FairParRootFileIo* output=new FairParRootFileIo(kParameterMerged);
37  output->open(ParOutputFile.Data());
38  rtdb->setOutput(output);
39 
40  // Create and add detectors
41  //------------------------- CAVE -----------------------------
42  FairModule *Cave= new PndCave("CAVE");
43  Cave->SetGeometryFileName("pndcave.geo");
44  fRun->AddModule(Cave);
45  //------------------------- Magnet -------------------------------
46  //FairModule *Magnet= new PndMagnet("MAGNET");
47  //Magnet->SetGeometryFileName("FullSolenoid_V842.root");
48  //Magnet->SetGeometryFileName("FullSuperconductingSolenoid_v831.root");
49  //fRun->AddModule(Magnet);
50  //FairModule *Dipole= new PndMagnet("MAGNET");
51  //Dipole->SetGeometryFileName("dipole.geo");
52  //fRun->AddModule(Dipole);
53  //------------------------- Pipe --------------------------
54  //FairModule *Pipe= new PndPipe("PIPE");
55  //Pipe->SetGeometryFileName("beampipe_201309.root");
56  //fRun->AddModule(Pipe);
57  // FairModule *Pipe= new PndPipe("PIPE");
58  // Pipe->SetGeometryFileName("pipebeamtarget.geo");
59  // fRun->AddModule(Pipe);
60  //------------------------- STT -----------------
61  //FairDetector *Stt= new PndStt("STT", kTRUE);
62  //Stt->SetGeometryFileName("straws_skewed_blocks_35cm_pipe.geo");
63  //fRun->AddModule(Stt);
64  //------------------------- MVD -----------------
65  //FairDetector *Mvd = new PndMvdDetector("MVD", kTRUE);
66  //Mvd->SetGeometryFileName("Mvd-2.1_FullVersion.root");
67  //fRun->AddModule(Mvd);
68  //-------------------------GEM----------------------------------
69  FairDetector *Gem = new PndGemDetector("GEM", kTRUE);
70  Gem->SetGeometryFileName("gem_3Stations_realistic_v1.root");
71  Gem->SetVerboseLevel(0);
72  fRun->AddModule(Gem);
73 
74  //----------------- Event generator-----------------------------
75  FairPrimaryGenerator* primGen = new FairPrimaryGenerator();
76  fRun->SetGenerator(primGen);
77 
78  FairBoxGenerator* boxGen = new FairBoxGenerator(211,1);
79  //FairBoxGenerator* boxGen = new FairBoxGenerator(13,1); /// 13 = muon
80  //FairBoxGenerator* boxGen = new FairBoxGenerator(22,1); /// 22 = gamma
81  //FairBoxGenerator* boxGen = new FairBoxGenerator(0,1); /// 0 = geantino
82  boxGen->SetThetaRange(0.0001,45);
83  boxGen->SetPhiRange (0.,360.);
84  boxGen->SetPRange (0.05,10.0);
85  primGen->AddGenerator(boxGen);
86 
87  // PndDpmDirect* dpmGen = new PndDpmDirect(momentum,1); //0. - only inelastic, 2 only elastic, 1 both
88  // primGen->AddGenerator(dpmGen);
89 
90  //---------------------Create and Set the Field(s)--------------
91  PndMultiField *fField= new PndMultiField("FULL");
92  fRun->SetField(fField);
93  //-----------end of Bfield stuff---------------------------------
94 
95  //--------------- support event display?------------------------
96  fRun->SetStoreTraj(kFALSE);
97  fRun->SetRadLenRegister(kTRUE);
98  fRun->SetRadMapRegister(kTRUE);
99 
100  fRun->Init();
101  //------------------- Transport nEvents----------------------------
102  fRun->Run(nEvents);
103  rtdb->saveOutput();
104  rtdb->print();
105 
106  timer.Stop();
107  Double_t rtime = timer.RealTime();
108  Double_t ctime = timer.CpuTime();
109  printf("RealTime=%f seconds, CpuTime=%f seconds\n",rtime,ctime);
110 
111  return 1;
112 }
113 
PndMultiField * fField
Definition: sim_emc_apd.C:97
printf("RealTime=%f seconds, CpuTime=%f seconds\n", rtime, ctime)
int verboseLevel
Definition: Lars/runMvdSim.C:7
Bool_t kParameterMerged
Definition: sim_emc_apd.C:113
Int_t gem_material_sim1(Double_t momentum=15., Int_t nEvents=1000000, int verboseLevel=0)
double BeamMomentum
Definition: sim_ftof_stof.C:17
FairParRootFileIo * output
Definition: sim_emc_apd.C:120
FairPrimaryGenerator * primGen
Definition: sim_emc_apd.C:81
FairRunAna * fRun
Definition: hit_dirc.C:58
FairDetector * Gem
Definition: runJohan.C:71
Double_t
FairModule * Cave
Definition: sim_emc_apd.C:32
Int_t nEvents
Definition: hit_dirc.C:11
TStopwatch timer
Definition: hit_dirc.C:51
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
gDebug
Definition: sim_emc_apd.C:6
Double_t ctime
Definition: hit_dirc.C:114
FairBoxGenerator * boxGen
Definition: sim_emc_apd.C:85
Double_t rtime
Definition: hit_dirc.C:113
Definition: PndCave.h:8