File indexing completed on 2025-12-16 09:18:03
0001
0002
0003
0004
0005
0006
0007
0008 #ifndef InttSeedTracking_h
0009 #define InttSeedTracking_h
0010
0011 #include <TROOT.h>
0012 #include <TChain.h>
0013 #include <TFile.h>
0014
0015 #include <TMath.h>
0016 #include "Fit/Fitter.h"
0017 #include <Math/Functor.h>
0018
0019 #include <TH1.h>
0020 #include <TH2.h>
0021 #include <TH3.h>
0022 #include <TF1.h>
0023
0024 #include "TArc.h"
0025
0026 #include "SPHTracKuma.h"
0027
0028
0029
0030 #include <vector>
0031
0032 class InttSeedTracking {
0033 public :
0034 InttSeedTracking(std::vector<tracKuma>& tracks,\
0035 std::vector<hitStruct > vFMvtxHits,\
0036 std::vector<hitStruct > vSMvtxHits, std::vector<hitStruct > vTMvtxHits,\
0037 std::vector<hitStruct > vIInttHits, std::vector<hitStruct > vOInttHits,\
0038 std::vector<hitStruct > vEmcalHits,\
0039 std::vector<hitStruct > vIHCalHits, std::vector<hitStruct > vOHcalHits);
0040 InttSeedTracking(){};
0041 virtual ~InttSeedTracking();
0042
0043 void HitMatching(std::vector<tracKuma>& tracks, std::vector<hitStruct > vFMvtxHits,\
0044 std::vector<hitStruct > vSMvtxHits, std::vector<hitStruct > vTMvtxHits,\
0045 std::vector<hitStruct > vIInttHits, std::vector<hitStruct > vOInttHits,\
0046 std::vector<hitStruct > vEmcalHits,\
0047 std::vector<hitStruct > vIHCalHits, std::vector<hitStruct > vOHcalHits);
0048
0049 Double_t AddHCalE(Double_t emcalPhi, Double_t emcalE,\
0050 std::vector<hitStruct > vIHCalHits, std::vector<hitStruct > vOHcalHits);
0051
0052 void AddMvtxHits(tracKuma& trk, std::vector<hitStruct > vFMvtxHits, \
0053 std::vector<hitStruct > vSMvtxHits, std::vector<hitStruct > vTMvtxHits, Double_t dRThre, Int_t& seleFMvtxId, Int_t& seleSMvtxId, Int_t& seleTMvtxId);
0054
0055
0056 void TrackPropertiesEstimation(tracKuma& trk, std::vector<hitStruct > vFMvtxHits,\
0057 std::vector<hitStruct > vSMvtxHits, std::vector<hitStruct > vTMvtxHits);
0058 Double_t TrackPtEstimation(tracKuma& trk, std::vector<hitStruct > vFMvtxHits,\
0059 std::vector<hitStruct > vSMvtxHits, std::vector<hitStruct > vTMvtxHits);
0060
0061
0062 void RoughEstiSagittaCenter3Point(Double_t& sagittaR, \
0063 Double_t& sagittaCenterX, Double_t& sagittaCenterY, Double_t HitsXY[3][2]);
0064
0065 Double_t AccuratePtEstimation(Double_t sagittaR, Double_t centerX, Double_t centerY, tracKuma& trk, std::vector<hitStruct > vFMvtxHits,\
0066 std::vector<hitStruct > vSMvtxHits, std::vector<hitStruct > vTMvtxHits);
0067
0068 bool RecoTracksInttSeed2(std::vector<tracKuma>& tracks, std::vector<hitStruct > vFMvtxHits,\
0069 std::vector<hitStruct > vSMvtxHits, std::vector<hitStruct > vTMvtxHits,\
0070 std::vector<hitStruct > vIInttHits, std::vector<hitStruct > vOInttHits,\
0071 std::vector<hitStruct > vEmcalHits,\
0072 std::vector<hitStruct > vIHCalHits, std::vector<hitStruct > vOHcalHits);
0073 bool InttSeedMatching(std::vector<tracKuma>& tracks, Int_t inttId,\
0074 hitStruct baseInttHit,\
0075 std::vector<hitStruct >& vFMvtxHits,\
0076 std::vector<hitStruct >& vSMvtxHits, std::vector<hitStruct >& vTMvtxHits,\
0077 std::vector<hitStruct >& vIInttHits,\
0078 std::vector<hitStruct > vEmcalHits,\
0079 std::vector<hitStruct > vIHCalHits, std::vector<hitStruct > vOHCalHits);
0080 Double_t EstiChiTrkOnCircle(tracKuma trk, Double_t cX, Double_t cY, Double_t sagittaR);
0081 bool CheckTrkRequirements(tracKuma trk);
0082 void RefindCalHit(tracKuma trk, std::vector<hitStruct > vEmcalHits,\
0083 std::vector<hitStruct > vIHCalHits, std::vector<hitStruct > vOHcalHits);
0084 void CalESumAndCorrPosi(tracKuma trk, std::vector<hitStruct > vEmcalHits,\
0085 std::vector<hitStruct > vIHCalHits, std::vector<hitStruct > vOHcalHits);
0086 void TrackPropertiesEstimation2(tracKuma& trk);
0087
0088 Int_t TempINTTIOMatching(Double_t oINTTPhi, std::vector<hitStruct > vIInttHits);
0089 Double_t TempCalcdPhidR(Int_t iInttID, Int_t oInttID, \
0090 std::vector<hitStruct > vIInttHits, std::vector<hitStruct > vOInttHits);
0091 Int_t TempInttCalMatch(Int_t iInttID, Int_t oInttID,\
0092 std::vector<hitStruct > vIInttHits, std::vector<hitStruct > vOInttHits,\
0093 std::vector<hitStruct > vEmcalHits);
0094
0095 void SagittaRByCircleFit(Double_t& cernterX, Double_t& centerY, Double_t& sagittaR,
0096 std::vector<Double_t > r, std::vector<Double_t > phi , Double_t oInttPhi, Double_t emcalPhi);
0097
0098 bool FindHitsOnCircle(Int_t& hitMatchId, Double_t& hitX, Double_t& hitY, Double_t sagittaR, \
0099 Double_t centerX, Double_t centerY, std::vector<hitStruct > vHits, Double_t dRThre);
0100
0101 Double_t EstimateRecoTheta(tracKuma trk, Int_t type);
0102 void EstiVertex(tracKuma trk);
0103
0104 void ParticleIdentify(tracKuma trk);
0105
0106 bool ReturnHitsRPhiVect(std::vector<Double_t >& hitR, std::vector<Double_t >& hitPhi,\
0107 std::vector<Int_t > subDetSet, tracKuma trk);
0108
0109 void SetHitParaInTrack(tracKuma& trk, Int_t detId, hitStruct hitSt);
0110
0111
0112
0113 inline Double_t CalcSagittaPt(Double_t sagittaR){
0114 Double_t pT = 0.3*1.4*sagittaR*0.01;
0115 return pT;
0116 }
0117
0118 inline Double_t dPhiOInttEmcalEsti(tracKuma trk){
0119 Double_t xIIntt = trk.getHitR(4)*std::cos(trk.getHitPhi(4));
0120 Double_t yIIntt = trk.getHitR(4)*std::sin(trk.getHitPhi(4));
0121 Double_t xOIntt = trk.getHitR(5)*std::cos(trk.getHitPhi(5));
0122 Double_t yOIntt = trk.getHitR(5)*std::sin(trk.getHitPhi(5));
0123 Double_t xEmcal = trk.getHitR(6)*std::cos(trk.getHitPhi(6));
0124 Double_t yEmcal = trk.getHitR(6)*std::sin(trk.getHitPhi(6));
0125
0126 Double_t phiIInttOIntt = std::atan((yOIntt-yIIntt)/(xOIntt-xIIntt));
0127 Double_t phiOInttEmcal = std::atan((yEmcal-yOIntt)/(xEmcal-xOIntt));
0128
0129 Double_t dPhiOInttEmcal = phiOInttEmcal - phiIInttOIntt;
0130
0131 return dPhiOInttEmcal;
0132 }
0133
0134 inline Double_t dPhiVtxIInttEmcalEsti(tracKuma trk){
0135 Double_t xVtx = trk.getHitR(0)*std::cos(trk.getHitPhi(0));
0136 Double_t yVtx = trk.getHitR(0)*std::sin(trk.getHitPhi(0));
0137 Double_t xIIntt = trk.getHitR(4)*std::cos(trk.getHitPhi(4));
0138 Double_t yIIntt = trk.getHitR(4)*std::sin(trk.getHitPhi(4));
0139 Double_t xEmcal = trk.getHitR(6)*std::cos(trk.getHitPhi(6));
0140 Double_t yEmcal = trk.getHitR(6)*std::sin(trk.getHitPhi(6));
0141
0142 Double_t phiVtxIIntt = std::atan((yIIntt-yVtx)/(xIIntt-xVtx));
0143 Double_t phiIInttEmcal = std::atan((yEmcal-yIIntt)/(xEmcal-xIIntt));
0144
0145 Double_t dPhiIInttEmcal = phiIInttEmcal - phiVtxIIntt;
0146
0147 return dPhiIInttEmcal;
0148 }
0149
0150 inline Double_t dPhiVtxOInttEmcalEsti(tracKuma trk){
0151 Double_t xVtx = trk.getHitR(0)*std::cos(trk.getHitPhi(0));
0152 Double_t yVtx = trk.getHitR(0)*std::sin(trk.getHitPhi(0));
0153 Double_t xOIntt = trk.getHitR(5)*std::cos(trk.getHitPhi(5));
0154 Double_t yOIntt = trk.getHitR(5)*std::sin(trk.getHitPhi(5));
0155 Double_t xEmcal = trk.getHitR(6)*std::cos(trk.getHitPhi(6));
0156 Double_t yEmcal = trk.getHitR(6)*std::sin(trk.getHitPhi(6));
0157
0158 Double_t phiVtxOIntt = std::atan((yOIntt-yVtx)/(xOIntt-xVtx));
0159 Double_t phiOInttEmcal = std::atan((yEmcal-yOIntt)/(xEmcal-xOIntt));
0160
0161 Double_t dPhiOInttEmcal = phiOInttEmcal - phiVtxOIntt;
0162
0163 return dPhiOInttEmcal;
0164 }
0165
0166 inline Double_t dPhiVtxInttEmcalEsti(tracKuma trk){
0167 Double_t xVtx = trk.getHitR(0)*std::cos(trk.getHitPhi(0));
0168 Double_t yVtx = trk.getHitR(0)*std::sin(trk.getHitPhi(0));
0169 Double_t xIIntt = trk.getHitR(4)*std::cos(trk.getHitPhi(4));
0170 Double_t yIIntt = trk.getHitR(4)*std::sin(trk.getHitPhi(4));
0171 Double_t xOIntt = trk.getHitR(5)*std::cos(trk.getHitPhi(5));
0172 Double_t yOIntt = trk.getHitR(5)*std::sin(trk.getHitPhi(5));
0173 Double_t xEmcal = trk.getHitR(6)*std::cos(trk.getHitPhi(6));
0174 Double_t yEmcal = trk.getHitR(6)*std::sin(trk.getHitPhi(6));
0175
0176 Double_t phiVtxIIntt = std::atan((yIIntt-yVtx)/(xIIntt-xVtx));
0177 Double_t phiVtxOIntt = std::atan((yOIntt-yVtx)/(xOIntt-xVtx));
0178 Double_t phiVtxIntt = (phiVtxOIntt - phiVtxIIntt)/2;
0179 Double_t phiInttEmcal = std::atan((yEmcal-yOIntt)/(xEmcal-xOIntt));
0180
0181 Double_t dPhiInttEmcal = phiInttEmcal - phiVtxIntt;
0182
0183 return dPhiInttEmcal;
0184 }
0185
0186
0187 inline Double_t dPhiMvtxIInttEmcalEsti(tracKuma trk){
0188 Int_t numOfMvtx = 0;
0189 Double_t xMvtx = 0.;
0190 Double_t yMvtx = 0.;
0191 for(Int_t i = 1; i < 4; i++){
0192 if(trk.getHitIs(i)){
0193 xMvtx += trk.getHitR(i)*std::cos(trk.getHitPhi(i));
0194 yMvtx += trk.getHitR(i)*std::sin(trk.getHitPhi(i));
0195 numOfMvtx++;
0196 }
0197 }
0198 xMvtx /= numOfMvtx;
0199 yMvtx /= numOfMvtx;
0200
0201 Double_t xIIntt = trk.getHitR(4)*std::cos(trk.getHitPhi(4));
0202 Double_t yIIntt = trk.getHitR(4)*std::sin(trk.getHitPhi(4));
0203
0204 Double_t xEmcal = trk.getHitR(6)*std::cos(trk.getHitPhi(6));
0205 Double_t yEmcal = trk.getHitR(6)*std::sin(trk.getHitPhi(6));
0206
0207 Double_t phiMvtxIntt = std::atan((yIIntt-yMvtx)/(xIIntt-xMvtx));
0208 Double_t phiInttEmcal = std::atan((yEmcal-yIIntt)/(xEmcal-xIIntt));
0209
0210 Double_t dPhiIInttEmcal = phiInttEmcal - phiMvtxIntt;
0211
0212 return dPhiIInttEmcal;
0213 }
0214
0215 inline Double_t dPhiMvtxOInttEmcalEsti(tracKuma trk){
0216 Int_t numOfMvtx = 0;
0217 Double_t xMvtx = 0.;
0218 Double_t yMvtx = 0.;
0219 for(Int_t i = 1; i < 4; i++){
0220 if(trk.getHitIs(i)){
0221 xMvtx += trk.getHitR(i)*std::cos(trk.getHitPhi(i));
0222 yMvtx += trk.getHitR(i)*std::sin(trk.getHitPhi(i));
0223 numOfMvtx++;
0224 }
0225 }
0226 xMvtx /= numOfMvtx;
0227 yMvtx /= numOfMvtx;
0228
0229 Double_t xOIntt = trk.getHitR(5)*std::cos(trk.getHitPhi(5));
0230 Double_t yOIntt = trk.getHitR(5)*std::sin(trk.getHitPhi(5));
0231
0232 Double_t xEmcal = trk.getHitR(6)*std::cos(trk.getHitPhi(6));
0233 Double_t yEmcal = trk.getHitR(6)*std::sin(trk.getHitPhi(6));
0234
0235 Double_t phiMvtxOIntt = std::atan((yOIntt-yMvtx)/(xOIntt-xMvtx));
0236 Double_t phiOInttEmcal = std::atan((yEmcal-yOIntt)/(xEmcal-xOIntt));
0237
0238 Double_t dPhiOInttEmcal = phiOInttEmcal - phiMvtxOIntt;
0239
0240 return dPhiOInttEmcal;
0241 }
0242
0243 inline Double_t dPhiMvtxInttEmcalEsti(tracKuma trk){
0244 Int_t numOfMvtx = 0;
0245 Double_t xMvtx = 0.;
0246 Double_t yMvtx = 0.;
0247 for(Int_t i = 1; i < 4; i++){
0248 if(trk.getHitIs(i)){
0249 xMvtx += trk.getHitR(i)*std::cos(trk.getHitPhi(i));
0250 yMvtx += trk.getHitR(i)*std::sin(trk.getHitPhi(i));
0251 numOfMvtx++;
0252 }
0253 }
0254 xMvtx /= numOfMvtx;
0255 yMvtx /= numOfMvtx;
0256
0257 Double_t xIIntt = trk.getHitR(4)*std::cos(trk.getHitPhi(4));
0258 Double_t yIIntt = trk.getHitR(4)*std::sin(trk.getHitPhi(4));
0259 Double_t xOIntt = trk.getHitR(5)*std::cos(trk.getHitPhi(5));
0260 Double_t yOIntt = trk.getHitR(5)*std::sin(trk.getHitPhi(5));
0261 Double_t xIntt = (xIIntt + xOIntt)/2;
0262 Double_t yIntt = (yIIntt + yOIntt)/2;
0263
0264 Double_t xEmcal = trk.getHitR(6)*std::cos(trk.getHitPhi(6));
0265 Double_t yEmcal = trk.getHitR(6)*std::sin(trk.getHitPhi(6));
0266
0267 Double_t phiMvtxIntt = std::atan((yIntt-yMvtx)/(xIntt-xMvtx));
0268 Double_t phiInttEmcal = std::atan((yEmcal-yIntt)/(xEmcal-xIntt));
0269
0270 Double_t dPhiInttEmcal = phiInttEmcal - phiMvtxIntt;
0271
0272 return dPhiInttEmcal;
0273 }
0274
0275
0276 inline Double_t FitFunctionPt(Double_t dPhi){
0277 Double_t pT = 0.0249291 + 0.232799/dPhi -0.000752825/(dPhi*dPhi);
0278 return pT;
0279 }
0280
0281 inline Double_t FitFunctionPt_VtxIInttEmcal(Double_t dPhi){
0282 Double_t pT = -0.075123 + 0.241657 /dPhi + 7.81643e-05/(dPhi*dPhi);
0283 return pT;
0284 }
0285
0286 inline Double_t FitFunctionPt_VtxOInttEmcal(Double_t dPhi){
0287 Double_t pT = 0.0592938 + 0.210859/dPhi - 0.00157569/(dPhi*dPhi);
0288 return pT;
0289 }
0290
0291 inline Double_t FitFunctionPt_VtxInttEmcal(Double_t dPhi){
0292 Double_t pT = 0.089305 + 0.202456/dPhi + 0.0018411/(dPhi*dPhi);
0293 return pT;
0294 }
0295
0296 inline Double_t FitFunctionPt_MVtxIInttEmcal(Double_t dPhi){
0297 Double_t pT = -0.125412 + 0.293345/dPhi - 0.00233657/(dPhi*dPhi);
0298 return pT;
0299 }
0300
0301 inline Double_t FitFunctionPt_MVtxOInttEmcal(Double_t dPhi){
0302 Double_t pT = -0.0708943 + 0.269396/dPhi - 0.00158606/(dPhi*dPhi);
0303 return pT;
0304 }
0305
0306 inline Double_t FitFunctionPt_MVtxInttEmcal(Double_t dPhi){
0307 Double_t pT = -0.0645267 + 0.268115/dPhi -0.00156531/(dPhi*dPhi);
0308 return pT;
0309 }
0310
0311
0312 inline void Set3PointsXY(Double_t (&HitXY)[3][2], tracKuma trk, Int_t type)
0313 {
0314 if(type==0)
0315 {
0316 HitXY[0][0] = 0.;
0317 HitXY[0][1] = 0.;
0318 HitXY[1][0] = trk.getHitR(4)*cos(trk.getHitPhi(4));
0319 HitXY[1][1] = trk.getHitR(4)*sin(trk.getHitPhi(4));
0320 HitXY[2][0] = trk.getHitR(5)*cos(trk.getHitPhi(5));
0321 HitXY[2][1] = trk.getHitR(5)*sin(trk.getHitPhi(5));
0322 }
0323 else if(type==1)
0324 {
0325 HitXY[0][0] = trk.getHitR(4)*cos(trk.getHitPhi(4));
0326 HitXY[0][1] = trk.getHitR(4)*sin(trk.getHitPhi(4));
0327 HitXY[1][0] = trk.getHitR(5)*cos(trk.getHitPhi(5));
0328 HitXY[1][1] = trk.getHitR(5)*sin(trk.getHitPhi(5));
0329 HitXY[2][0] = trk.getHitR(6)*cos(trk.getHitPhi(6));
0330 HitXY[2][1] = trk.getHitR(6)*sin(trk.getHitPhi(6));
0331 }
0332 else if(type==2)
0333 {
0334 HitXY[0][0] = 0.;
0335 HitXY[0][1] = 0.;
0336 HitXY[1][0] = trk.getHitR(5)*cos(trk.getHitPhi(5));
0337 HitXY[1][1] = trk.getHitR(5)*sin(trk.getHitPhi(5));
0338 HitXY[2][0] = trk.getHitR(6)*cos(trk.getHitPhi(6));
0339 HitXY[2][1] = trk.getHitR(6)*sin(trk.getHitPhi(6));
0340 }
0341 }
0342
0343 inline void CalcPerpBis(Double_t& slope, Double_t& sector,
0344 Double_t HitXY1[2], Double_t HitXY2[2]){
0345
0346 Double_t tempSlope = (HitXY2[1] - HitXY1[1])/(HitXY2[0] - HitXY1[0]);
0347 Double_t tempSector = HitXY1[1] - tempSlope*HitXY1[0];
0348
0349 Double_t centerX = (HitXY2[0] + HitXY1[0])/2;
0350 Double_t centerY = (HitXY2[1] + HitXY1[1])/2;
0351
0352 slope = -1/tempSlope;
0353 sector = centerY - slope*centerX;
0354
0355 return;
0356 }
0357
0358 inline void CalcLineEq(Double_t& slope, Double_t& section,\
0359 Double_t x1, Double_t y1, Double_t x2, Double_t y2){
0360 slope = (y2 - y1)/(x2 - x1);
0361 section = y1 - slope*x1;
0362 }
0363
0364 inline void CalcCrossPoint(Double_t& crossX, Double_t& crossY,\
0365 Double_t origSlope, Double_t origSection, Double_t nextPointX, Double_t nextPointY){
0366 Double_t invSlope = - (1/origSlope);
0367 Double_t invSection = nextPointY - invSlope*nextPointX;
0368 crossX = - (invSection - origSection)/(invSlope - origSlope);
0369 crossY = origSlope*crossX + origSection;
0370 }
0371
0372 inline Double_t CalcSection(Double_t x, Double_t y, Double_t slope){
0373 Double_t section = y - slope*x;
0374 return section;
0375 }
0376
0377 inline void CrossLineCircle(Double_t& xv, Double_t& yv,\
0378 Double_t xc, Double_t yc, Double_t R){
0379
0380
0381
0382
0383 Double_t slope = yc/xc;
0384
0385 Double_t b = (xc - slope*yc)/(slope*slope + 1);
0386 Double_t c = (yc*yc - R*R + xc*xc)/(slope*slope + 1);
0387
0388 Double_t xv1 = -b + std::sqrt(b*b - c);
0389 Double_t xv2 = -b - std::sqrt(b*b - c);
0390 Double_t yv1 = slope*xv1;
0391 Double_t yv2 = slope*xv2;
0392
0393 if((xv1*xv1+yv1*yv1) < (xv2*xv2+yv2*yv2)) xv = xv1;
0394 else xv = xv2;
0395
0396 yv = slope*xv;
0397
0398 return;
0399 }
0400
0401 inline Double_t CrossCircleCircle(Double_t& xcal, Double_t& ycal, Double_t xc, Double_t yc, Double_t R, Double_t phiIntt)
0402 {
0403
0404 Double_t slope = -1*(xc/yc);
0405 Double_t section = (xc*xc+yc*yc-R*R)/(2*yc);
0406
0407 Double_t coeff1 = (slope*section)/(slope*slope + 1);
0408 Double_t coeff2 = (section*section - m_ECalR*m_ECalR)/(slope*slope + 1);
0409
0410 Double_t x1 = -1 * coeff1 + std::sqrt(coeff1*coeff1 - coeff2);
0411 Double_t x2 = -1 * coeff1 - std::sqrt(coeff1*coeff1 - coeff2);
0412
0413 Double_t y1 = slope*x1 + section;
0414 Double_t y2 = slope*x2 + section;
0415
0416 Double_t phi1 = std::atan(y1/x1);
0417 if((phi1 < 0)&&(x1 < 0)) phi1 += TMath::Pi();
0418 else if((phi1 > 0)&&(x1 < 0)) phi1 -= TMath::Pi();
0419
0420 Double_t phi2 = std::atan(y2/x2);
0421 if((phi2 < 0)&&(x2 < 0)) phi2 += TMath::Pi();
0422 else if((phi2 > 0)&&(x2 < 0)) phi2 -= TMath::Pi();
0423
0424 Double_t calPhi = 0.;
0425 if((std::abs(phi1-phiIntt)) < (std::abs(phi2-phiIntt)))
0426 {
0427 xcal = x1;
0428 ycal = y1;
0429 calPhi = phi1;
0430 }
0431 else
0432 {
0433 xcal = x2;
0434 ycal = y2;
0435 calPhi = phi2;
0436 }
0437
0438 return calPhi;
0439 }
0440
0441
0442
0443 private:
0444 Double_t m_ECalR = 93.5;
0445
0446 Double_t m_iHCalR = 127.502;
0447 Double_t m_oHCalR = 225.87;
0448
0449 Double_t m_InttMatchPhiMin = -0.05;
0450 Double_t m_InttMatchPhiMax = 0.05;
0451
0452 Double_t m_EcalEThre = 0.1;
0453 Double_t m_iHcalEThre = 5;
0454 Double_t m_oHcalEThre = 5;
0455
0456 Double_t m_EcalSearchPhiMin = 0.;
0457 Double_t m_EcalSearchPhiMax = 0.;
0458
0459 Double_t m_EOverPElectron = 0.2;
0460
0461 std::vector<hitStruct > m_fMvtxHits;
0462 std::vector<hitStruct > m_sMvtxHits;
0463 std::vector<hitStruct > m_tMvtxHits;
0464
0465 std::vector<hitStruct > m_iInttHits;
0466 std::vector<hitStruct > m_oInttHits;
0467
0468 std::vector<hitStruct > m_emcalHits;
0469 std::vector<hitStruct > m_iHCalHits;
0470 std::vector<hitStruct > m_oHCalHits;
0471
0472 };
0473
0474
0475
0476 #endif