StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StEmbeddingQA.h
1 //****************************************************************************************************
2 // Class StEmbeddingQA
3 //
4 // See http://drupal.star.bnl.gov/STAR/comp/embedding/embedding-procedures/qa-documentation
5 // for instructions
6 //****************************************************************************************************
7 /****************************************************************************************************
8  * $Id: StEmbeddingQA.h,v 1.12 2016/10/27 15:50:06 zhux Exp $
9  * $Log: StEmbeddingQA.h,v $
10  * Revision 1.12 2016/10/27 15:50:06 zhux
11  * added an option to set the maximum pT cut, by Zachariah Miller
12  *
13  * Revision 1.11 2012/03/05 10:32:29 cpowell
14  * Functions added to cut on refMult
15  *
16  * Revision 1.10 2011/04/01 05:05:47 hmasui
17  * Track selections by StEmbeddingQAUtilities. Added 1/pt(RC)-1/pt(MC) vs pt, and pt dependent Ncommon vs NhitFit histograms
18  *
19  * Revision 1.9 2011/02/11 03:55:44 hmasui
20  * Change geantid type to integer
21  *
22  * Revision 1.8 2011/01/31 21:32:10 hmasui
23  * Modify histogram keys to TString to take into account parent geantid
24  *
25  * Revision 1.7 2010/07/12 21:29:40 hmasui
26  * Move isGeantIdOk() function into StEmbeddingQAUtilities
27  *
28  * Revision 1.6 2010/05/14 19:50:11 hmasui
29  * Add rapidity and trigger cuts.
30  *
31  * Revision 1.5 2010/04/24 20:21:18 hmasui
32  * Add geant process check for contaminated pairs
33  *
34  * Revision 1.4 2010/02/16 02:13:34 hmasui
35  * Add parent-parent geant id in the histogram name
36  *
37  * Revision 1.3 2010/01/28 21:50:37 hmasui
38  * Add Vx vs Vz and Vy vs Vz histograms.
39  *
40  * Revision 1.2 2010/01/26 17:46:31 hmasui
41  * Add histograms for eventid, runnumber, and number of particles
42  *
43  * Revision 1.1 2009/12/22 21:41:17 hmasui
44  * Change class name from StEmbeddingQAMaker to StEmbeddingQA
45  *
46  ****************************************************************************************************/
47 
48 #ifndef __StEmbeddingQA_h__
49 #define __StEmbeddingQA_h__
50 
51 #include <map>
52 #include <vector>
53 
54 #include "TMath.h"
55 #include "TString.h"
56 
57 #include "StEmbeddingQAUtilities.h"
58 
59 class TFile ;
60 class TH1 ;
61 class TH2 ;
62 class TH3 ;
63 class TObject ;
64 //class TTree ;
65 
66 class StContamPair ;
67 class StEmbeddingQAPair ;
68 class StEmbeddingQATrack ;
69 class StMiniMcEvent ;
70 class StMuDstMaker ;
71 class StMuEvent ;
72 class StMuTrack ;
73 class StTinyMcTrack ;
74 
75 //____________________________________________________________________________________________________
76 // Analyze either minimc trees for embedding QA or microDST for real data QA
78  public:
80  StEmbeddingQA() ;
81 
83  StEmbeddingQA(const Int_t year, const TString production, const Bool_t isSimulation = kTRUE);
84 
86  virtual ~StEmbeddingQA();
87 
88  void init() ;
89 
94  Bool_t book(const TString outputFileName = "");
95 
97  Bool_t make(const TString inputFileName, const Bool_t isSimulation = kTRUE);
98 
100  Bool_t run(const TString inputFileList) ;
101 
103  Bool_t runRealData(const TString inputFileList) ;
104 
106  Bool_t runEmbedding(const TString inputFileList) ;
107 
109  Bool_t end() const;
110 
113  void setZVertexCut(const Float_t vz) ;
114  void setRefMultMinCut(const Int_t refMultMin) ;
115  void setRefMultMaxCut(const Int_t refMultMax) ;
116 
119  void addTriggerIdCut(const UInt_t id) ;
120 
123  void setRapidityCut(const Float_t ycut) ;
124 
126  void setPtMax(Float_t ptmax) ;
127 
128  private:
129  const Int_t mYear ;
130  const TString mProduction ;
131  const Bool_t mIsSimulation ;
132  Float_t mPtMax ;
133 
134  void clear() ;
135 
136  // Fill histograms
137  Bool_t fillEmbedding(const TString inputFileName) ; // Fill embedding histograms
138  Bool_t fillRealData(const TString inputFileName) ; // Fill real data histograms
139 
141  void fillEmbeddingTracks(const StMiniMcEvent& mcevent, const Int_t categoryid, const Int_t itrk) ;
142 
144  StEmbeddingQATrack* getEmbeddingQATrack(const StMiniMcEvent& mcevent, const Int_t categoryid, const Int_t itrk) ;
145 
147  void fillRealTracks(const StMuTrack& track, const Int_t categoryid, const Int_t itrk);
148 
150  void fillHistograms(const StEmbeddingQATrack& track, const Int_t categoryid);
151 
153  void expandHistograms(const Int_t categoryid, const Int_t geantid, const Int_t parentid,
154  const Int_t parentparentid, const Int_t geantprocess);
155 
157  Bool_t pushBackGeantId(const Int_t categoryid, const Int_t geantid, const Int_t parentid,
158  const Int_t parentparentid, const Int_t geantprocess) ;
159 
161  Bool_t isZVertexOk(const StMiniMcEvent& mcevent) const ;
162 
164  Bool_t isRefMultOk(const StMiniMcEvent& mcevent) const ;
165 
167  Bool_t isTriggerOk(StMuEvent* event) const ;
168 
170  Int_t getNtrack(const Int_t categoryid, const StMiniMcEvent& mcevent) const ;
171 
174  TString getIdCollection(const Int_t geantid, const Int_t parentid, const Int_t parentparentid) const ;
175 
176  StMuDstMaker* mMuDstMaker ;
177 
178  TFile* mOutput ;
179 
180  // Event-wise histograms
181  Double_t mVz ;
182  TH1* mhVz ;
183  TH1* mhVzAccepted ;
184  TH1* mhRef ;
185  TH1* mhRefAccepted ;
186 
187  TH2* mhVyVx ;
188  TH2* mhVxVz ;
189  TH2* mhVyVz ;
190  TH1* mhdVx ;
191  TH1* mhdVy ;
192  TH1* mhdVz ;
193 
194  TH1* mhEventId ;
195  TH1* mhRunNumber ;
196  TH1* mhNParticles[StEmbeddingQAConst::mNCategory] ;
197 
200  //
201  // - Use MC momentum instead of reconstructed momentum (Update on Nov/13/2009)
202  // - Add p (reco) vs p (MC) (Update on Nov/13/2009)
203  std::vector<Int_t> mGeantId[StEmbeddingQAConst::mNCategory] ;
204  std::vector<TString> mGeantIdCollection ;
205 
207  TH1* mhGeantId[StEmbeddingQAConst::mNCategory];
208  std::map<TString, TH3*> mhNHit[StEmbeddingQAConst::mNCategory] ;
209  std::map<TString, TH3*> mhDca[StEmbeddingQAConst::mNCategory] ;
210  std::map<TString, TH2*> mhPtVsEta[StEmbeddingQAConst::mNCategory] ;
211  std::map<TString, TH2*> mhPtVsY[StEmbeddingQAConst::mNCategory] ;
212  std::map<TString, TH2*> mhPtVsPhi[StEmbeddingQAConst::mNCategory] ;
213  std::map<TString, TH2*> mhPtVsMom[StEmbeddingQAConst::mNCategory] ;
214  std::map<TString, TH2*> mhdPtVsPt[StEmbeddingQAConst::mNCategory] ;
215  std::map<TString, TH2*> mhdInvPtVsPt[StEmbeddingQAConst::mNCategory] ;
216  std::map<TString, TH2*> mhMomVsEta[StEmbeddingQAConst::mNCategory] ;
217  std::map<TString, TH2*> mhdEdxVsMomMc[StEmbeddingQAConst::mNCategory] ;
218  std::map<TString, TH2*> mhdEdxVsMomMcPidCut[StEmbeddingQAConst::mNCategory] ;
219  std::map<TString, TH2*> mhdEdxVsMomReco[StEmbeddingQAConst::mNCategory] ;
220  std::map<TString, TH2*> mhdEdxVsMomRecoPidCut[StEmbeddingQAConst::mNCategory] ;
221  std::map<TString, TH2*> mhRecoPVsMcP[StEmbeddingQAConst::mNCategory] ;
222  std::map<TString, TH3*> mhNCommonHitVsNHit[StEmbeddingQAConst::mNCategory] ;
223 
224  std::map<TString, TH2*> mhEtaVsPhi[StEmbeddingQAConst::mNCategory] ;
225  std::map<TString, TH2*> mhEtaVsVz[StEmbeddingQAConst::mNCategory] ;
226  std::map<TString, TH2*> mhYVsVz[StEmbeddingQAConst::mNCategory] ;
227 
228  ClassDef(StEmbeddingQA, 1);
229 };
230 
231 #endif
232 
233 
234 
void setRapidityCut(const Float_t ycut)
StEmbeddingQA()
Default constructor (default argument is 2007, P08ic)
Bool_t run(const TString inputFileList)
Either RunRealData or RunEmbedding according to the kIsSimulation flag.
void setPtMax(Float_t ptmax)
Set Maximum Range of pT histograms; binning = 10*ptmax.
Bool_t runRealData(const TString inputFileList)
Analyzer real data.
Bool_t book(const TString outputFileName="")
Initialization.
virtual ~StEmbeddingQA()
Destructor.
Bool_t runEmbedding(const TString inputFileList)
Analyzer embedding data.
Bool_t make(const TString inputFileName, const Bool_t isSimulation=kTRUE)
Either fillEmbedding or fillRealData according to the isSimulation flag.
void addTriggerIdCut(const UInt_t id)
Bool_t end() const
Close output file.
void setZVertexCut(const Float_t vz)