Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DCdEdxSelector.C
Go to the documentation of this file.
1 #define DCdEdxSelector_cxx
2 // The class definition in DCdEdxSelector.h has been generated automatically
3 // by the ROOT utility TTree::MakeSelector(). This class is derived
4 // from the ROOT class TSelector. For more information on the TSelector
5 // framework see $ROOTSYS/README/README.SELECTOR or the ROOT User Manual.
6 
7 // The following methods are defined in this file:
8 // Begin(): called every time a loop on the tree starts,
9 // a convenient place to create your histograms.
10 // SlaveBegin(): called after Begin(), when on PROOF called only on the
11 // slave servers.
12 // Process(): called for each event, in this function you decide what
13 // to read and fill your histograms.
14 // SlaveTerminate: called at the end of the loop on the tree, when on PROOF
15 // called only on the slave servers.
16 // Terminate(): called at the end of the loop on the tree,
17 // a convenient place to draw/fit your histograms.
18 //
19 // To use this file, try the following session on your Tree T:
20 //
21 // Root > T->Process("DCdEdxSelector.C")
22 // Root > T->Process("DCdEdxSelector.C","some options")
23 // Root > T->Process("DCdEdxSelector.C+")
24 //
25 
26 #include "DCdEdxSelector.h"
27 #include <TH2.h>
28 #include <TStyle.h>
29 
30 void DCdEdxSelector::Begin(TTree * /*tree*/)
31 {
32  // The Begin() function is called at the start of the query.
33  // When running with PROOF Begin() is only called on the client.
34  // The tree argument is deprecated (on PROOF 0 is passed).
35 
36  TString option = GetOption();
37  TF1 *locFunc;
38 
39  dRhoZoverA_CDC = 1.29409e-07/0.1535E-3;
40  dRhoZoverA_FDC = 1.3277e-07/0.1535E-3;
41 
42  dCalcFOMManuallyFlag = true;
43 // dParticleName = "Proton";
44 // dParticleName = "PiPlus";
45  dParticleName = "KPlus";
46 
47  dOutputFile = new TFile("dh_DCdEdxHists.root", "RECREATE");
48 
49  ddEdxMeanFunc_FDC_Proton = new TF1("dPID_dEdxMeanFunc_FDC_Proton", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
50  ddEdxMeanFunc_FDC_Proton->SetParameters(47.5481, -4.70207, 2.22134, -0.357701);
51 
52  ddEdxMeanFunc_CDC_Proton = new TF1("dPID_dEdxMeanFunc_CDC_Proton", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
53  ddEdxMeanFunc_CDC_Proton->SetParameters(40.3217, -4.33611, 2.36034, -0.381628);
54 
55  ddEdxMeanFunc_FDC_KPlus = new TF1("dPID_dEdxMeanFunc_FDC_KPlus", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
56  ddEdxMeanFunc_FDC_KPlus->SetParameters(11.1221, -2.65063, 1.50373, -0.0202496);
57 
58  ddEdxMeanFunc_CDC_KPlus = new TF1("dPID_dEdxMeanFunc_CDC_KPlus", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
59  ddEdxMeanFunc_CDC_KPlus->SetParameters(12.3025, -2.51467, 1.53035, -0.0120788);
60 
61  dSigmaFuncArray_FDC_Proton = new TObjArray(4);
62  locFunc = new TF1("dPID_dEdxSigmaFunc_FDC_Proton1", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
63  locFunc->SetParameters(15.0407, -4.60489, 0.501318, -0.0913598);
64  dSigmaFuncArray_FDC_Proton->AddAt(locFunc, 0);
65  locFunc = new TF1("dPID_dEdxSigmaFunc_FDC_Proton2", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
66  locFunc->SetParameters(26.6243, -5.92033, 0.384732, -0.0820965);
67  dSigmaFuncArray_FDC_Proton->AddAt(locFunc, 1);
68  locFunc = new TF1("dPID_dEdxSigmaFunc_FDC_Proton3", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
69  locFunc->SetParameters(10.0415, -4.956, 0.29667, -0.0551989);
70  dSigmaFuncArray_FDC_Proton->AddAt(locFunc, 2);
71  locFunc = new TF1("dPID_dEdxSigmaFunc_FDC_Proton4", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
72  locFunc->SetParameters(14.5991, -5.9726, 0.274239, -0.0546458);
73  dSigmaFuncArray_FDC_Proton->AddAt(locFunc, 3);
75  for(unsigned int loc_i = 0; loc_i < 4; loc_i++)
76  ddEdxSigmaNumHitsVector_FDC_Proton[loc_i] = 3 + 3*loc_i;
77 
78  dSigmaFuncArray_CDC_Proton = new TObjArray(6);
79  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_Proton1", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
80  locFunc->SetParameters(14.263, -3.91119, 0.471623, -0.0867709);
81  dSigmaFuncArray_CDC_Proton->AddAt(locFunc, 0);
82  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_Proton2", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
83  locFunc->SetParameters(10.5091, -4.2701, 0.638759, -0.191789);
84  dSigmaFuncArray_CDC_Proton->AddAt(locFunc, 1);
85  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_Proton3", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
86  locFunc->SetParameters(5.74486, -3.46166, 0.353456, -0.0594251);
87  dSigmaFuncArray_CDC_Proton->AddAt(locFunc, 2);
88  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_Proton4", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
89  locFunc->SetParameters(6.93564, -4.55521, 0.372744, -0.0809587);
90  dSigmaFuncArray_CDC_Proton->AddAt(locFunc, 3);
91  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_Proton5", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
92  locFunc->SetParameters(11.3186, -6.17783, 0.269247, -0.0525885);
93  dSigmaFuncArray_CDC_Proton->AddAt(locFunc, 4);
94  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_Proton6", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
95  locFunc->SetParameters(8.41072, -5.90604, 0.23806, -0.0406181);
96  dSigmaFuncArray_CDC_Proton->AddAt(locFunc, 5);
98  for(unsigned int loc_i = 0; loc_i < 6; loc_i++)
99  ddEdxSigmaNumHitsVector_CDC_Proton[loc_i] = 4 + 2*loc_i;
100 
101  dSigmaFuncArray_CDC_PiPlus = new TObjArray(6);
102  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_PiPlus1", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
103  locFunc->SetParameters(1.8644, -1.92051, 0.401207, -0.00510045);
104  dSigmaFuncArray_CDC_PiPlus->AddAt(locFunc, 0);
105  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_PiPlus2", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
106  locFunc->SetParameters(3.08596, -2.55591, 0.326127, -0.00439305);
107  dSigmaFuncArray_CDC_PiPlus->AddAt(locFunc, 1);
108  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_PiPlus3", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
109  locFunc->SetParameters(1.63183, -2.18999, 0.276187, -0.00277214);
110  dSigmaFuncArray_CDC_PiPlus->AddAt(locFunc, 2);
111  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_PiPlus4", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
112  locFunc->SetParameters(0.42353, -1.12885, 0.218266, 0.000475549);
113  dSigmaFuncArray_CDC_PiPlus->AddAt(locFunc, 3);
114  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_PiPlus5", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
115  locFunc->SetParameters(0.599551, -1.27373, 0.15483, 0.00125788);
116  dSigmaFuncArray_CDC_PiPlus->AddAt(locFunc, 4);
117  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_PiPlus6", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
118  locFunc->SetParameters(1.02319, -2.15814, 0.146108, 0.00126295);
119  dSigmaFuncArray_CDC_PiPlus->AddAt(locFunc, 5);
120  for(unsigned int loc_i = 0; loc_i < 6; loc_i++)
121  ddEdxSigmaNumHitsVector_CDC_PiPlus.push_back(4 + 2*loc_i);
122 
123  dSigmaFuncArray_FDC_PiPlus = new TObjArray(4);
124  locFunc = new TF1("dPID_dEdxSigmaFunc_FDC_PiPlus1", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
125  locFunc->SetParameters(0.972821, -1.13778, 0.28177, 0.00414764);
126  dSigmaFuncArray_FDC_PiPlus->AddAt(locFunc, 0);
127  locFunc = new TF1("dPID_dEdxSigmaFunc_FDC_PiPlus2", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
128  locFunc->SetParameters(0.72011, -1.13622, 0.205834, 0.00248226);
129  dSigmaFuncArray_FDC_PiPlus->AddAt(locFunc, 1);
130  locFunc = new TF1("dPID_dEdxSigmaFunc_FDC_PiPlus3", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
131  locFunc->SetParameters(0.783602, -1.58954, 0.175932, 0.00139313);
132  dSigmaFuncArray_FDC_PiPlus->AddAt(locFunc, 2);
133  locFunc = new TF1("dPID_dEdxSigmaFunc_FDC_PiPlus4", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
134  locFunc->SetParameters(0.666107, -1.50611, 0.144197, 0.0019479);
135  dSigmaFuncArray_FDC_PiPlus->AddAt(locFunc, 3);
136  for(unsigned int loc_i = 0; loc_i < 4; loc_i++)
137  ddEdxSigmaNumHitsVector_FDC_PiPlus.push_back(3 + 3*loc_i);
138 
139  //sigmas, k+
140  dSigmaFuncArray_CDC_KPlus = new TObjArray(6);
141  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_KPlus1", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
142  locFunc->SetParameters(10.3149, -3.48486, 0.590082, -0.0757512);
143  dSigmaFuncArray_CDC_KPlus->AddAt(locFunc, 0);
144  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_KPlus2", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
145  locFunc->SetParameters(10.9306, -3.49101, 0.443805, -0.0527794);
146  dSigmaFuncArray_CDC_KPlus->AddAt(locFunc, 1);
147  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_KPlus3", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
148  locFunc->SetParameters(7.71149, -3.19237, 0.319949, -0.0248138);
149  dSigmaFuncArray_CDC_KPlus->AddAt(locFunc, 2);
150  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_KPlus4", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
151  locFunc->SetParameters(2.70214, -2.51235, 0.224439, -0.00564491);
152  dSigmaFuncArray_CDC_KPlus->AddAt(locFunc, 3);
153  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_KPlus5", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
154  locFunc->SetParameters(1.12165, -2.16051, 0.155717, 0.000195505);
155  dSigmaFuncArray_CDC_KPlus->AddAt(locFunc, 4);
156  locFunc = new TF1("dPID_dEdxSigmaFunc_CDC_KPlus6", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
157  locFunc->SetParameters(1.60267, -2.62908, 0.153127, -0.000775554);
158  dSigmaFuncArray_CDC_KPlus->AddAt(locFunc, 5);
159  for(unsigned int loc_i = 0; loc_i < 6; loc_i++)
160  ddEdxSigmaNumHitsVector_CDC_KPlus.push_back(4 + 2*loc_i);
161 
162  dSigmaFuncArray_FDC_KPlus = new TObjArray(4);
163  locFunc = new TF1("dPID_dEdxSigmaFunc_FDC_KPlus1", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
164  locFunc->SetParameters(3.4624, -2.51387, 0.328651, -0.00906703);
165  dSigmaFuncArray_FDC_KPlus->AddAt(locFunc, 0);
166  locFunc = new TF1("dPID_dEdxSigmaFunc_FDC_KPlus2", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
167  locFunc->SetParameters(8.77765, -3.16076, 0.20651, 0.00175279);
168  dSigmaFuncArray_FDC_KPlus->AddAt(locFunc, 1);
169  locFunc = new TF1("dPID_dEdxSigmaFunc_FDC_KPlus3", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
170  locFunc->SetParameters(6.44965, -3.01797, 0.175224, -0.000185835);
171  dSigmaFuncArray_FDC_KPlus->AddAt(locFunc, 2);
172  locFunc = new TF1("dPID_dEdxSigmaFunc_FDC_KPlus4", "[0]*exp([1]*x) + [2] + [3]*x", 0.0, 100.0);
173  locFunc->SetParameters(1.37552, -2.21704, 0.148445, 0.0014621);
174  dSigmaFuncArray_FDC_KPlus->AddAt(locFunc, 3);
175  for(unsigned int loc_i = 0; loc_i < 4; loc_i++)
176  ddEdxSigmaNumHitsVector_FDC_KPlus.push_back(3 + 3*loc_i);
177 
178  unsigned int locNumBetaBins = 240;
179  unsigned int locNumdEdxBins = 5000;
180  unsigned int locNumMomentumBins = 500;
181  unsigned int locNumBetaGammaBins = 1500;
182  unsigned int locNumHitsBins = 20;
183  unsigned int locThetaBins = 600;
184  unsigned int locNumdxBins = 500;
185  float locBetaMin = 0.0, locBetaMax = 1.0;
186  float locdEdxMin = 0.0, locdEdxMax = 100.0;
187  float locMomentumMin = 0.0, locMomentumMax = 3.0;
188  float locBetaGammaMin = 0.0, locBetaGammaMax = 15.0;
189  float locNumHitsMin = -0.5, locNumHitsMax = 19.5;
190  float locThetaMin = 0.0, locThetaMax = 150.0;
191  float locdxMin = 0.0, locdxMax = 50.0;
192 
193  //no dx bins, no #hits cutoff
194  //also do: dx bins (probably very small dependence), different #hits cutoff
195  dSelectorHist_dEdxVsBeta_HitsCutoff_FDC = new TH2F("dSelectorHist_dEdxVsBeta_HitsCutoff_FDC", ";#beta;FDC #frac{dE}{dx} (keV/cm)", locNumBetaBins, locBetaMin, locBetaMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
196  dSelectorHist_dEdxVsP_HitsCutoff_FDC = new TH2F("dSelectorHist_dEdxVsP_HitsCutoff_FDC", ";p (GeV/c);FDC #frac{dE}{dx} (keV/cm)", locNumMomentumBins, locMomentumMin, locMomentumMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
197  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC = new TH2F("dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC", ";#beta#gamma;FDC #frac{dE}{dx} (keV/cm)", locNumBetaGammaBins, locBetaGammaMin, locBetaGammaMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
198  dSelectorHist_dEdxVsBeta_HitsCutoff_CDC = new TH2F("dSelectorHist_dEdxVsBeta_HitsCutoff_CDC", ";#beta;CDC #frac{dE}{dx} (keV/cm)", locNumBetaBins, locBetaMin, locBetaMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
199  dSelectorHist_dEdxVsP_HitsCutoff_CDC = new TH2F("dSelectorHist_dEdxVsP_HitsCutoff_CDC", ";p (GeV/c);CDC #frac{dE}{dx} (keV/cm)", locNumMomentumBins, locMomentumMin, locMomentumMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
200  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC = new TH2F("dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC", ";#beta#gamma;CDC #frac{dE}{dx} (keV/cm)", locNumBetaGammaBins, locBetaGammaMin, locBetaGammaMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
201 
202  dSelectorHist_NumHitsVsTheta_CDC = new TH2F("dSelectorHist_NumHitsVsTheta_CDC", ";#theta#circ;#CDC Hits / 2", locThetaBins, locThetaMin, locThetaMax, locNumHitsBins, locNumHitsMin, locNumHitsMax);
203  dSelectorHist_dxVsNumHits_CDC = new TH2F("dSelectorHist_dxVsNumHits_CDC", ";#CDC Hits / 2;dx (cm)", locNumHitsBins, locNumHitsMin, locNumHitsMax, locNumdxBins, locdxMin, locdxMax);
204  dSelectorHist_dxVsTheta_CDC = new TH2F("dSelectorHist_dxVsTheta_CDC", ";#theta#circ;CDC dx (cm)", locThetaBins, locThetaMin, locThetaMax, locNumdxBins, locdxMin, locdxMax);
205  dSelectorHist_NumHitsVsTheta_FDC = new TH2F("dSelectorHist_NumHitsVsTheta_FDC", ";#theta#circ;#FDC Hits / 2", locThetaBins, locThetaMin, locThetaMax, locNumHitsBins, locNumHitsMin, locNumHitsMax);
206  dSelectorHist_dxVsNumHits_FDC = new TH2F("dSelectorHist_dxVsNumHits_FDC", ";#FDC Hits / 2;dx (cm)", locNumHitsBins, locNumHitsMin, locNumHitsMax, locNumdxBins, locdxMin, locdxMax);
207  dSelectorHist_dxVsTheta_FDC = new TH2F("dSelectorHist_dxVsTheta_FDC", ";#theta#circ;FDC dx (cm)", locThetaBins, locThetaMin, locThetaMax, locNumdxBins, locdxMin, locdxMax);
208 
209  dSelectorHist_NumHitsFDCVsNumHitsCDC = new TH2F("dSelectorHist_NumHitsFDCVsNumHitsCDC", ";#FDC Hits / 2;#CDC Hits / 2", locNumHitsBins, locNumHitsMin, locNumHitsMax, locNumHitsBins, locNumHitsMin, locNumHitsMax);
210 
211 
212  dSelectorHist_NotEnoughHitsInTheta = new TH1F("dSelectorHist_NotEnoughHitsInTheta", "Not enough dE/dx Info;#theta#circ", locThetaBins, locThetaMin, locThetaMax);
213  dSelectorHist_ThetaDistribution = new TH1F("dSelectorHist_ThetaDistribution", "All Tracks;#theta#circ", locThetaBins, locThetaMin, locThetaMax);
214  dSelectorHist_PercentageNotEnoughHitsInTheta = new TH1F("dSelectorHist_PercentageNotEnoughHitsInTheta", "% Not enough dE/dx Info;#theta#circ", locThetaBins, locThetaMin, locThetaMax);
215 
216  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC_3Hits = new TH2F("dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC_3Hits", "# FDC Hits Used = 3;#beta#gamma;FDC #frac{dE}{dx} (keV/cm)", locNumBetaGammaBins, locBetaGammaMin, locBetaGammaMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
217  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC_6Hits = new TH2F("dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC_6Hits", "# FDC Hits Used = 6;#beta#gamma;FDC #frac{dE}{dx} (keV/cm)", locNumBetaGammaBins, locBetaGammaMin, locBetaGammaMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
218  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC_9Hits = new TH2F("dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC_9Hits", "# FDC Hits Used = 9;#beta#gamma;FDC #frac{dE}{dx} (keV/cm)", locNumBetaGammaBins, locBetaGammaMin, locBetaGammaMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
219  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC_12Hits = new TH2F("dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC_12Hits", "# FDC Hits Used = 12;#beta#gamma;FDC #frac{dE}{dx} (keV/cm)", locNumBetaGammaBins, locBetaGammaMin, locBetaGammaMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
220 
221  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_2Hits = new TH2F("dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_2Hits", "# CDC Hits Used = 2;#beta#gamma;CDC #frac{dE}{dx} (keV/cm)", locNumBetaGammaBins, locBetaGammaMin, locBetaGammaMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
222  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_4Hits = new TH2F("dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_4Hits", "# CDC Hits Used = 4;#beta#gamma;CDC #frac{dE}{dx} (keV/cm)", locNumBetaGammaBins, locBetaGammaMin, locBetaGammaMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
223  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_6Hits = new TH2F("dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_6Hits", "# CDC Hits Used = 6;#beta#gamma;CDC #frac{dE}{dx} (keV/cm)", locNumBetaGammaBins, locBetaGammaMin, locBetaGammaMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
224  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_8Hits = new TH2F("dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_8Hits", "# CDC Hits Used = 8;#beta#gamma;CDC #frac{dE}{dx} (keV/cm)", locNumBetaGammaBins, locBetaGammaMin, locBetaGammaMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
225  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_10Hits = new TH2F("dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_10Hits", "# CDC Hits Used = 10;#beta#gamma;CDC #frac{dE}{dx} (keV/cm)", locNumBetaGammaBins, locBetaGammaMin, locBetaGammaMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
226  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_12Hits = new TH2F("dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_12Hits", "# CDC Hits Used = 12;#beta#gamma;CDC #frac{dE}{dx} (keV/cm)", locNumBetaGammaBins, locBetaGammaMin, locBetaGammaMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
227  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_14Hits = new TH2F("dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_14Hits", "# CDC Hits Used = 14;#beta#gamma;CDC #frac{dE}{dx} (keV/cm)", locNumBetaGammaBins, locBetaGammaMin, locBetaGammaMax, locNumdEdxBins, locdEdxMin, locdEdxMax);
228 
229  //confidence level
230  unsigned int locNumConfidenceLevelBins = 400;
231  dSelectorHist_ConfidenceLevel_FDC = new TH1F("dSelectorHist_ConfidenceLevel_FDC", ";Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
232  dSelectorHist_ConfidenceLevel_CDC = new TH1F("dSelectorHist_ConfidenceLevel_CDC", ";Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
233  dSelectorHist_ConfidenceLevel_Both = new TH1F("dSelectorHist_ConfidenceLevel_Both", ";Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
234 
235  dSelectorHist_ConfidenceLevel_FDC_3Hits = new TH1F("dSelectorHist_ConfidenceLevel_FDC_3Hits", "# FDC Hits Used = 3;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
236  dSelectorHist_ConfidenceLevel_FDC_6Hits = new TH1F("dSelectorHist_ConfidenceLevel_FDC_6Hits", "# FDC Hits Used = 6;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
237  dSelectorHist_ConfidenceLevel_FDC_9Hits = new TH1F("dSelectorHist_ConfidenceLevel_FDC_9Hits", "# FDC Hits Used = 9;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
238  dSelectorHist_ConfidenceLevel_FDC_12Hits = new TH1F("dSelectorHist_ConfidenceLevel_FDC_12Hits", "# FDC Hits Used = 12;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
239 
240  dSelectorHist_ConfidenceLevel_CDC_4Hits = new TH1F("dSelectorHist_ConfidenceLevel_CDC_4Hits", "# CDC Hits Used = 4;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
241  dSelectorHist_ConfidenceLevel_CDC_6Hits = new TH1F("dSelectorHist_ConfidenceLevel_CDC_6Hits", "# CDC Hits Used = 6;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
242  dSelectorHist_ConfidenceLevel_CDC_8Hits = new TH1F("dSelectorHist_ConfidenceLevel_CDC_8Hits", "# CDC Hits Used = 8;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
243  dSelectorHist_ConfidenceLevel_CDC_10Hits = new TH1F("dSelectorHist_ConfidenceLevel_CDC_10Hits", "# CDC Hits Used = 10;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
244  dSelectorHist_ConfidenceLevel_CDC_12Hits = new TH1F("dSelectorHist_ConfidenceLevel_CDC_12Hits", "# CDC Hits Used = 12;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
245  dSelectorHist_ConfidenceLevel_CDC_14Hits = new TH1F("dSelectorHist_ConfidenceLevel_CDC_14Hits", "# CDC Hits Used = 14;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
246 
247  dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin1 = new TH1F("dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin1", "0.3 < #beta#gamma < 0.6;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
248  dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin2 = new TH1F("dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin2", "0.6 < #beta#gamma < 0.9;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
249  dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin3 = new TH1F("dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin3", "0.9 < #beta#gamma < 1.2;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
250  dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin4 = new TH1F("dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin4", "1.2 < #beta#gamma < 1.5;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
251  dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin5 = new TH1F("dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin5", "1.5 < #beta#gamma < 1.8;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
252  dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin6 = new TH1F("dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin6", "1.8 < #beta#gamma < 2.1;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
253  dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin7 = new TH1F("dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin7", "2.1 < #beta#gamma < 4.1;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
254  dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin8 = new TH1F("dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin8", "4.1 < #beta#gamma < 6.1;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
255  dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin9 = new TH1F("dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin9", "6.1 < #beta#gamma < 8.1;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
256  dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin10 = new TH1F("dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin10", "8.1 < #beta#gamma < 10.1;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
257  dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin11 = new TH1F("dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin11", "10.1 < #beta#gamma < 12.1;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
258  dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin12 = new TH1F("dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin12", "12.1 < #beta#gamma < 14.1;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
259 
260  dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin1 = new TH1F("dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin1", "0.3 < #beta#gamma < 0.6;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
261  dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin2 = new TH1F("dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin2", "0.6 < #beta#gamma < 0.9;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
262  dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin3 = new TH1F("dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin3", "0.9 < #beta#gamma < 1.2;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
263  dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin4 = new TH1F("dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin4", "1.2 < #beta#gamma < 1.5;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
264  dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin5 = new TH1F("dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin5", "1.5 < #beta#gamma < 1.8;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
265  dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin6 = new TH1F("dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin6", "1.8 < #beta#gamma < 2.1;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
266  dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin7 = new TH1F("dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin7", "2.1 < #beta#gamma < 4.1;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
267  dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin8 = new TH1F("dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin8", "4.1 < #beta#gamma < 6.1;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
268  dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin9 = new TH1F("dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin9", "6.1 < #beta#gamma < 8.1;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
269  dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin10 = new TH1F("dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin10", "8.1 < #beta#gamma < 10.1;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
270  dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin11 = new TH1F("dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin11", "10.1 < #beta#gamma < 12.1;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
271  dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin12 = new TH1F("dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin12", "12.1 < #beta#gamma < 14.1;Confidence Level", locNumConfidenceLevelBins, 0.0, 1.0);
272 }
273 
274 void DCdEdxSelector::SlaveBegin(TTree * /*tree*/)
275 {
276  // The SlaveBegin() function is called after the Begin() function.
277  // When running with PROOF SlaveBegin() is called on each slave server.
278  // The tree argument is deprecated (on PROOF 0 is passed).
279 
280  TString option = GetOption();
281 
282 }
283 
284 Bool_t DCdEdxSelector::Process(Long64_t entry)
285 {
286  // The Process() function is called for each entry in the tree (or possibly
287  // keyed object in the case of PROOF) to be processed. The entry argument
288  // specifies which entry in the currently loaded tree is to be processed.
289  // It can be passed to either DCdEdxSelector::GetEntry() or TBranch::GetEntry()
290  // to read either all or the required parts of the data. When processing
291  // keyed objects with PROOF, the object is already loaded and is available
292  // via the fObject pointer.
293  //
294  // This function should contain the "body" of the analysis. It can contain
295  // simple or elaborate selection criteria, run algorithms on the data
296  // of the event and typically fill histograms.
297  //
298  // The processing can be stopped by calling Abort().
299  //
300  // Use fStatus to set the return value of TTree::Process().
301  //
302  // The return value is currently not used.
303 
304  GetEntry(entry);
305  if((entry % 1000) == 0)
306  cout << "entry = " << entry << endl;
307 
308  double ddEdxkeV_FDC = 1000000.0*ddEdx_FDC; //convert GeV to keV (gas!)
309  double ddEdxkeV_CDC = 1000000.0*ddEdx_CDC; //convert GeV to keV (gas!)
310 
311  float locBetaGamma = dBeta/sqrt(1.0 - dBeta*dBeta);
312  float locTheta_Degrees = dTheta*180.0/3.141592654;
313 
314  unsigned int locMinimumNumHitsUsed = 3; //minimum num hits = 2x this
315 
318  dSelectorHist_dxVsTheta_CDC->Fill(locTheta_Degrees, ddx_CDC);
321  dSelectorHist_dxVsTheta_FDC->Fill(locTheta_Degrees, ddx_FDC);
322 
324 
325  if((dNumHitsUsedFordEdx_FDC < locMinimumNumHitsUsed) && (dNumHitsUsedFordEdx_CDC < locMinimumNumHitsUsed)) //not enough hits (6+) in either!
326  dSelectorHist_NotEnoughHitsInTheta->Fill(locTheta_Degrees);
327  dSelectorHist_ThetaDistribution->Fill(locTheta_Degrees);
328 
329  if(dNumHitsUsedFordEdx_FDC >= locMinimumNumHitsUsed){ //6+ original hits
330  dSelectorHist_dEdxVsBeta_HitsCutoff_FDC->Fill(dBeta, ddEdxkeV_FDC);
331  dSelectorHist_dEdxVsP_HitsCutoff_FDC->Fill(dMomentum, ddEdxkeV_FDC);
332  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC->Fill(locBetaGamma, ddEdxkeV_FDC);
333  }
334 
335  if(dNumHitsUsedFordEdx_CDC >= locMinimumNumHitsUsed){ //6+ original hits
336  dSelectorHist_dEdxVsBeta_HitsCutoff_CDC->Fill(dBeta, ddEdxkeV_CDC);
337  dSelectorHist_dEdxVsP_HitsCutoff_CDC->Fill(dMomentum, ddEdxkeV_CDC);
338  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC->Fill(locBetaGamma, ddEdxkeV_CDC);
339  }
340 
341  if(dNumHitsUsedFordEdx_FDC == 3)
342  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC_3Hits->Fill(locBetaGamma, ddEdxkeV_FDC);
343  if(dNumHitsUsedFordEdx_FDC == 6)
344  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC_6Hits->Fill(locBetaGamma, ddEdxkeV_FDC);
345  if(dNumHitsUsedFordEdx_FDC == 9)
346  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC_9Hits->Fill(locBetaGamma, ddEdxkeV_FDC);
347  if(dNumHitsUsedFordEdx_FDC == 12)
348  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC_12Hits->Fill(locBetaGamma, ddEdxkeV_FDC);
349 
350  if(dNumHitsUsedFordEdx_CDC == 2)
351  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_2Hits->Fill(locBetaGamma, ddEdxkeV_CDC);
352  if(dNumHitsUsedFordEdx_CDC == 4)
353  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_4Hits->Fill(locBetaGamma, ddEdxkeV_CDC);
354  if(dNumHitsUsedFordEdx_CDC == 6)
355  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_6Hits->Fill(locBetaGamma, ddEdxkeV_CDC);
356  if(dNumHitsUsedFordEdx_CDC == 8)
357  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_8Hits->Fill(locBetaGamma, ddEdxkeV_CDC);
358  if(dNumHitsUsedFordEdx_CDC == 10)
359  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_10Hits->Fill(locBetaGamma, ddEdxkeV_CDC);
360  if(dNumHitsUsedFordEdx_CDC == 12)
361  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_12Hits->Fill(locBetaGamma, ddEdxkeV_CDC);
362  if(dNumHitsUsedFordEdx_CDC == 14)
363  dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_14Hits->Fill(locBetaGamma, ddEdxkeV_CDC);
364 
365  double locFOM;
366  bool locCalcStatus = Calc_FOM(locFOM);
367 /*
368  if(fabs(dFOM - locFOM) > 0.001){
369  cout << "orig fom, new fom, diff = " << dFOM << ", " << locFOM << ", " << (dFOM - locFOM) << endl;
370  cout << "#hits CDC, #hits FDC, dE/dx CDC, dE/dx FDC = " << dNumHitsUsedFordEdx_CDC << ", " << dNumHitsUsedFordEdx_FDC << ", " << ddEdx_CDC << ", " << ddEdx_FDC << endl;
371  }
372 */
373  if(dCalcFOMManuallyFlag == true){
374  if(locCalcStatus == false)
375  return kTRUE;
376  dFOM = locFOM;
377  }
378 
379 
380  if((dNumHitsUsedFordEdx_FDC >= 3) && (dNumHitsUsedFordEdx_CDC < 3)){ //FDC only
382  if(dNumHitsUsedFordEdx_FDC == 3)
384  if(dNumHitsUsedFordEdx_FDC == 6)
386  if(dNumHitsUsedFordEdx_FDC == 9)
388  if(dNumHitsUsedFordEdx_FDC == 12)
390 
391  if((locBetaGamma >= 0.3) && (locBetaGamma < 0.6))
393  if((locBetaGamma >= 0.6) && (locBetaGamma < 0.9))
395  if((locBetaGamma >= 0.9) && (locBetaGamma < 1.2))
397  if((locBetaGamma >= 1.2) && (locBetaGamma < 1.5))
399  if((locBetaGamma >= 1.5) && (locBetaGamma < 1.8))
401  if((locBetaGamma >= 1.8) && (locBetaGamma < 2.1))
403  if((locBetaGamma >= 2.1) && (locBetaGamma < 4.1))
405  if((locBetaGamma >= 4.1) && (locBetaGamma < 6.1))
407  if((locBetaGamma >= 6.1) && (locBetaGamma < 8.1))
409  if((locBetaGamma >= 8.1) && (locBetaGamma < 10.1))
411  if((locBetaGamma >= 10.1) && (locBetaGamma < 12.1))
413  if((locBetaGamma >= 12.1) && (locBetaGamma < 14.1))
415  }
416  if((dNumHitsUsedFordEdx_CDC >= 3) && (dNumHitsUsedFordEdx_FDC < 3)){ //CDC only
418  if(dNumHitsUsedFordEdx_CDC == 4)
420  if(dNumHitsUsedFordEdx_CDC == 6)
422  if(dNumHitsUsedFordEdx_CDC == 8)
424  if(dNumHitsUsedFordEdx_CDC == 10)
426  if(dNumHitsUsedFordEdx_CDC == 12)
428  if(dNumHitsUsedFordEdx_CDC == 14)
430 
431  if((locBetaGamma >= 0.3) && (locBetaGamma < 0.6))
433  if((locBetaGamma >= 0.6) && (locBetaGamma < 0.9))
435  if((locBetaGamma >= 0.9) && (locBetaGamma < 1.2))
437  if((locBetaGamma >= 1.2) && (locBetaGamma < 1.5))
439  if((locBetaGamma >= 1.5) && (locBetaGamma < 1.8))
441  if((locBetaGamma >= 1.8) && (locBetaGamma < 2.1))
443  if((locBetaGamma >= 2.1) && (locBetaGamma < 4.1))
445  if((locBetaGamma >= 4.1) && (locBetaGamma < 6.1))
447  if((locBetaGamma >= 6.1) && (locBetaGamma < 8.1))
449  if((locBetaGamma >= 8.1) && (locBetaGamma < 10.1))
451  if((locBetaGamma >= 10.1) && (locBetaGamma < 12.1))
453  if((locBetaGamma >= 12.1) && (locBetaGamma < 14.1))
455 
456 
457  }
458 
459  if((dNumHitsUsedFordEdx_CDC >= 3) && (dNumHitsUsedFordEdx_FDC >= 3)){ //Both
461 
462  double locMeandEdx_CDC, locMeandEdx_FDC;
463  if(!GetdEdxMean_CDC(dBeta, dNumHitsUsedFordEdx_CDC, locMeandEdx_CDC))
464  return kTRUE;
465  if(!GetdEdxMean_FDC(dBeta, dNumHitsUsedFordEdx_FDC, locMeandEdx_FDC))
466  return kTRUE;
467 
468  //Correlation coefficient:
469  //http://en.wikipedia.org/wiki/Pearson_product-moment_correlation_coefficient#Geometric_interpretation
470 
471 
472  }
473 
474  return kTRUE;
475 }
476 
478 {
479  // The SlaveTerminate() function is called after all entries or objects
480  // have been processed. When running with PROOF SlaveTerminate() is called
481  // on each slave server.
482 
483 }
484 
486 {
487  // The Terminate() function is the last function to be called during
488  // a query. It always runs on the client, it can be used to present
489  // the results graphically or save the results to file.
490 
491  float locNumNotEnoughHits, locNumGenerated, locNotEnoughHitsRatio, locNotEnoughHitsRatioUncertainty;
492  for(unsigned int loc_i = 1; loc_i <= dSelectorHist_ThetaDistribution->GetNbinsX(); loc_i++){
493  locNumNotEnoughHits = float(dSelectorHist_NotEnoughHitsInTheta->GetBinContent(loc_i));
494  locNumGenerated = float(dSelectorHist_ThetaDistribution->GetBinContent(loc_i));
495  locNotEnoughHitsRatio = locNumNotEnoughHits/locNumGenerated;
496  locNotEnoughHitsRatioUncertainty = sqrt(locNumNotEnoughHits*(1.0 - locNotEnoughHitsRatio))/locNumGenerated;
497  dSelectorHist_PercentageNotEnoughHitsInTheta->SetBinContent(loc_i, locNotEnoughHitsRatio);
498  dSelectorHist_PercentageNotEnoughHitsInTheta->SetBinError(loc_i, locNotEnoughHitsRatioUncertainty);
499  }
500 
501  dOutputFile->Write();
502 }
503 
504 bool DCdEdxSelector::GetdEdxMean_CDC(double locBeta, unsigned int locNumHitsUsedFordEdx, double& locMeandEdx){
505  double locBetaGammaValue = locBeta/sqrt(1.0 - locBeta*locBeta);
506  if(dParticleName == "Proton"){
507  if(locBetaGammaValue < 0.3)
508  return false;
509  if(locBetaGammaValue > 2.15)
510  return false;
511 
512  locMeandEdx = (ddEdxMeanFunc_CDC_Proton->Eval(locBetaGammaValue))/1000000.0;
513  return true;
514  }
515  if(dParticleName == "KPlus"){
516  if(locBetaGammaValue < 0.65)
517  return false; //K+ very likely to decay, don't compute chisq!!
518  if(locBetaGammaValue > 4.1)
519  return false;
520 
521  locMeandEdx = (ddEdxMeanFunc_CDC_KPlus->Eval(locBetaGammaValue))/1000000.0;
522  return true;
523  }
524 
525  double locBetaGamma[44];
526  double locdEdxMean[44];
527  locBetaGamma[0] = 0.48; locBetaGamma[1] = 0.8; locBetaGamma[2] = 1.12; locBetaGamma[3] = 1.44; locBetaGamma[4] = 1.76; locBetaGamma[5] = 2.08; locBetaGamma[6] = 2.4; locBetaGamma[7] = 2.72; locBetaGamma[8] = 3.04; locBetaGamma[9] = 3.36; locBetaGamma[10] = 3.68; locBetaGamma[11] = 4; locBetaGamma[12] = 4.32; locBetaGamma[13] = 4.64; locBetaGamma[14] = 4.96; locBetaGamma[15] = 5.28; locBetaGamma[16] = 5.6; locBetaGamma[17] = 5.92; locBetaGamma[18] = 6.24; locBetaGamma[19] = 6.56; locBetaGamma[20] = 6.88; locBetaGamma[21] = 7.2; locBetaGamma[22] = 7.52; locBetaGamma[23] = 7.84; locBetaGamma[24] = 8.16; locBetaGamma[25] = 8.48; locBetaGamma[26] = 8.8; locBetaGamma[27] = 9.12; locBetaGamma[28] = 9.44; locBetaGamma[29] = 9.76; locBetaGamma[30] = 10.08; locBetaGamma[31] = 10.4; locBetaGamma[32] = 10.72; locBetaGamma[33] = 11.04; locBetaGamma[34] = 11.36; locBetaGamma[35] = 11.68; locBetaGamma[36] = 12; locBetaGamma[37] = 12.32; locBetaGamma[38] = 12.64; locBetaGamma[39] = 12.96; locBetaGamma[40] = 13.28; locBetaGamma[41] = 13.6; locBetaGamma[42] = 13.92; locBetaGamma[43] = 14.24;
528  locdEdxMean[0] = 2.40705; locdEdxMean[1] = 2.64016; locdEdxMean[2] = 1.97719; locdEdxMean[3] = 1.74434; locdEdxMean[4] = 1.60673; locdEdxMean[5] = 1.52877; locdEdxMean[6] = 1.47792; locdEdxMean[7] = 1.44914; locdEdxMean[8] = 1.44456; locdEdxMean[9] = 1.44089; locdEdxMean[10] = 1.44077; locdEdxMean[11] = 1.44643; locdEdxMean[12] = 1.46221; locdEdxMean[13] = 1.47146; locdEdxMean[14] = 1.48451; locdEdxMean[15] = 1.49671; locdEdxMean[16] = 1.50718; locdEdxMean[17] = 1.51892; locdEdxMean[18] = 1.53035; locdEdxMean[19] = 1.54224; locdEdxMean[20] = 1.5531; locdEdxMean[21] = 1.56325; locdEdxMean[22] = 1.57347; locdEdxMean[23] = 1.58422; locdEdxMean[24] = 1.59239; locdEdxMean[25] = 1.60057; locdEdxMean[26] = 1.61224; locdEdxMean[27] = 1.62025; locdEdxMean[28] = 1.6288; locdEdxMean[29] = 1.63748; locdEdxMean[30] = 1.64582; locdEdxMean[31] = 1.65254; locdEdxMean[32] = 1.65864; locdEdxMean[33] = 1.66693; locdEdxMean[34] = 1.67503; locdEdxMean[35] = 1.68328; locdEdxMean[36] = 1.68964; locdEdxMean[37] = 1.69732; locdEdxMean[38] = 1.70258; locdEdxMean[39] = 1.70849; locdEdxMean[40] = 1.71469; locdEdxMean[41] = 1.72359; locdEdxMean[42] = 1.73032; locdEdxMean[43] = 1.73672;
529 
530  if(locBetaGammaValue < 0.81)
531  return false;
532  if(locBetaGammaValue > 14.2)
533  return false;
534 
535  double locSlope, locIntercept;
536  for(unsigned int loc_i = 0; loc_i < 43; loc_i++){
537  if(locBetaGammaValue > locBetaGamma[loc_i + 1])
538  continue;
539  locSlope = (locdEdxMean[loc_i + 1] - locdEdxMean[loc_i])/(locBetaGamma[loc_i + 1] - locBetaGamma[loc_i]);
540  locIntercept = locdEdxMean[loc_i] - locSlope*locBetaGamma[loc_i]; //y = mx + b, b = y - mx
541  locMeandEdx = (locSlope*locBetaGammaValue + locIntercept)/1000000.0;
542  return true;
543  }
544 
545  return false;
546 }
547 
548 bool DCdEdxSelector::GetdEdxMean_FDC(double locBeta, unsigned int locNumHitsUsedFordEdx, double& locMeandEdx){
549  double locBetaGammaValue = locBeta/sqrt(1.0 - locBeta*locBeta);
550 
551  if(dParticleName == "Proton"){
552  if(locBetaGammaValue < 0.3)
553  return false;
554  if(locBetaGammaValue > 2.15)
555  return false;
556 
557  locMeandEdx = (ddEdxMeanFunc_FDC_Proton->Eval(locBetaGammaValue))/1000000.0;
558  return true;
559  }
560  if(dParticleName == "KPlus"){
561  if(locBetaGammaValue < 0.9)
562  return false; //K+ very likely to decay, don't compute chisq!!
563  if(locBetaGammaValue > 4.1)
564  return false;
565 
566  locMeandEdx = (ddEdxMeanFunc_FDC_KPlus->Eval(locBetaGammaValue))/1000000.0;
567  return true;
568  }
569 
570  double locBetaGamma[44];
571  double locdEdxMean[44];
572  locBetaGamma[0] = 0.48; locBetaGamma[1] = 0.8; locBetaGamma[2] = 1.12; locBetaGamma[3] = 1.44; locBetaGamma[4] = 1.76; locBetaGamma[5] = 2.08; locBetaGamma[6] = 2.4; locBetaGamma[7] = 2.72; locBetaGamma[8] = 3.04; locBetaGamma[9] = 3.36; locBetaGamma[10] = 3.68; locBetaGamma[11] = 4; locBetaGamma[12] = 4.32; locBetaGamma[13] = 4.64; locBetaGamma[14] = 4.96; locBetaGamma[15] = 5.28; locBetaGamma[16] = 5.6; locBetaGamma[17] = 5.92; locBetaGamma[18] = 6.24; locBetaGamma[19] = 6.56; locBetaGamma[20] = 6.88; locBetaGamma[21] = 7.2; locBetaGamma[22] = 7.52; locBetaGamma[23] = 7.84; locBetaGamma[24] = 8.16; locBetaGamma[25] = 8.48; locBetaGamma[26] = 8.8; locBetaGamma[27] = 9.12; locBetaGamma[28] = 9.44; locBetaGamma[29] = 9.76; locBetaGamma[30] = 10.08; locBetaGamma[31] = 10.4; locBetaGamma[32] = 10.72; locBetaGamma[33] = 11.04; locBetaGamma[34] = 11.36; locBetaGamma[35] = 11.68; locBetaGamma[36] = 12; locBetaGamma[37] = 12.32; locBetaGamma[38] = 12.64; locBetaGamma[39] = 12.96; locBetaGamma[40] = 13.28; locBetaGamma[41] = 13.6; locBetaGamma[42] = 13.92; locBetaGamma[43] = 14.24;
573  locdEdxMean[0] = 1.32105; locdEdxMean[1] = 1.50904; locdEdxMean[2] = 1.50873; locdEdxMean[3] = 1.57426; locdEdxMean[4] = 1.52213; locdEdxMean[5] = 1.45914; locdEdxMean[6] = 1.4271; locdEdxMean[7] = 1.42623; locdEdxMean[8] = 1.42268; locdEdxMean[9] = 1.42632; locdEdxMean[10] = 1.43233; locdEdxMean[11] = 1.43742; locdEdxMean[12] = 1.4467; locdEdxMean[13] = 1.45642; locdEdxMean[14] = 1.46337; locdEdxMean[15] = 1.4719; locdEdxMean[16] = 1.48252; locdEdxMean[17] = 1.49443; locdEdxMean[18] = 1.50234; locdEdxMean[19] = 1.513; locdEdxMean[20] = 1.52394; locdEdxMean[21] = 1.5289; locdEdxMean[22] = 1.53878; locdEdxMean[23] = 1.54396; locdEdxMean[24] = 1.55427; locdEdxMean[25] = 1.56105; locdEdxMean[26] = 1.5674; locdEdxMean[27] = 1.57221; locdEdxMean[28] = 1.58027; locdEdxMean[29] = 1.5865; locdEdxMean[30] = 1.59519; locdEdxMean[31] = 1.60211; locdEdxMean[32] = 1.60751; locdEdxMean[33] = 1.61522; locdEdxMean[34] = 1.62066; locdEdxMean[35] = 1.6277; locdEdxMean[36] = 1.63227; locdEdxMean[37] = 1.63661; locdEdxMean[38] = 1.64457; locdEdxMean[39] = 1.65043; locdEdxMean[40] = 1.6557; locdEdxMean[41] = 1.66085; locdEdxMean[42] = 1.66486; locdEdxMean[43] = 1.66916;
574 
575  if(locBetaGammaValue < 1.77)
576  return false;
577  if(locBetaGammaValue > 14.2)
578  return false;
579 
580  double locSlope, locIntercept;
581  for(unsigned int loc_i = 0; loc_i < 43; loc_i++){
582  if(locBetaGammaValue > locBetaGamma[loc_i + 1])
583  continue;
584  locSlope = (locdEdxMean[loc_i + 1] - locdEdxMean[loc_i])/(locBetaGamma[loc_i + 1] - locBetaGamma[loc_i]);
585  locIntercept = locdEdxMean[loc_i] - locSlope*locBetaGamma[loc_i]; //y = mx + b, b = y - mx
586  locMeandEdx = (locSlope*locBetaGammaValue + locIntercept)/1000000.0;
587  return true;
588  }
589  return false;
590 }
591 
592 
593 bool DCdEdxSelector::GetdEdxSigma_FDC(double locBeta, unsigned int locNumHitsUsedFordEdx, double& locSigmadEdx){
594  double locBetaGamma = locBeta/sqrt(1.0 - locBeta*locBeta);
595  double locSigmadEdx_LowSide, locSigmadEdx_HighSide; //for linear interpolation/extrapolation
596 
597  vector<int> locNumHitsVector;
598  TObjArray *locFuncArray;
599  if(dParticleName == "Proton"){
600  locNumHitsVector = ddEdxSigmaNumHitsVector_FDC_Proton;
601  locFuncArray = dSigmaFuncArray_FDC_Proton;
602  }
603  if(dParticleName == "KPlus"){
604  locNumHitsVector = ddEdxSigmaNumHitsVector_FDC_KPlus;
605  locFuncArray = dSigmaFuncArray_FDC_KPlus;
606  }
607  if(dParticleName == "PiPlus"){
608  locNumHitsVector = ddEdxSigmaNumHitsVector_FDC_PiPlus;
609  locFuncArray = dSigmaFuncArray_FDC_PiPlus;
610  }
611 
612  double locSlope, locIntercept;
613 
614  for(unsigned int loc_i = 0; loc_i < locNumHitsVector.size(); loc_i++){
615  if(locNumHitsUsedFordEdx == locNumHitsVector[loc_i]){
616  locSigmadEdx = (((TF1*)(*locFuncArray)[loc_i])->Eval(locBetaGamma))/1000000.0;
617  return true;
618  }
619  if(loc_i == (locNumHitsVector.size() - 1)){
620 //extrapolate!!
621  locSigmadEdx_LowSide = ((TF1*)(*locFuncArray)[loc_i - 1])->Eval(locBetaGamma);
622  locSigmadEdx_HighSide = ((TF1*)(*locFuncArray)[loc_i])->Eval(locBetaGamma);
623 
624  locSlope = (locSigmadEdx_HighSide - locSigmadEdx_LowSide)/(locNumHitsVector[loc_i] - locNumHitsVector[loc_i - 1]);
625  locIntercept = locSigmadEdx_HighSide - locSlope*locNumHitsVector[loc_i]; //y = mx + b, b = y - mx
626  locSigmadEdx = (locSlope*locNumHitsUsedFordEdx + locIntercept)/1000000.0;
627 
628  return true;
629  }
630  if((locNumHitsUsedFordEdx > locNumHitsVector[loc_i]) && (locNumHitsUsedFordEdx < locNumHitsVector[loc_i + 1])){
631  locSigmadEdx_LowSide = ((TF1*)(*locFuncArray)[loc_i])->Eval(locBetaGamma);
632  locSigmadEdx_HighSide = ((TF1*)(*locFuncArray)[loc_i + 1])->Eval(locBetaGamma);
633 
634  locSlope = (locSigmadEdx_HighSide - locSigmadEdx_LowSide)/(locNumHitsVector[loc_i + 1] - locNumHitsVector[loc_i]);
635  locIntercept = locSigmadEdx_HighSide - locSlope*locNumHitsVector[loc_i + 1]; //y = mx + b, b = y - mx
636  locSigmadEdx = (locSlope*locNumHitsUsedFordEdx + locIntercept)/1000000.0;
637 
638  return true;
639  }
640  }
641  return true;
642 }
643 
644 bool DCdEdxSelector::GetdEdxSigma_CDC(double locBeta, unsigned int locNumHitsUsedFordEdx, double& locSigmadEdx){
645  double locBetaGamma = locBeta/sqrt(1.0 - locBeta*locBeta);
646  double locSigmadEdx_LowSide, locSigmadEdx_HighSide; //for linear interpolation/extrapolation
647 
648  vector<int> locNumHitsVector;
649  TObjArray *locFuncArray;
650  if(dParticleName == "Proton"){
651  locNumHitsVector = ddEdxSigmaNumHitsVector_CDC_Proton;
652  locFuncArray = dSigmaFuncArray_CDC_Proton;
653  }
654  if(dParticleName == "KPlus"){
655  locNumHitsVector = ddEdxSigmaNumHitsVector_CDC_KPlus;
656  locFuncArray = dSigmaFuncArray_CDC_KPlus;
657  }
658  if(dParticleName == "PiPlus"){
659  locNumHitsVector = ddEdxSigmaNumHitsVector_CDC_PiPlus;
660  locFuncArray = dSigmaFuncArray_CDC_PiPlus;
661  }
662 
663  double locSlope, locIntercept;
664 
665  for(unsigned int loc_i = 0; loc_i < locNumHitsVector.size(); loc_i++){
666  if(locNumHitsUsedFordEdx == locNumHitsVector[loc_i]){
667  locSigmadEdx = (((TF1*)(*locFuncArray)[loc_i])->Eval(locBetaGamma))/1000000.0;
668  return true;
669  }
670  if(loc_i == (locNumHitsVector.size() - 1)){
671 //extrapolate!!
672  locSigmadEdx_LowSide = ((TF1*)(*locFuncArray)[loc_i - 1])->Eval(locBetaGamma);
673  locSigmadEdx_HighSide = ((TF1*)(*locFuncArray)[loc_i])->Eval(locBetaGamma);
674 
675  locSlope = (locSigmadEdx_HighSide - locSigmadEdx_LowSide)/(locNumHitsVector[loc_i] - locNumHitsVector[loc_i - 1]);
676  locIntercept = locSigmadEdx_HighSide - locSlope*locNumHitsVector[loc_i]; //y = mx + b, b = y - mx
677  locSigmadEdx = (locSlope*locNumHitsUsedFordEdx + locIntercept)/1000000.0;
678 
679  return true;
680  }
681  if((locNumHitsUsedFordEdx > locNumHitsVector[loc_i]) && (locNumHitsUsedFordEdx < locNumHitsVector[loc_i + 1])){
682  locSigmadEdx_LowSide = ((TF1*)(*locFuncArray)[loc_i])->Eval(locBetaGamma);
683  locSigmadEdx_HighSide = ((TF1*)(*locFuncArray)[loc_i + 1])->Eval(locBetaGamma);
684 
685  locSlope = (locSigmadEdx_HighSide - locSigmadEdx_LowSide)/(locNumHitsVector[loc_i + 1] - locNumHitsVector[loc_i]);
686  locIntercept = locSigmadEdx_HighSide - locSlope*locNumHitsVector[loc_i + 1]; //y = mx + b, b = y - mx
687  locSigmadEdx = (locSlope*locNumHitsUsedFordEdx + locIntercept)/1000000.0;
688 
689  return true;
690  }
691  }
692 
693  return true;
694 }
695 
696 
697 bool DCdEdxSelector::Calc_FOM(double& locFOM){
698  double locDCdEdx, locChiSq;
699  unsigned int locNDF;
700  unsigned int locMinimumNumberUsedHitsForConfidence = 3; //dE/dx is landau-distributed, so to approximate Gaussian must remove hits with largest dE/dx //3 means 6 or more hits originally
701 
702  bool locUseCDCHitsFlag = (dNumHitsUsedFordEdx_CDC > locMinimumNumberUsedHitsForConfidence) ? true : false;
703  bool locUseFDCHitsFlag = (dNumHitsUsedFordEdx_FDC > locMinimumNumberUsedHitsForConfidence) ? true : false;
704 
705 //cout << "#hits cdc, fdc = " << locNumHitsUsedFordEdx_CDC << ", " << locNumHitsUsedFordEdx_FDC << endl;
706  if((locUseCDCHitsFlag == false) && (locUseFDCHitsFlag == false))
707  return false; //not enough hits, use other sources of information for PID
708 
709  double locMeandEdx_FDC, locMeandEdx_CDC, locSigmadEdx_FDC, locSigmadEdx_CDC;
710  double locDeltadEdx_CDC = 0.0, locDeltadEdx_FDC = 0.0;
711 
712  if(locUseCDCHitsFlag == true){
713  if(GetdEdxMean_CDC(dBeta, dNumHitsUsedFordEdx_CDC, locMeandEdx_CDC) != true)
714  return false;
715  if(GetdEdxSigma_CDC(dBeta, dNumHitsUsedFordEdx_CDC, locSigmadEdx_CDC) != true)
716  return false;
717  locDeltadEdx_CDC = ddEdx_CDC - locMeandEdx_CDC;
718 //cout << "CDC: #hits, beta, actual, mean, sigma = " << dNumHitsUsedFordEdx_CDC << ", " << dBeta << ", " << ddEdx_CDC << ", " << locMeandEdx_CDC << ", " << locSigmadEdx_CDC << endl;
719  }
720 
721  if(locUseFDCHitsFlag == true){
722  if(GetdEdxMean_FDC(dBeta, dNumHitsUsedFordEdx_FDC, locMeandEdx_FDC) != true)
723  return false;
724  if(GetdEdxSigma_FDC(dBeta, dNumHitsUsedFordEdx_FDC, locSigmadEdx_FDC) != true)
725  return false;
726  locDeltadEdx_FDC = ddEdx_FDC - locMeandEdx_FDC;
727 //cout << "FDC: #hits, beta, actual, mean, sigma = " << dNumHitsUsedFordEdx_FDC << ", " << dBeta << ", " << ddEdx_FDC << ", " << locMeandEdx_FDC << ", " << locSigmadEdx_FDC << endl;
728  }
729 
730  if((locUseCDCHitsFlag == true) && ((locUseFDCHitsFlag == false) || (dNumHitsUsedFordEdx_CDC >= dNumHitsUsedFordEdx_FDC))){
731  locDCdEdx = ddEdx_CDC;
732  locChiSq = locDeltadEdx_CDC/locSigmadEdx_CDC;
733  locChiSq *= locChiSq;
734  locNDF = 1;
735  locFOM = TMath::Prob(locChiSq, locNDF);
736  return true;
737  }
738 
739  if((locUseFDCHitsFlag == true) && ((locUseCDCHitsFlag == false) || (dNumHitsUsedFordEdx_FDC > dNumHitsUsedFordEdx_CDC)) ){
740  locDCdEdx = ddEdx_FDC;
741  locChiSq = locDeltadEdx_FDC/locSigmadEdx_FDC;
742  locChiSq *= locChiSq;
743  locNDF = 1;
744  locFOM = TMath::Prob(locChiSq, locNDF);
745  return true;
746  }
747  return true;
748 }
749 
vector< int > ddEdxSigmaNumHitsVector_CDC_PiPlus
TH1F * dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin8
TH1F * dSelectorHist_ConfidenceLevel_Both
TH1F * dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin3
TH2F * dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_12Hits
TH2F * dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_4Hits
TH2F * dSelectorHist_NumHitsVsTheta_CDC
TH1F * dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin1
vector< int > ddEdxSigmaNumHitsVector_CDC_KPlus
TH2F * dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC
bool GetdEdxSigma_CDC(double locBeta, unsigned int locNumHitsUsedFordEdx, double &locSigmadEdx)
TH2F * dSelectorHist_dEdxVsP_HitsCutoff_FDC
TH2F * dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC
TH1F * dSelectorHist_PercentageNotEnoughHitsInTheta
bool GetdEdxMean_FDC(double locBeta, unsigned int locNumHitsUsedFordEdx, double &locMeandEdx)
TH1F * dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin2
TH2F * dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC_12Hits
bool GetdEdxMean_CDC(double locBeta, unsigned int locNumHitsUsedFordEdx, double &locMeandEdx)
TH1F * dSelectorHist_ConfidenceLevel_CDC_14Hits
bool Calc_FOM(double &locFOM)
TF1 * ddEdxMeanFunc_FDC_KPlus
TH1F * dSelectorHist_ConfidenceLevel_CDC_4Hits
bool dCalcFOMManuallyFlag
TH2F * dSelectorHist_dEdxVsP_HitsCutoff_CDC
TF1 * ddEdxMeanFunc_FDC_Proton
TH1F * dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin1
TFile * dOutputFile
UInt_t dNumHitsUsedFordEdx_FDC
TH2F * dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_14Hits
TH1F * dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin5
virtual Int_t GetEntry(Long64_t entry, Int_t getall=0)
TH1F * dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin11
TH2F * dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC_3Hits
TH1F * dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin2
virtual Bool_t Process(Long64_t entry)
virtual void SlaveBegin(TTree *tree)
Double_t ddEdx_FDC
TH1F * dSelectorHist_ConfidenceLevel_FDC_6Hits
TObjArray * dSigmaFuncArray_FDC_PiPlus
TObjArray * dSigmaFuncArray_CDC_Proton
Double_t dMomentum
TF1 * ddEdxMeanFunc_CDC_KPlus
TH1F * dSelectorHist_ConfidenceLevel_CDC_10Hits
vector< int > ddEdxSigmaNumHitsVector_FDC_PiPlus
string dParticleName
TH2F * dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_8Hits
Double_t ddEdx_CDC
UInt_t dNumHitsUsedFordEdx_CDC
Double_t ddx_CDC
vector< int > ddEdxSigmaNumHitsVector_FDC_KPlus
TH1F * dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin7
bool GetdEdxSigma_FDC(double locBeta, unsigned int locNumHitsUsedFordEdx, double &locSigmadEdx)
TH2F * dSelectorHist_dxVsNumHits_FDC
TH1F * dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin12
TH1F * dSelectorHist_ConfidenceLevel_CDC_6Hits
TH2F * dSelectorHist_dxVsNumHits_CDC
vector< int > ddEdxSigmaNumHitsVector_CDC_Proton
double dRhoZoverA_FDC
TH2F * dSelectorHist_dxVsTheta_FDC
TH2F * dSelectorHist_NumHitsVsTheta_FDC
TH1F * dSelectorHist_ConfidenceLevel_CDC_12Hits
TH1F * dSelectorHist_ConfidenceLevel_FDC_3Hits
TH1F * dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin6
TH1F * dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin4
TH1F * dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin7
double sqrt(double)
TH1F * dSelectorHist_ThetaDistribution
double dRhoZoverA_CDC
TH1F * dSelectorHist_ConfidenceLevel_FDC
TH2F * dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC_9Hits
TH2F * dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_6Hits
TH1F * dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin9
TH1F * dSelectorHist_ConfidenceLevel_FDC_9Hits
Double_t ddx_FDC
TObjArray * dSigmaFuncArray_CDC_PiPlus
TH1F * dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin9
TH1F * dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin12
TH2F * dSelectorHist_dxVsTheta_CDC
TH1F * dSelectorHist_ConfidenceLevel_CDC
TH2F * dSelectorHist_dEdxVsBeta_HitsCutoff_FDC
TH1F * dSelectorHist_ConfidenceLevel_FDC_12Hits
TObjArray * dSigmaFuncArray_CDC_KPlus
TH2F * dSelectorHist_dEdxVsBetaGamma_HitsCutoff_FDC_6Hits
TH1F * dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin3
TH1F * dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin10
virtual void Terminate()
TH1F * dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin5
TH1F * dSelectorHist_ConfidenceLevel_CDC_8Hits
TH2F * dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_2Hits
TH1F * dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin8
TH2F * dSelectorHist_dEdxVsBetaGamma_HitsCutoff_CDC_10Hits
TH2F * dSelectorHist_dEdxVsBeta_HitsCutoff_CDC
TH1F * dSelectorHist_NotEnoughHitsInTheta
virtual void Begin(TTree *tree)
TH2F * dSelectorHist_NumHitsFDCVsNumHitsCDC
TF1 * ddEdxMeanFunc_CDC_Proton
TH1F * dSelectorHist_ConfidenceLevel_FDC_BetaGammaBin10
vector< int > ddEdxSigmaNumHitsVector_FDC_Proton
TObjArray * dSigmaFuncArray_FDC_Proton
TH1F * dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin4
TH1F * dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin11
virtual void SlaveTerminate()
TH1F * dSelectorHist_ConfidenceLevel_CDC_BetaGammaBin6
TObjArray * dSigmaFuncArray_FDC_KPlus