EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
CbmRichTestSim.cxx
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file CbmRichTestSim.cxx
1 
8 #include "CbmRichTestSim.h"
9 
10 #include "CbmRichPoint.h"
11 #include "CbmGeoRichPar.h"
12 
13 #include "FairTask.h"
14 #include "FairRootManager.h"
15 #include "CbmMCTrack.h"
16 #include "FairRunAna.h"
17 #include "FairRuntimeDb.h"
18 #include "FairBaseParSet.h"
19 #include "FairGeoNode.h"
20 
21 #include "TClonesArray.h"
22 #include "TVector3.h"
23 #include "TH2D.h"
24 #include "TH1D.h"
25 
26 #include <map>
27 #include <iostream>
28 
29 using std::cout;
30 using std::endl;
31 using std::map;
32 
33 
35  FairTask("CbmRichTestSim"),
36  fNEvents(0),
37  fMCRichPointArray(NULL),
38  fMCTrackArray(NULL),
39 
40  fh_Det1ev(NULL),
41  fh_Det1ev_zoom(NULL),
42  fh_Detall(NULL),
43  fh_n_vs_p(NULL),
44  fh_v_el(NULL),
45 
46  fh_Nall(NULL),
47  fh_Nel(NULL),
48  fh_Nelprim(NULL),
49  fh_Npi(NULL),
50  fh_Nk(NULL),
51  fh_Nhad(NULL),
52 
53  fSensNodes(NULL),
54  fPar(NULL),
55  fDetZ(0.)
56 {
57  fh_Det1ev = new TH2D("fh_Det1ev","points in detector plane (1 event)",170,-170,170,250,-250,250);
58  fh_Det1ev_zoom = new TH2D("fh_Det1ev_zoom","points in detector plane (1 event, zoom in)",100,10,60,100,100,150);
59  fh_Detall = new TH2D("fh_Detall","points in detector plane (all events)",170,-170,170,250,-250,250);
60  fh_n_vs_p = new TH2D("fh_n_vs_p","Npoints versus momentum",150,0,15,100,0,400);
61  fh_v_el = new TH2D("fh_v_el","(y,z) of production vertex of electrons",230,0,460,300,-300,300);
62 
63  fh_Nall = new TH1D("fh_Nall","Number of all rings in RICH",150,0,150);
64  fh_Nel = new TH1D("fh_Nel","Number of electron rings in RICH",150,0,150);
65  fh_Nelprim = new TH1D("fh_Nelprim","Number of electron (STS>6) rings in RICH",150,0,150);
66  fh_Npi = new TH1D("fh_Npi","Number of pi rings in RICH",150,0,150);
67  fh_Nk = new TH1D("fh_Nk","Number of K rings in RICH",150,0,150);
68  fh_Nhad = new TH1D("fh_Nhad","Number of hadrons crossing PMT plane",50,0,50);
69 }
70 
72 {
73  // write histograms to a file
74  fh_Det1ev->Write();
75  fh_Det1ev_zoom->Write();
76  fh_Detall->Write();
77  fh_n_vs_p->Write();
78  fh_v_el->Write();
79 
80  fh_Nall->Write();
81  fh_Nel->Write();
82  fh_Nelprim->Write();
83  fh_Npi->Write();
84  fh_Nk->Write();
85  fh_Nhad->Write();
86 }
87 
89 {
90  // Get Base Container
92  FairRuntimeDb* rtdb=ana->GetRuntimeDb();
93  fPar=(CbmGeoRichPar*)(rtdb->getContainer("CbmGeoRichPar"));
94 }
95 
97 {
99  // get detector position (z):
100  FairGeoNode *det= dynamic_cast<FairGeoNode*> (fSensNodes->FindObject("rich1d#1"));
101  if ( NULL == det) Fatal("CbmRichTestSim::Init", "no RICH Geo Node found!");
102  FairGeoTransform* detTr = det->getLabTransform(); // detector position in labsystem
103  FairGeoVector detPosLab = detTr->getTranslation(); // ... in cm
104  FairGeoTransform detCen = det->getCenterPosition(); // center in Detector system
105  FairGeoVector detPosCen = detCen.getTranslation();
106  fDetZ = detPosLab.Z() + detPosCen.Z(); // z coordinate of photodetector (Labsystem, cm)
107 
109  if (NULL == ioman) { Fatal("CbmRichTestSim::Init()", "RootManager not instantised!");}
110 
111  fMCTrackArray = (TClonesArray*) ioman->GetObject("MCTrack");
112  if ( NULL == fMCTrackArray) { Fatal("CbmRichTestSim::Init()", "No MCTrack array!");}
113 
114  fMCRichPointArray = (TClonesArray*) ioman->GetObject("RichPoint");
115  if ( NULL == fMCRichPointArray) { Fatal("CbmRichTestSim::Init()", "No RichPoint array!");}
116 
117  return kSUCCESS;
118 }
119 
121  Option_t* option)
122 {
123  fNEvents++;
124 
125  // define some variables
126  TVector3 position;
127  TVector3 momentum;
128  TVector3 vertex;
129 
130  Int_t Nhad, Nring, Nel, Nel_prim, Npi, Nk;
131  Nhad = Nring = Nel = Nel_prim = Npi = Nk = 0;
132 
133  map<Int_t,Int_t> pointMap;
134 
135  Int_t nPoints = fMCRichPointArray->GetEntriesFast();
136  Int_t nMCTracks = fMCTrackArray->GetEntriesFast();
137  for(Int_t iPoint=0; iPoint<nPoints; iPoint++){
138  CbmRichPoint* pPoint= (CbmRichPoint*)fMCRichPointArray->At(iPoint);
139  if( NULL == pPoint ) continue;
140  pPoint->Position(position);
141  fh_Detall->Fill(position.X(),position.Y());
142 
143  if (fNEvents == 1) {
144  fh_Det1ev->Fill(position.X(),position.Y());
145  fh_Det1ev_zoom->Fill(position.X(),position.Y());
146  }
147  if (pPoint->GetTrackID() < 0 ) continue;
148  CbmMCTrack* pTrack = (CbmMCTrack*)fMCTrackArray->At(pPoint->GetTrackID());
149  Int_t gcode = pTrack->GetPdgCode();
150  if (gcode != 50000050) Nhad++;
151 
152  if (gcode == 50000050){
153  pTrack->GetStartVertex(vertex);
154  Int_t motherID = pTrack->GetMotherId();
155  if (motherID == -1) continue;
156  pTrack = (CbmMCTrack*)fMCTrackArray->At(motherID);
157  pointMap[motherID]++;
158  } // select cherenkov photons
159  } // iPoint
160 
161  for (Int_t iMCTrack=0; iMCTrack<nMCTracks; iMCTrack++) {
162  CbmMCTrack* pTrack = (CbmMCTrack*) fMCTrackArray->At(iMCTrack);
163  if ( NULL == pTrack ) continue;
164  Int_t Npoints = pointMap[iMCTrack];
165  if (Npoints){
166  Nring++;
167  pTrack->GetMomentum(momentum);
168  fh_n_vs_p->Fill(momentum.Mag(),Npoints);
169  pTrack->GetStartVertex(vertex);
170  Int_t gcode = pTrack->GetPdgCode();
171  if (TMath::Abs(gcode) == 11) {
172  Nel++;
173  fh_v_el->Fill(vertex.Z(),vertex.Y());
174  if (pTrack->GetNPoints(kSTS) > 5) Nel_prim++;
175  } // electrons
176  if (TMath::Abs(gcode) == 211) Npi++;
177  if (TMath::Abs(gcode) == 321) Nk++;
178  } // tracks with points in Rich
179  } // iMCTrack
180 
181  fh_Nall->Fill(Nring);
182  fh_Nel->Fill(Nel);
183  fh_Nelprim->Fill(Nel_prim);
184  fh_Npi->Fill(Npi);
185  fh_Nk->Fill(Nk);
186  fh_Nhad->Fill(Nhad);
187 
188  if (fVerbose){
189  cout << "--------------------------------------------------------------------------" << endl;
190  cout << "----------- Test Rich Simulation -----------" << endl;
191  cout << endl;
192  cout << " Number of particles in RICH detector ----- event number " << fNEvents << endl;
193  cout << " hadrons in RICH (no Cherenkov photons) = " << Nhad << endl;
194  cout << " total rings = " << Nring << endl;
195  cout << " electrons = " << Nel << endl;
196  cout << " electrons (STS>=6) = " << Nel_prim << endl;
197  cout << " pions = " << Npi << endl;
198  cout << " Kaons = " << Nk << endl;
199  cout << endl;
200  cout << "--------------------------------------------------------------------------" << endl;
201  }
202 }
203 
205 {
206  fMCTrackArray->Clear();
207  fMCRichPointArray->Clear();
208 }
209 
211