00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055 #include <iostream>
00056 #include "StEventTypes.h"
00057 #include "Stypes.h"
00058 #include "StThreeVectorD.hh"
00059 #include "StThreeVectorF.hh"
00060 #include "StHelix.hh"
00061 #include "StTrackGeometry.h"
00062 #include "StEventUtilities/StuRefMult.hh"
00063 #include "PhysicalConstants.h"
00064 #include "StPhysicalHelixD.hh"
00065 #include "StTofUtil/StTofGeometry.h"
00066 #include "StTofUtil/StTofSlatCollection.h"
00067 #include "tables/St_pvpdStrobeDef_Table.h"
00068 #include "TH1.h"
00069 #include "TH2.h"
00070 #include "TFile.h"
00071 #include "TOrdCollection.h"
00072 #include "StMessMgr.h"
00073 #include "StMemoryInfo.hh"
00074 #include "StTimer.hh"
00075 #include "StTofUtil/tofPathLength.hh"
00076 #include "StTofpMatchMaker.h"
00077
00078
00079 ClassImp(StTofpMatchMaker)
00080
00081
00083 StTofpMatchMaker::StTofpMatchMaker(const Char_t *name): StMaker(name){
00084
00085 mEventCounter = 0;
00086 mAcceptedEventCounter = 0;
00087 mTofEventCounter = 0;
00088 mTofStrobeEventCounter = 0;
00089 mAcceptAndStrobe = 0;
00090 mAcceptAndBeam = 0;
00091 mTofGeom = 0;
00092
00093
00094 setValidAdcRange(0,1024);
00095 setValidTdcRange(30,1200);
00096 setOuterTrackGeometry();
00097 setMinHitsPerTrack(0);
00098 setMinFitPointsPerTrack(0);
00099 setMaxDCA(9999.);
00100
00101 createHistograms(kTRUE);
00102 doPrintMemoryInfo = kFALSE;
00103 doPrintCpuInfo = kFALSE;
00104 }
00105
00106 StTofpMatchMaker::~StTofpMatchMaker(){ }
00107
00108 void StTofpMatchMaker::Clear(Option_t *opt){StMaker::Clear();}
00109
00110
00112 Int_t StTofpMatchMaker::Init(){
00113 gMessMgr->Info("StTofpMatchMaker -- initializing ...","OS");
00114 gMessMgr->Info("","OST") << "Valid TDC range: " << mMinValidTdc << " " << mMaxValidTdc << endm;
00115 gMessMgr->Info("","OST") << "Valid ADC range: " << mMinValidAdc << " " << mMaxValidAdc << endm;
00116 gMessMgr->Info("","OST") << "Minimum hits per track: " << mMinHitsPerTrack << endm;
00117 gMessMgr->Info("","OST") << "Minimum fitpoints per track: " << mMinFitPointsPerTrack << endm;
00118 gMessMgr->Info("","OST") << "Maximum DCA: " << mMaxDCA << endm;
00119 if (!mOuterTrackGeometry)
00120 gMessMgr->Warning("Warning: using standard trackgeometry()","OST");
00121
00122 if(m_Mode) {
00123 setHistoFileName("tofana.root");
00124 } else {
00125 setHistoFileName("");
00126 }
00127
00128 if (mHisto){
00129 bookHistograms();
00130 gMessMgr->Info("Histograms are booked","OST");
00131 if (mHistoFileName!="")
00132 gMessMgr->Info("","OST") << "Histograms will be stored in " << mHistoFileName.c_str() << endm;
00133 }
00134
00135
00136 mEventCounter = 0;
00137 mAcceptedEventCounter = 0;
00138 mTofEventCounter = 0;
00139 mTofStrobeEventCounter = 0;
00140 mAcceptAndStrobe = 0;
00141 mAcceptAndBeam = 0;
00142
00143 return kStOK;
00144 }
00145
00146
00147
00149 Int_t StTofpMatchMaker::InitRun(int runnumber){
00150
00151
00152 mYear2 = (runnumber<4000000);
00153 mYear3 = (runnumber>4000000&&runnumber<5000000);
00154 mYear4 = (runnumber>5000000&&runnumber<6000000);
00155 mYear5 = (runnumber>6000000);
00156
00157 gMessMgr->Info("StTofpMatchMaker -- reinitializing TofGeometry (InitRun)","OS" );
00158 mTofGeom = new StTofGeometry();
00159 mTofGeom->init(this);
00160
00161 gMessMgr->Info(" -- retrieving run parameters","OS");
00162 TDataSet *mDbDataSet = GetDataBase("Calibrations/tof");
00163 if (!mDbDataSet){
00164 gMessMgr->Error("unable to get TOF run parameters","OS");
00165
00166 return kStErr;
00167 }
00168 St_pvpdStrobeDef* pvpdStrobeDef = static_cast<St_pvpdStrobeDef*>(mDbDataSet->Find("pvpdStrobeDef"));
00169 if (!pvpdStrobeDef){
00170 gMessMgr->Error("unable to find TOF run param table","OS");
00171
00172 return kStErr;
00173 }
00174 pvpdStrobeDef_st *strobeDef = static_cast<pvpdStrobeDef_st*>(pvpdStrobeDef->GetArray());
00175 int numRows = pvpdStrobeDef->GetNRows();
00176 if (NPVPD != numRows) gMessMgr->Warning("#tubes inconsistency in dbase");
00177 for (int i=0;i<NPVPD;i++){
00178 int ii = strobeDef[i].id - 1;
00179 mStrobeTdcMin[ii] = strobeDef[i].strobeTdcMin;
00180 mStrobeTdcMax[ii] = strobeDef[i].strobeTdcMax;
00181 if (Debug())
00182 LOG_INFO << "tube " << strobeDef[i].id << " min:"<< strobeDef[i].strobeTdcMin
00183 <<" max:"<< strobeDef[i].strobeTdcMax<< endm;
00184 }
00185
00186
00187
00188 return kStOK;
00189 }
00190
00192 Int_t StTofpMatchMaker::FinishRun(int runnumber){
00193 gMessMgr->Info("StTofpMatchMaker -- cleaning up geometry (FinishRun)","OS" );
00194 if (mTofGeom) delete mTofGeom;
00195 mTofGeom=0;
00196 return kStOK;
00197 }
00198
00199
00200
00202 Int_t StTofpMatchMaker::Finish(){
00203 gMessMgr->Info("","OS") << "StTofpMatchMaker ----- RUN SUMMARY ----- (Finish)\n"
00204 << "\tProcessed " << mEventCounter << " events."
00205 << " Accepted " << mAcceptedEventCounter << " events."
00206 << " Rejected " << mEventCounter - mAcceptedEventCounter << " events\n"
00207 << "\tTOF events " << mTofEventCounter
00208 << ". Beam " << mTofEventCounter - mTofStrobeEventCounter
00209 << " Strobe " << mTofStrobeEventCounter
00210 << "\n\t Accept & Strobe " << mAcceptAndStrobe << " events\n"
00211 << "\t Accept & Beam " << mAcceptAndBeam << " events" << endm;
00212
00213 if (mHistoFileName!="") writeHistogramsToFile();
00214 return kStOK;
00215 }
00216
00217
00218
00220 Int_t StTofpMatchMaker::Make(){
00221 gMessMgr->Info("StTofpMatchMaker -- welcome","OS");
00222
00223 if(mYear5) {
00224 gMessMgr->Info("StTofpMatchMaker -- no TOFp in and after Run 5","OS");
00225 return kStOK;
00226 }
00227
00228
00229 StEvent *event = (StEvent *) GetInputDS("StEvent");
00230 if (!validEvent(event)){
00231 gMessMgr->Info("StTofpMatchMaker -- nothing to do ... bye-bye","OST");
00232 return kStOK;
00233 }
00234
00235
00236 StTimer timer;
00237 if (doPrintCpuInfo) timer.start();
00238 if (doPrintMemoryInfo) StMemoryInfo::instance()->snapshot();
00239
00240
00241
00242 StTofCollection *theTof = event->tofCollection();
00243 getTofData(theTof);
00244
00245
00246
00247 float sumAdcTofp(0.); int nAdcTofp(0), nTdcTofp(0), nAdcTdcTofp(0);
00248 for (int i=0;i<NTOFP;i++){
00249 sumAdcTofp += mTofpAdc[i];
00250 bool tdcValid = validTdc(mTofpTdc[i]);
00251 bool adcValid = validAdc(mTofpAdc[i]);
00252 if (tdcValid) nTdcTofp++;
00253 if (adcValid) nAdcTofp++;
00254 if (adcValid && tdcValid) nAdcTdcTofp++;
00255 }
00256 if (Debug())
00257 LOG_INFO << " TOFp #Adc:" << nAdcTofp << " #Tdc:" << nTdcTofp << endm;
00258
00259
00260 float sumAdcPvpd=0; int nAdcPvpd=0, nTdcPvpd=0;
00261 for (int i=0;i<NPVPD;i++){
00262 sumAdcPvpd += mPvpdAdc[i];
00263 bool tdcValid = validTdc(mPvpdTdc[i]);
00264 bool adcValid = validAdc(mPvpdAdc[i]);
00265 if (tdcValid) nTdcPvpd++;
00266 if (adcValid) nAdcPvpd++;
00267 }
00268 if (Debug())
00269 LOG_INFO << " pVPD #Adc:" << nAdcPvpd << " #Tdc:" << nTdcPvpd << endm;
00270
00271
00272
00273 if (Debug()){
00274
00275
00276 int refmult(0);
00277 bool richTofMuDST = (event->summary()->numberOfExoticTracks() == -999);
00278 if (richTofMuDST)
00279 refmult = event->summary()->numberOfGoodTracks();
00280 else
00281 refmult = uncorrectedNumberOfPrimaries(*event);
00282
00283 LOG_INFO << " #Tracks :" << event->summary()->numberOfTracks()
00284 << "\n #goodPrimaryTracks:" << event->summary()->numberOfGoodPrimaryTracks()
00285 << "\n #uncorr.prim.tracks :" << refmult << endm;
00286 if (!richTofMuDST)
00287 LOG_INFO << " #goodTracks (global):" << event->summary()->numberOfGoodTracks() << endm;
00288 }
00289
00290
00291
00292
00293
00294 idVector validSlatIdVec;
00295 for (int i=0;i<NTOFP;i++){
00296 unsigned short slatId = mTofGeom->daqToSlatId(i);
00297 float rawAdc = mTofpAdc[i];
00298 float rawTdc = mTofpTdc[i];
00299 if (mHisto) hTofpSlatIdA0->Fill(i+1);
00300
00301 if (validAdc(rawAdc) && validTdc(rawTdc)){
00302 validSlatIdVec.push_back(slatId);
00303 if (mHisto) hTofpSlatIdA1->Fill(i+1);
00304 }
00305 }
00306 gMessMgr->Info("","OST") << "A: #valid slats: " << validSlatIdVec.size() << endm;
00307
00308 if(!validSlatIdVec.size()) return kStOK;
00309
00310
00311
00312
00313 tofSlatHitVector allSlatsHitVec;
00314 allSlatsHitVec.clear();
00315 StSPtrVecTrackNode& nodes = event->trackNodes();
00316 int nAllTracks=0;
00317 for (unsigned int iNode=0; iNode<nodes.size(); iNode++){
00318 tofSlatHitVector slatHitVec;
00319 slatHitVec.clear();
00320 StTrack *theTrack = nodes[iNode]->track(global);
00321
00322
00323 if (validTrack(theTrack)){
00324 nAllTracks++;
00325 StPhysicalHelixD theHelix = trackGeometry(theTrack)->helix();
00326
00327 idVector projTrayVec;
00328 if(!mTofGeom->projTrayVector(theHelix, projTrayVec)) continue;
00329
00330 Bool_t hitTofp = kFALSE;
00331 for(size_t ii = 0; ii<projTrayVec.size(); ii++) {
00332 if(projTrayVec[ii]==mTofpTrayId) hitTofp = kTRUE;
00333 }
00334 if(!hitTofp) continue;
00335
00336 slatHitVec = mTofGeom->tofHelixToArray(theHelix, validSlatIdVec);
00337 if (slatHitVec.size()>0 && mHisto) hTofpSlatHitVecSize->Fill(slatHitVec.size());
00338
00339 for (size_t ii = 0; ii < slatHitVec.size(); ii++) {
00340 slatHitVec[ii].trackIdVec.push_back(iNode);
00341 allSlatsHitVec.push_back(slatHitVec[ii]);
00342 if (mHisto){
00343 int id = mTofGeom->slatIdToDaq(slatHitVec[ii].slatIndex);
00344 float xhit = slatHitVec[ii].hitPosition.x();
00345 float yhit = slatHitVec[ii].hitPosition.y();
00346 float zhit = slatHitVec[ii].hitPosition.z();
00347 float phihit = atan2(yhit,xhit);
00348 hTofpSlatIdB1->Fill(id);
00349 hTofpHitMap1->Fill(zhit,phihit);
00350 }
00351
00352 if (Debug()){
00353 LOG_INFO << "B: trackid=";
00354 idVectorIter ij = slatHitVec[ii].trackIdVec.begin();
00355 while (ij != slatHitVec[ii].trackIdVec.end()) {LOG_INFO << " " << *ij; ij++;}
00356 LOG_INFO << "\tind=" << mTofGeom->slatIdToDaq(slatHitVec[ii].slatIndex)
00357 << "\thitprof="<< slatHitVec[ii].hitProfile
00358 << "\ts="<<slatHitVec[ii].s << "\tthxy="<<slatHitVec[ii].theta_xy
00359 << "\tthzr="<<slatHitVec[ii].theta_zr;
00360 if (slatHitVec.size()>1) LOG_INFO << " M" << endm;
00361 else LOG_INFO << endm;
00362 }
00363 }
00364 }
00365 }
00366
00367 gMessMgr->Info("","OST") << "B: #matched/#avail/#total tracknodes: "
00368 <<allSlatsHitVec.size() << "/" <<nAllTracks
00369 << "/" << nodes.size() << endm;
00370
00371 if(!allSlatsHitVec.size()) return kStOK;
00372
00373
00374
00375
00376 if (mHisto){
00377 for (tofSlatHitVectorIter ij = allSlatsHitVec.begin();ij!=allSlatsHitVec.end(); ij++){
00378 int slatId=ij->slatIndex;
00379 int daqId= mTofGeom->slatIdToDaq(slatId);
00380
00381 idVector slatNeighbours = mTofGeom->slatNeighboursWide(slatId);
00382 idVector slatCloseNeighbours = mTofGeom->slatNeighbours(slatId);
00383
00384
00385 bool matchedNeighbours(false);
00386 for (idVectorIter kk=slatNeighbours.begin();kk!=slatNeighbours.end();kk++){
00387
00388 for (tofSlatHitVectorIter jj = allSlatsHitVec.begin();jj!=allSlatsHitVec.end();jj++){
00389 if (jj->slatIndex == *kk) matchedNeighbours=true;
00390 }
00391 }
00392
00393
00394 if (!matchedNeighbours){
00395 bool tdcOK = validTdc(mTofpTdc[daqId-1]);
00396 int nHitNeighbours(0),nNoHitNeighbours(0);
00397
00398 for (idVectorIter kk=slatCloseNeighbours.begin();kk!=slatCloseNeighbours.end();kk++){
00399 int neighborDaqId=mTofGeom->slatIdToDaq(*kk);
00400 int neighborSlatId=*kk;
00401 if (validTdc(mTofpTdc[neighborDaqId-1])){
00402 int iEta = mTofGeom->tofSlat(neighborSlatId).ieta - mTofGeom->tofSlat(slatId).ieta;
00403 int iPhi = mTofGeom->tofSlat(neighborSlatId).iphi - mTofGeom->tofSlat(slatId).iphi;
00404
00405
00406 if (tdcOK){
00407 hTofpMatchHit[daqId-1]->Fill(iEta,iPhi);
00408 nHitNeighbours++;
00409 }
00410 else {
00411 hTofpMatchNoHit[daqId-1]->Fill(iEta,iPhi);
00412 nNoHitNeighbours++;
00413 }
00414 }
00415 }
00416
00417 if (tdcOK){
00418 hTofpMatchHit[daqId-1]->Fill(0.,0.,nHitNeighbours);
00419 }
00420 else{
00421 hTofpMatchNoHit[daqId-1]->Fill(0.,0.,nNoHitNeighbours);
00422 }
00423
00424 }
00425 }
00426 }
00427
00428
00429
00430
00431
00432
00433
00434 int nSingleHitSlats(0);
00435 tofSlatHitVector singleHitSlatsVec;
00436 StructSlatHit slatHit;
00437
00438 tofSlatHitVector tempVec = allSlatsHitVec;
00439 tofSlatHitVector erasedVec = tempVec;
00440 while (tempVec.size() != 0) {
00441 int nTracks = 0;
00442 vector<StThreeVectorD> vPosition;
00443 vector<vector<StThreeVectorD> > vLayerHitPositions;
00444 vector<Int_t> vHitProfile;
00445 vector<Float_t> vS, vTheta_xy, vTheta_zr;
00446 idVector trackIdVec;
00447
00448 tofSlatHitVectorIter tempIter=tempVec.begin();
00449 tofSlatHitVectorIter erasedIter=erasedVec.begin();
00450 while(erasedIter!= erasedVec.end()) {
00451 if(tempIter->slatIndex == erasedIter->slatIndex) {
00452 nTracks++;
00453
00454 trackIdVec.push_back(erasedIter->trackIdVec.back());
00455 vPosition.push_back(erasedIter->hitPosition);
00456 vLayerHitPositions.push_back(erasedIter->layerHitPositions);
00457 vHitProfile.push_back(erasedIter->hitProfile);
00458 vS.push_back(erasedIter->s);
00459 vTheta_xy.push_back(erasedIter->theta_xy);
00460 vTheta_zr.push_back(erasedIter->theta_zr);
00461
00462 if (mHisto){
00463 float xhit = erasedIter->hitPosition.x();
00464 float yhit = erasedIter->hitPosition.y();
00465 float zhit = erasedIter->hitPosition.z();
00466 float phihit = atan2(yhit,xhit);
00467 hTofpHitMap2->Fill(zhit,phihit);
00468 }
00469
00470 erasedVec.erase(erasedIter);
00471 erasedIter--;
00472 }
00473 erasedIter++;
00474 }
00475
00476 if (mHisto){
00477 int daqId = mTofGeom->slatIdToDaq(tempIter->slatIndex);
00478 hTofpSlatIdD1->Fill(daqId);
00479 }
00480
00481 if (nTracks==1){
00482 nSingleHitSlats++;
00483
00484 slatHit.slatIndex = tempIter->slatIndex;
00485 slatHit.hitPosition = vPosition[0];
00486 slatHit.layerHitPositions = vLayerHitPositions[0];
00487 slatHit.trackIdVec = trackIdVec;
00488 slatHit.hitProfile = vHitProfile[0];
00489 slatHit.s = vS[0];
00490 slatHit.theta_xy = vTheta_xy[0];
00491 slatHit.theta_zr = vTheta_zr[0];
00492
00493 singleHitSlatsVec.push_back(slatHit);
00494
00495 if (mHisto){
00496 int daqId = mTofGeom->slatIdToDaq(tempIter->slatIndex);
00497 float xhit = slatHit.hitPosition.x();
00498 float yhit = slatHit.hitPosition.y();
00499 float zhit = slatHit.hitPosition.z();
00500 float phihit = atan2(yhit,xhit);
00501 hTofpSlatIdD2->Fill(daqId);
00502 hTofpHitMap3->Fill(zhit,phihit);
00503 }
00504
00505
00506 if (Debug()) {
00507 LOG_INFO << "D: ind=" << mTofGeom->slatIdToDaq(slatHit.slatIndex)
00508 << "\thitprof="<< slatHit.hitProfile << "\ts="<<slatHit.s
00509 << "\tthxy="<<slatHit.theta_xy << "\tthzr="<<slatHit.theta_zr << "\ttrackid:";
00510 idVectorIter ij=trackIdVec.begin();
00511 while (ij != trackIdVec.end()) { LOG_INFO << " " << *ij; ij++; }
00512 LOG_INFO <<endm;
00513 }
00514 }
00515 else if (nTracks>1){
00516
00517
00518 } else
00519 gMessMgr->Warning("","OST") << "D: no tracks extrapolate to matched slat ... should not happen!" << endm;
00520
00521 tempVec = erasedVec;
00522 }
00523 gMessMgr->Info("","OST") << "D: #before/#after: " << allSlatsHitVec.size()
00524 << "/" << singleHitSlatsVec.size() << endm;
00525
00526
00527
00528
00529
00530
00531
00532 tofSlatHitVector allMatchedSlatsVec;
00533 tempVec = singleHitSlatsVec;
00534 erasedVec = tempVec;
00535 while (tempVec.size() != 0) {
00536 StructSlatHit slatHit;
00537 int nSlats = 0;
00538 vector<StThreeVectorD> vPosition;
00539 vector< vector<StThreeVectorD> > vLayerHitPositions;
00540 vector<Int_t> vHitProfile;
00541 vector<Float_t> vS, vTheta_xy, vTheta_zr;
00542 idVector vTrackId;
00543 vector<Int_t> slatIndex;
00544
00545 tofSlatHitVectorIter tempIter=tempVec.begin();
00546 tofSlatHitVectorIter erasedIter=erasedVec.begin();
00547 while(erasedIter!= erasedVec.end()) {
00548 if(tempIter->trackIdVec.back() == erasedIter->trackIdVec.back()) {
00549 nSlats++;
00550
00551 slatIndex.push_back(erasedIter->slatIndex);
00552 vTrackId.push_back(erasedIter->trackIdVec.back());
00553 vPosition.push_back(erasedIter->hitPosition);
00554 vLayerHitPositions.push_back(erasedIter->layerHitPositions);
00555 vHitProfile.push_back(erasedIter->hitProfile);
00556 vS.push_back(erasedIter->s);
00557 vTheta_xy.push_back(erasedIter->theta_xy);
00558 vTheta_zr.push_back(erasedIter->theta_zr);
00559
00560 erasedVec.erase(erasedIter);
00561 erasedIter--;
00562 }
00563 erasedIter++;
00564 }
00565
00566
00567 if (nSlats==1){
00568
00569 slatHit.slatIndex = slatIndex[0];
00570 slatHit.hitPosition = vPosition[0];
00571 slatHit.layerHitPositions = vLayerHitPositions[0];
00572 slatHit.trackIdVec.push_back(vTrackId[0]);
00573 slatHit.hitProfile = vHitProfile[0];
00574 slatHit.s = vS[0];
00575 slatHit.theta_xy = vTheta_xy[0];
00576 slatHit.theta_zr = vTheta_zr[0];
00577 slatHit.matchFlag = 0;
00578
00579 allMatchedSlatsVec.push_back(slatHit);
00580
00581 if (mHisto){
00582 int daqId = mTofGeom->slatIdToDaq(slatIndex[0]);
00583 hTofpSlatIdE1->Fill(daqId);
00584 }
00585
00586
00587 if (Debug()) {
00588 LOG_INFO << "E: ind=" << mTofGeom->slatIdToDaq(slatHit.slatIndex)
00589 << "\thitprof="<< slatHit.hitProfile << "\ts="<<slatHit.s
00590 << "\tthxy="<<slatHit.theta_xy << "\tthzr="<<slatHit.theta_zr << "\ttrackid:";
00591 idVectorIter ij=vTrackId.begin();
00592 while (ij != vTrackId.end()) { LOG_INFO << " " << *ij; ij++; }
00593 LOG_INFO <<endm;
00594 }
00595 }
00596 else if (nSlats>1){
00597 int thiscandidate(-99);
00598 int thisMatchFlag(0);
00599
00600
00601 int weight(0);
00602 vector<int> weightCandidates;
00603 thisMatchFlag = 1;
00604 if (Debug()) LOG_INFO << "E: find ... weight ";
00605 for (int i=0;i<nSlats;i++){
00606 int hitWeight = vLayerHitPositions[i].size();
00607 if (Debug()) LOG_INFO << mTofGeom->slatIdToDaq(slatIndex[i]) << "("<<hitWeight<<")"<<" ";
00608 if (hitWeight>weight) {
00609 weight=hitWeight;
00610 weightCandidates.clear();
00611 weightCandidates.push_back(i);
00612 } else if (hitWeight == weight)
00613 weightCandidates.push_back(i);
00614 }
00615 if (weightCandidates.size()==1){
00616 thiscandidate = weightCandidates[0];
00617 int daqId = mTofGeom->slatIdToDaq(slatIndex[thiscandidate]);
00618 if (mHisto) hTofpSlatIdE2->Fill(daqId);
00619 if (Debug()) LOG_INFO << "candidate =" << daqId << endm;
00620 }
00621
00622
00623 if (weightCandidates.size()>1){
00624 Float_t ss(0);
00625 vector<int> ssCandidates;
00626 thisMatchFlag = 2;
00627 if (Debug()) LOG_INFO << " ss ";
00628 for (unsigned int i=0;i<weightCandidates.size();i++){
00629 int ii=weightCandidates[i];
00630 if (Debug()) LOG_INFO << mTofGeom->slatIdToDaq(slatIndex[ii]) << " ";
00631 if (vS[ii]>ss){
00632 ss = vS[ii];
00633 ssCandidates.clear();
00634 ssCandidates.push_back(ii);
00635 }else if (vS[ii]==ss)
00636 ssCandidates.push_back(ii);
00637 }
00638 if (ssCandidates.size()==1){
00639 thiscandidate = ssCandidates[0];
00640 int daqId = mTofGeom->slatIdToDaq(slatIndex[thiscandidate]);
00641 if (mHisto) hTofpSlatIdE3->Fill(daqId);
00642 if (Debug()) LOG_INFO << "candidate =" << daqId << endm;
00643 }
00644
00645
00646 if (ssCandidates.size()>1){
00647 Int_t hitprof(0);
00648 vector<int> profileCandidates;
00649 thisMatchFlag = 3;
00650 if (Debug()) LOG_INFO << " hprof ";
00651 for (unsigned int i=0;i<ssCandidates.size();i++){
00652 int ii=ssCandidates[i];
00653 if (Debug()) LOG_INFO << mTofGeom->slatIdToDaq(slatIndex[ii]) << " ";
00654 if (vHitProfile[ii]>hitprof){
00655 hitprof = vHitProfile[ii];
00656 profileCandidates.clear();
00657 profileCandidates.push_back(ii);
00658 }else if (vHitProfile[ii]==hitprof)
00659 profileCandidates.push_back(ii);
00660 }
00661 if (profileCandidates.size()==1){
00662 thiscandidate = profileCandidates[0];
00663 int daqId = mTofGeom->slatIdToDaq(slatIndex[thiscandidate]);
00664 if (mHisto) hTofpSlatIdE4->Fill(daqId);
00665 if (Debug()) LOG_INFO << "candidate =" << daqId << endm;
00666 }
00667 else
00668 if (Debug()) LOG_INFO << "none" << endm;
00669 }
00670
00671
00672
00673 if (thiscandidate == -99 && Debug()){
00674 LOG_INFO << "E: ind=";
00675 for (unsigned int ii=0;ii<slatIndex.size();ii++)
00676 LOG_INFO << mTofGeom->slatIdToDaq(slatIndex[ii]) << " ";
00677 LOG_INFO << "\ttrkid:" << vTrackId[0] << " Unable to decide. ";
00678 LOG_INFO << "(hitprofs:";
00679 for (unsigned int ii=0;ii<slatIndex.size();ii++)
00680 LOG_INFO << vHitProfile[ii] << " ";
00681 LOG_INFO << " ss:";
00682 for (unsigned int ii=0;ii<slatIndex.size();ii++)
00683 LOG_INFO << vS[ii] << " ";
00684 LOG_INFO << ")" << endm;
00685 }
00686
00687 }
00688
00689
00690 if (thiscandidate>=0){
00691 slatHit.slatIndex = slatIndex[thiscandidate];
00692 slatHit.hitPosition = vPosition[thiscandidate];
00693 slatHit.layerHitPositions = vLayerHitPositions[thiscandidate];
00694
00695 slatHit.trackIdVec.push_back(vTrackId[thiscandidate]);
00696 slatHit.hitProfile = vHitProfile[thiscandidate];
00697 slatHit.s = vS[thiscandidate];
00698 slatHit.theta_xy = vTheta_xy[thiscandidate];
00699 slatHit.theta_zr = vTheta_zr[thiscandidate];
00700 slatHit.matchFlag = thisMatchFlag;
00701
00702 allMatchedSlatsVec.push_back(slatHit);
00703
00704 if (mHisto){
00705 int daqId = mTofGeom->slatIdToDaq(slatIndex[thiscandidate]);
00706 hTofpSlatIdE5->Fill(daqId);
00707 }
00708
00709
00710 if (Debug()) {
00711 LOG_INFO << "E: ind=" << mTofGeom->slatIdToDaq(slatHit.slatIndex)
00712 << "\thitprof="<< slatHit.hitProfile << "\ts="<<slatHit.s
00713 << "\tthxy="<<slatHit.theta_xy << "\tthzr="<<slatHit.theta_zr << "\ttrackid:"
00714 << vTrackId[thiscandidate] << endm;
00715 }
00716 }
00717 } else
00718 gMessMgr->Warning("","OS") << "E: no slats belong to this track ... should not happen!" << endm;
00719
00720 tempVec = erasedVec;
00721 }
00722 gMessMgr->Info("","OST") << "E: #before/#after: " << singleHitSlatsVec.size()
00723 << "/" << allMatchedSlatsVec.size() << endm;
00724
00725
00726
00727
00728
00729
00730
00731
00732 StTofSlatCollection *mSlatCollection = new StTofSlatCollection;
00733 int nValidSingleHitSlats(0), nValidSinglePrimHitSlats(0);
00734
00735 for (size_t ii=0; ii < allMatchedSlatsVec.size(); ii++){
00736 int daqId = mTofGeom->slatIdToDaq(allMatchedSlatsVec[ii].slatIndex);
00737 int jj = daqId-1;
00738
00739 if (allMatchedSlatsVec[ii].trackIdVec.size()!=1)
00740 gMessMgr->Warning("","OST") << "F: WHAT!?! mult.matched slat in single slat list " << daqId << endm;
00741
00742
00743 if (validTdc(mTofpTdc[jj])) nValidSingleHitSlats++;
00744
00745
00746 unsigned int trackNode = allMatchedSlatsVec[ii].trackIdVec[0];
00747 StTrack *theTrack = nodes[trackNode]->track(primary);
00748
00749
00750 if (validTofTrack(theTrack)){
00751 nValidSinglePrimHitSlats++;
00752 if (mHisto){
00753 float xhit = allMatchedSlatsVec[ii].hitPosition.x();
00754 float yhit = allMatchedSlatsVec[ii].hitPosition.y();
00755 float zhit = allMatchedSlatsVec[ii].hitPosition.z();
00756 float phihit = atan2(yhit,xhit);
00757 hTofpHitMap4->Fill(zhit,phihit);
00758 hTofpSlatIdF1->Fill(daqId);
00759 }
00760
00761
00762 int nHitsPerTrack = theTrack->topologyMap().numberOfHits(kTpcId);
00763 if(mHisto) hTofpNumberOfTrackHits->Fill(nHitsPerTrack);
00764
00765
00766 StTrackGeometry *theTrackGeometry = trackGeometry(theTrack);
00767
00768
00769 const StThreeVectorF momentum = theTrackGeometry->momentum();
00770 if (mHisto) hTofpPtTrack->Fill(momentum.perp());
00771
00772
00773 double pathLength = tofPathLength(&event->primaryVertex()->position(),
00774 &allMatchedSlatsVec[ii].hitPosition,
00775 theTrackGeometry->helix().curvature());
00776
00777
00778
00779
00780
00781 StThreeVectorD *pInnerLayer, *pOuterLayer;
00782 pInnerLayer = &(*(allMatchedSlatsVec[ii].layerHitPositions.begin()));
00783 pOuterLayer = &(*(allMatchedSlatsVec[ii].layerHitPositions.end() - 1));
00784
00785
00786 float dedx(0.), cherang(0);
00787 int dedx_np(0), cherang_nph(0);
00788 if (Debug()){
00789 StSPtrVecTrackPidTraits& traits = theTrack->pidTraits();
00790 for (unsigned int it=0;it<traits.size();it++){
00791 if (traits[it]->detector() == kTpcId){
00792 StDedxPidTraits *pid = dynamic_cast<StDedxPidTraits*>(traits[it]);
00793 if (pid && pid->method() ==kTruncatedMeanId){
00794 dedx = pid->mean();
00795 dedx_np = pid->numberOfPoints();
00796 }
00797 } else if (traits[it]->detector() == kRichId){
00798 StRichPidTraits *pid = dynamic_cast<StRichPidTraits*>(traits[it]);
00799 if (pid){
00800 StRichSpectra* pidinfo = pid->getRichSpectra();
00801 if (pidinfo && pidinfo->getCherenkovPhotons()>2){
00802 cherang = pidinfo->getCherenkovAngle();
00803 cherang_nph = pidinfo->getCherenkovPhotons();
00804 }
00805 }
00806 }
00807 }
00808 }
00809
00810
00811 float localHitPos = mTofGeom->slatHitPosition(&allMatchedSlatsVec[ii].hitPosition);
00812
00813
00814
00815 StTofSlat *tofSlat = new StTofSlat(jj+1,(int)mTofpAdc[jj],(int)mTofpTdc[jj],theTrack,
00816 localHitPos, allMatchedSlatsVec[ii].hitProfile,
00817 allMatchedSlatsVec[ii].matchFlag);
00818 tofSlat->setPosition(allMatchedSlatsVec[ii].hitPosition);
00819 mSlatCollection->push_back(tofSlat);
00820
00821
00822 if (Debug()){
00823 LOG_INFO << "F: ind=" << mTofGeom->slatIdToDaq(allMatchedSlatsVec[ii].slatIndex)
00824 << "\ttrackid:";
00825 idVectorIter ij=allMatchedSlatsVec[ii].trackIdVec.begin();
00826 while (ij != allMatchedSlatsVec[ii].trackIdVec.end()) { LOG_INFO << " " << *ij; ij++; }
00827 LOG_INFO << "\tR=" << 1/(theTrackGeometry->helix().curvature())
00828 << "\tpT=" << momentum.perp() << "\tp=" << momentum.mag()
00829 << "\thits="<< nHitsPerTrack << "\ts="<< pathLength
00830 << "\t#fitp=" <<theTrack->fitTraits().numberOfFitPoints(kTpcId)
00831 << "\t#trkp=" <<theTrack->detectorInfo()->numberOfPoints(kTpcId)
00832 << " \tdedx=" << dedx
00833 << " \tdca="<< theTrack->geometry()->helix().distance(event->primaryVertex()->position());
00834 if (cherang!=0) LOG_INFO << " \trich="<< cherang << " (" << cherang_nph << ")";
00835 LOG_INFO << endm;
00836 }
00837
00838 }
00839 }
00840
00841 storeMatchData(mSlatCollection,theTof);
00842 delete mSlatCollection;
00843
00844 gMessMgr->Info("","OST") << "F: #before/#after: " << allMatchedSlatsVec.size()
00845 << "/" <<nValidSinglePrimHitSlats << endm;
00846
00847
00848 if (theTof->dataPresent())
00849 gMessMgr->Info("- TofCollection: raw data container present","OST");
00850 if (theTof->slatsPresent()){
00851 gMessMgr->Info("- TofCollection: slat container present","OST");
00852 if (Debug()){
00853 StSPtrVecTofSlat& tmpSlatTofVec = theTof->tofSlats();
00854 for (size_t i = 0; i < tmpSlatTofVec.size(); i++) {
00855 StTofSlat* p = tmpSlatTofVec[i];
00856 LOG_INFO << p->slatIndex() << " " << p->adc() << " " << p->tdc()
00857 << " " << p->associatedTrack() << endm;
00858 }
00859 }
00860 }
00861
00862
00863
00864
00865
00866
00867 if (mHisto){
00868 hTofpNumberOfValidAdc->Fill(nAdcTofp);
00869 hTofpNumberOfValidTdc->Fill(nTdcTofp);
00870 hTofpNumberOfValidSlats->Fill(nAdcTdcTofp);
00871 hTofpNumberOfGlobalTracks->Fill(allSlatsHitVec.size());
00872 hTofpNumberOfHitSlats->Fill(allMatchedSlatsVec.size());
00873 hTofpNumberOfSingleHitTracks->Fill(nSingleHitSlats);
00874 hTofpNumberOfSingleValidHitTracks->Fill(nValidSingleHitSlats);
00875 }
00876 gMessMgr->Info("","OST") << "#(slat tracks): " << allSlatsHitVec.size()
00877 << " #(hit slats): " << allMatchedSlatsVec.size()
00878 << " #slats (valid tdc): " << nTdcTofp
00879 << " #(single hits): " << nSingleHitSlats
00880 << " #(single valid hits): " << nValidSingleHitSlats
00881 << " #(single prim valid hits): " << nValidSinglePrimHitSlats
00882 << endm;
00883
00884
00885
00886 if (doPrintMemoryInfo) {
00887 StMemoryInfo::instance()->snapshot();
00888 StMemoryInfo::instance()->print();
00889 }
00890 if (doPrintCpuInfo) {
00891 timer.stop();
00892 gMessMgr->Info("","OST") << "CPU time for StTofpMatchMaker::Make(): "
00893 << timer.elapsedTime() << " sec\n" << endm;
00894 }
00895
00896 gMessMgr->Info("StTofpMatchMaker -- bye-bye","OS");
00897 return kStOK;
00898 }
00899
00900
00901
00903 Int_t StTofpMatchMaker::storeMatchData(StTofSlatCollection *slatCollection,
00904 StTofCollection* tofCollection){
00905 if(!tofCollection){
00906 gMessMgr->Error("No TofCollection -- returning","OS");
00907 return kStErr;
00908 }
00909
00910 for (size_t j=0;j<slatCollection->size();j++){
00911 tofCollection->addSlat(slatCollection->getSlat(j));
00912 if (Debug())
00913 LOG_INFO << "storing " << j << " "
00914 << slatCollection->getSlat(j)->slatIndex() << endm;
00915 }
00916 return kStOK;
00917 }
00918
00919
00920
00922 Int_t StTofpMatchMaker::getTofData(StTofCollection* tofCollection){
00923 if (!tofCollection) return kStERR;
00924 StSPtrVecTofData &tofData = tofCollection->tofData();
00925
00926
00927 bool dataOK(true);
00928 gMessMgr->Info("TOF raw data consistency test","OS");
00929 for (int i=0;i<48;i++){
00930 if (tofData[i]->dataIndex() != mTofGeom->daqToSlatId(i)) {
00931 dataOK = false;
00932 gMessMgr->Warning("","OST") << "===>WARNING: " << tofData[i]->dataIndex() << " " << mTofGeom->daqToSlatId(i) << endm;
00933 }
00934
00935 }
00936
00937
00938 for (int i=0;i<NTOFP;i++){
00939 mTofpAdc[i] = tofData[i]->adc();
00940 mTofpTdc[i] = tofData[i]->tdc();
00941 }
00942
00943 if (mYear2){
00944
00945 float tmp = mTofpAdc[3];
00946 mTofpAdc[3] = mTofpAdc[2];
00947 mTofpAdc[2] = tmp;
00948 }
00949
00950 for (int i=0;i<NPVPD;i++){
00951 mPvpdAdc[i] = tofData[42+i]->adc();
00952 mPvpdTdc[i] = tofData[42+i]->tdc();
00953 if (mYear3||mYear4)
00954 mPvpdAdcLoRes[i] = tofData[54+i]->adc();
00955 }
00956
00957 if (!dataOK) return kStWarn;
00958
00959 return kStOK;
00960 }
00961
00962
00963
00965 void StTofpMatchMaker::bookHistograms(void){
00966
00967
00968
00969
00970 mHitPosHistNames = new TOrdCollection;
00971 hTofpHitMap1 = new TH2D("tofpHitMap1","valid hit positions", 500,-250,0, 120,-1.23, -1.08);
00972 mHitPosHistNames->AddLast(hTofpHitMap1);
00973 hTofpHitMap2 = new TH2D("tofpHitMap2","valid hit positions", 500,-250,0, 120,-1.23, -1.08);
00974 mHitPosHistNames->AddLast(hTofpHitMap2);
00975 hTofpHitMap3 = new TH2D("tofpHitMap3","valid hit positions", 500,-250,0, 120,-1.23, -1.08);
00976 mHitPosHistNames->AddLast(hTofpHitMap3);
00977 hTofpHitMap4 = new TH2D("tofpHitMap4","valid hit positions", 500,-250,0, 120,-1.23, -1.08);
00978 mHitPosHistNames->AddLast(hTofpHitMap4);
00979 hTofpSlatHitVecSize = new TH1D("SlatMult","Slat Mult per Track",10,0,10);
00980 mHitPosHistNames->AddLast(hTofpSlatHitVecSize);
00981 hTofpSlatIdA0 = new TH1D("tofpSlatIdA0","events per slat",41,0.5,41.5);
00982 mHitPosHistNames->AddLast(hTofpSlatIdA0);
00983 hTofpSlatIdA1 = new TH1D("tofpSlatIdA1","valid slat",41,0.5,41.5);
00984 mHitPosHistNames->AddLast(hTofpSlatIdA1);
00985 hTofpSlatIdB1 = new TH1D("tofpSlatIdB1","#tracks match valid slat",41,0.5,41.5);
00986 mHitPosHistNames->AddLast(hTofpSlatIdB1);
00987 hTofpSlatIdD1 = new TH1D("tofpSlatIdD1","track match per valid slat",41,0.5,41.5);
00988 mHitPosHistNames->AddLast(hTofpSlatIdD1);
00989 hTofpSlatIdD2 = new TH1D("tofpSlatIdD2","single track match per slat",41,0.5,41.5);
00990 mHitPosHistNames->AddLast(hTofpSlatIdD2);
00991 hTofpSlatIdE1 = new TH1D("tofpSlatIdE1","one slat for one track match",41,0.5,41.5);
00992 mHitPosHistNames->AddLast(hTofpSlatIdE1);
00993 hTofpSlatIdE2 = new TH1D("tofpSlatIdE2","recovered from hitprof-weight",41,0.5,41.5);
00994 mHitPosHistNames->AddLast(hTofpSlatIdE2);
00995 hTofpSlatIdE3 = new TH1D("tofpSlatIdE3","recovered from ss",41,0.5,41.5);
00996 mHitPosHistNames->AddLast(hTofpSlatIdE3);
00997 hTofpSlatIdE4 = new TH1D("tofpSlatIdE4","recovered from closest hitplane",41,0.5,41.5);
00998 mHitPosHistNames->AddLast(hTofpSlatIdE4);
00999 hTofpSlatIdE5 = new TH1D("tofpSlatIdE5","total recovered slat per track match",41,0.5,41.5);
01000 mHitPosHistNames->AddLast(hTofpSlatIdE5);
01001 hTofpSlatIdF1 = new TH1D("tofpSlatIdF1","primary track match per slat",41,0.5,41.5);
01002 mHitPosHistNames->AddLast(hTofpSlatIdF1);
01003
01004
01005
01006
01007 mTrackHistNames = new TOrdCollection;
01008 hTofpNumberOfTrackHits = new TH1D("tofpNumberOfTrackHits","numberOfTrackHits",80,0,80);
01009 mTrackHistNames->AddLast(hTofpNumberOfTrackHits);
01010 hTofpPtTrack = new TH1D("tofpPtTrack","ptTrack",250,0.,10);
01011 mTrackHistNames->AddLast(hTofpPtTrack);
01012 hTofpDCATrackprimVertex = new TH1D("tofpDCATrackprimVertex","DCA distribution",6000,-30.,30.);
01013 mTrackHistNames->AddLast(hTofpDCATrackprimVertex);
01014
01015 mOccupancyHistNames = new TOrdCollection;
01016 hTofpNumberOfValidAdc = new TH1D("tofpNumberOfValidTdc","numberOfValidTdc",41,0,41);
01017 mOccupancyHistNames->AddLast(hTofpNumberOfValidAdc);
01018 hTofpNumberOfValidTdc = new TH1D("tofpNumberOfValidAdc","numberOfValidAdc",41,0,41);
01019 mOccupancyHistNames->AddLast(hTofpNumberOfValidTdc);
01020 hTofpNumberOfValidSlats = new TH1D("tofpNumberOfValidSlats","numberOfValidSlats",41,0,41);
01021 mOccupancyHistNames->AddLast(hTofpNumberOfValidSlats);
01022 hTofpNumberOfGlobalTracks = new TH1D("tofpNumberOfGlobalTracks","numberOfGlobalTracks",50,0,50);
01023 mOccupancyHistNames->AddLast(hTofpNumberOfGlobalTracks);
01024 hTofpNumberOfHitSlats = new TH1D("tofpNumberOfHitSlats","numberOfHitSlats",50,0,50);
01025 mOccupancyHistNames->AddLast(hTofpNumberOfHitSlats);
01026 hTofpNumberOfSingleHitTracks = new TH1D("tofpNumberOfSingleHitTracks","numberOfSingleHitTracks",50,0,50);
01027 mOccupancyHistNames->AddLast(hTofpNumberOfSingleHitTracks);
01028 hTofpNumberOfSingleValidHitTracks = new TH1D("tofpNumberOfSingleValidTracks","numberOfSingleValidHitTracks",50,0,50);
01029 mOccupancyHistNames->AddLast(hTofpNumberOfSingleValidHitTracks);
01030
01031 mMatchHistNames = new TOrdCollection;
01032 char buf[20];
01033 for (int i=0;i<NTOFP;i++){
01034 sprintf(buf,"tofpSlathit_%d",i+1);
01035 hTofpMatchHit[i] = new TH2D(buf,buf,5,-2.5,2.5,5,-2.5,2.5);
01036 hTofpMatchHit[i]->SetXTitle("iEta");
01037 hTofpMatchHit[i]->SetYTitle("iPhi");
01038 mMatchHistNames->AddLast(hTofpMatchHit[i]);
01039 sprintf(buf,"tofpSlatnohit_%d",i+1);
01040 hTofpMatchNoHit[i] = new TH2D(buf,buf,5,-2.5,2.5,5,-2.5,2.5);
01041 hTofpMatchNoHit[i]->SetXTitle("iEta");
01042 hTofpMatchNoHit[i]->SetYTitle("iPhi");
01043 mMatchHistNames->AddLast(hTofpMatchNoHit[i]);
01044 }
01045
01046 return;
01047 }
01048
01049
01050
01052 void StTofpMatchMaker::writeHistogramsToFile(){
01053
01054 TFile *theHistoFile = new TFile(mHistoFileName.c_str(), "RECREATE");
01055 gMessMgr->Info("","OST") << "StTofpMatchMaker::writeHistogramsToFile()"
01056 << " histogram file " << mHistoFileName << endm;
01057
01058 theHistoFile->mkdir("hitpos","hit position histograms");
01059 theHistoFile->cd("hitpos");
01060 mHitPosHistNames->Write();
01061 delete mHitPosHistNames; mHitPosHistNames=0;
01062
01063 theHistoFile->mkdir("track","Tracking Plots");
01064 theHistoFile->cd("track");
01065 mTrackHistNames->Write();
01066 delete mTrackHistNames; mTrackHistNames=0;
01067
01068 theHistoFile->mkdir("occup","Occupancy Plots");
01069 theHistoFile->cd("occup");
01070 mOccupancyHistNames->Write();
01071 delete mOccupancyHistNames; mOccupancyHistNames=0;
01072
01073 theHistoFile->mkdir("match","Matching Plots");
01074 theHistoFile->cd("match");
01075 mMatchHistNames->Write();
01076 delete mMatchHistNames; mMatchHistNames=0;
01077
01078 theHistoFile->Write();
01079 theHistoFile->Close();
01080
01081 return;
01082 }
01083
01084
01085
01087 bool StTofpMatchMaker::strobeEvent(StSPtrVecTofData& tofData){
01088
01089
01090 int nStrobedPvpdTdcs=0;
01091 for(int i=0;i<NPVPD;i++)
01092 if((tofData[42+i]->tdc()>mStrobeTdcMin[i]) &&
01093 (tofData[42+i]->tdc()<mStrobeTdcMax[i]))
01094 nStrobedPvpdTdcs++;
01095
01096 if (nStrobedPvpdTdcs==NPVPD) return true;
01097
01098 return false;
01099 }
01100
01101
01102
01104 bool StTofpMatchMaker::validEvent(StEvent *event){
01105 mEventCounter++;
01106
01107 if (!event) return false;
01108
01109
01110 if (!event->primaryVertex()) return false;
01111 mAcceptedEventCounter++;
01112
01113
01114 if (!event->tofCollection()){
01115 gMessMgr->Info("TOF is not present","OST");
01116 return false;
01117 }
01118
01119
01120 if (!(event->tofCollection()->dataPresent())){
01121 gMessMgr->Info("TOF is present but no Raw Data","OST");
01122 if (!(event->tofCollection()->slatsPresent())){
01123 gMessMgr->Info(" and no Slat Data","OST");
01124 return false;
01125 }
01126 }
01127 mTofEventCounter++;
01128
01129
01130
01131 StSPtrVecTofData &tofData = event->tofCollection()->tofData();
01132 if (strobeEvent(tofData)){
01133 mTofStrobeEventCounter++;
01134 if (event->primaryVertex()) mAcceptAndStrobe++;
01135 gMessMgr->Info("strobe event","OTS");
01136 return false;
01137 }
01138 mAcceptAndBeam++;
01139
01140
01141 gMessMgr->Info("TOF present ... and valid beam event","OTS");
01142
01143 return true;
01144 }
01145
01146
01147
01149 bool StTofpMatchMaker::validTrack(StTrack *track){
01150
01151 if (!track) return false;
01152
01153
01154 if (track->flag()<=0) return false;
01155
01156
01157 if (track->topologyMap().numberOfHits(kTpcId) < mMinHitsPerTrack) return false;
01158
01159 if (track->fitTraits().numberOfFitPoints(kTpcId) < mMinFitPointsPerTrack) return false;
01160
01161 return true;
01162 }
01163
01164
01165
01167 bool StTofpMatchMaker::validTofTrack(StTrack *track){
01168
01169
01170
01171 if (!track) return false;
01172
01173
01174 if (!dynamic_cast<StPrimaryTrack*>(track)) return false;
01175
01176
01177 double DCA= track->impactParameter();
01178 int charge = track->geometry()->charge();
01179 if (mHisto) hTofpDCATrackprimVertex->Fill(DCA*charge);
01180 if (DCA > mMaxDCA) {
01181 gMessMgr->Info("","OST") << "dca>max:" << DCA<< endm;
01182 return false;
01183 }
01184
01185 return true;
01186 }
01187
01188
01189
01191 StTrackGeometry* StTofpMatchMaker::trackGeometry(StTrack* track){
01192
01193 if (!track) return 0;
01194 StTrackGeometry *thisTrackGeometry;
01195 if (mOuterTrackGeometry)
01196 thisTrackGeometry = track->outerGeometry();
01197 else
01198 thisTrackGeometry = track->geometry();
01199 return thisTrackGeometry;
01200 }
01201
01202
01203