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
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172
00173
00174
00175
00176
00177
00178
00179
00180
00181
00182
00183
00184
00185
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218
00219
00220
00221
00222
00223
00224
00225
00226
00227
00228
00229
00230
00231
00232
00233
00234
00235
00236
00237
00238
00239
00240
00241
00242
00243
00244
00245
00246
00247
00248
00249
00250
00251
00252
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270
00271
00272
00273
00274
00275
00276
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315
00316
00317
00318
00319
00320
00321
00322
00323
00324
00325
00326
00327
00328
00329
00330
00331
00332
00333
00334
00335
00336
00337
00338
00339
00340
00341
00342
00343
00344
00345
00346
00347
00348
00349
00350
00351
00352
00353
00354
00355
00356
00357
00358
00359
00360
00361
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371
00372
00373
00374
00375
00376
00377
00378
00379
00380
00381
00382
00383
00384
00385
00386
00387
00388
00389
00390
00391
00392
00393
00394
00395
00396
00397
00398
00399
00400
00401
00402
00403
00404
00405
00406
00407
00408
00409
00410
00411
00412
00413
00414
00415
00416
00417
00418
00419
00420
00421
00422
00423
00424
00425
00426
00427
00428
00429
00430
00431
00432
00433
00434
00435
00436
00437
00438
00439
00440
00441
00442
00443
00444
00445
00446
00447
00448
00449
00450
00451
00452
00453
00454
00455
00456
00457
00458
00459 #include "RVersion.h"
00460 #if ROOT_VERSION_CODE < 331013
00461 #include "TCL.h"
00462 #else
00463 #include "TCernLib.h"
00464 #endif
00465
00466 #include "Stiostream.h"
00467 #include <algorithm>
00468 #include <stdexcept>
00469 using namespace std;
00470
00471
00472 #include "StPhysicalHelix.hh"
00473 #include "StThreeVector.hh"
00474 #include "StThreeVectorF.hh"
00475 #include "PhysicalConstants.h"
00476 #include "SystemOfUnits.h"
00477 #include "StTrackDefinitions.h"
00478 #include "StTrackMethod.h"
00479 #include "StDedxMethod.h"
00480
00481
00482 #include "StPrimaryVertex.h"
00483 #include "StEventTypes.h"
00484 #include "StDetectorId.h"
00485 #include "StHelix.hh"
00486 #include "StDcaGeometry.h"
00487 #include "StHit.h"
00488
00489
00490 #include "StEventUtilities/StEventHelper.h"
00491 #include "StEventUtilities/StuFixTopoMap.cxx"
00492
00493 #include "Sti/StiTrackContainer.h"
00494 #include "Sti/StiKalmanTrack.h"
00495 #include "StDetectorDbMaker/StiKalmanTrackFitterParameters.h"
00497 #include "StiUtilities/StiDebug.h"
00498 #include "StiUtilities/StiPullEvent.h"
00499
00500
00501 #include "StiMaker/StiStEventFiller.h"
00502 #include "TMath.h"
00503 #include "StTrack2FastDetectorMatcher.h"
00504 #define NICE(angle) StiKalmanTrackNode::nice((angle))
00505 map<StiKalmanTrack*, StTrackNode*> StiStEventFiller::mTrkNodeMap;
00506 map<StTrackNode*, StiKalmanTrack*> StiStEventFiller::mNodeTrkMap;
00507 StiStEventFiller *StiStEventFiller::fgStiStEventFiller = 0;
00508
00509
00510 StiStEventFiller::StiStEventFiller() : mEvent(0), mTrackStore(0), mFastDetectorMatcher(0)
00511 {
00512 fgStiStEventFiller = this;
00513 mUseAux = 0;
00514 mAux = 0;
00515 mGloPri = 0;
00516 mPullEvent=0;
00517
00518 originD = new StThreeVectorD(0,0,0);
00519 physicalHelix = new StPhysicalHelixD(0.,0.,0.,*originD,-1);
00520
00521
00522
00523
00524
00525
00526
00527
00528
00529
00530
00531
00532
00533
00534
00535 unsigned short bit = 1 << tpcOther;
00536 mStiEncoded = kITKalmanFitId + bit;
00537 mFastDetectorMatcher = new StTrack2FastDetectorMatcher();
00538 }
00539
00540 StiStEventFiller::~StiStEventFiller()
00541 {
00542 delete physicalHelix; physicalHelix=0;
00543 delete originD; originD =0;
00544 SafeDelete(mFastDetectorMatcher);
00545 cout <<"StiStEventFiller::~StiStEventFiller()"<<endl;
00546 }
00547
00599
00600 void StiStEventFiller::fillEvent(StEvent* e, StiTrackContainer* t)
00601 {
00602 mFastDetectorMatcher->Clear();
00603 mFastDetectorMatcher->fillArrays(e);
00604
00605 mGloPri=0;
00606 if (e==0 || t==0)
00607 {
00608 cout <<"StiStEventFiller::fillEvent(). ERROR:\t"
00609 <<"Null StEvent ("<<e<<") || StiTrackContainer ("<<t<<"). Exit"<<endl;
00610 return;
00611 }
00612 mEvent = e;
00613 StEventHelper::Remove(mEvent,"StSPtrVecTrackNode");
00614 StEventHelper::Remove(mEvent,"StSPtrVecPrimaryVertex");
00615
00616 if (mUseAux) { mAux = new StiAux; e->Add(mAux);}
00617 mTrackStore = t;
00618 mTrkNodeMap.clear();
00619 mNodeTrkMap.clear();
00620 StSPtrVecTrackNode& trNodeVec = mEvent->trackNodes();
00621 StSPtrVecTrackDetectorInfo& detInfoVec = mEvent->trackDetectorInfo();
00622 int errorCount=0;
00623
00624 int fillTrackCount1=0;
00625 int fillTrackCount2=0;
00626 int fillTrackCountG=0;
00627 StErrorHelper errh;
00628 mTrackNumber=0;
00629 for (vector<StiTrack*>::iterator trackIt = mTrackStore->begin(); trackIt!=mTrackStore->end();++trackIt)
00630 {
00631 StiKalmanTrack* kTrack = static_cast<StiKalmanTrack*>(*trackIt);
00632 if (!accept(kTrack)) continue;
00633 mTrackNumber++;
00634 StTrackDetectorInfo* detInfo = new StTrackDetectorInfo;
00635 fillDetectorInfo(detInfo,kTrack,true);
00636
00637 StTrackNode* trackNode = new StTrackNode;
00638
00639 StGlobalTrack* gTrack = new StGlobalTrack;
00640 try
00641 {
00642 fillTrackCount1++;
00643 fillTrack(gTrack,kTrack,detInfo);
00644
00645 detInfoVec.push_back(detInfo);
00646
00647 gTrack->setKey(kTrack->getId());
00648 trackNode->addTrack(gTrack);
00649 trNodeVec.push_back(trackNode);
00650
00651
00652
00653
00654
00655 StTrackNode* node = trNodeVec.back();
00656 mTrkNodeMap.insert(pair<StiKalmanTrack*,StTrackNode*> (kTrack,node) );
00657 mNodeTrkMap.insert(pair<StTrackNode*,StiKalmanTrack*> (node,kTrack) );
00658 if (trackNode->entries(global)<1)
00659 cout << "StiStEventFiller::fillEvent() -E- Track Node has no entries!! -------------------------" << endl;
00660 int ibad = gTrack->bad();
00661 errh.Add(ibad);
00662 if (ibad) {
00663
00664
00665 }
00666 fillTrackCount2++;
00667 fillPulls(kTrack,0);
00668 if (gTrack->numberOfPossiblePoints()<10) continue;
00669 if (gTrack->geometry()->momentum().mag()<0.1) continue;
00670 fillTrackCountG++;
00671
00672 }
00673 catch (runtime_error & rte )
00674 {
00675 cout << "StiStEventFiller::fillEvent() -W- runtime-e filling track"<<rte.what() << endl;
00676 delete trackNode;
00677 delete detInfo;
00678 delete gTrack;
00679 }
00680 catch (...)
00681 {
00682 cout << "StiStEventFiller::fillEvent() -W- Unknown exception filling track."<<endl;
00683 delete trackNode;
00684 delete detInfo;
00685 delete gTrack;
00686 }
00687 }
00688 if (errorCount>4)
00689 cout << "There were "<<errorCount<<"runtime_error while filling StEvent"<<endl;
00690
00691 cout <<"StiStEventFiller::fillEvent() -I- Number of filled as global(1):"<< fillTrackCount1<<endl;
00692 cout <<"StiStEventFiller::fillEvent() -I- Number of filled as global(2):"<< fillTrackCount2<<endl;
00693 cout <<"StiStEventFiller::fillEvent() -I- Number of filled GOOD globals:"<< fillTrackCountG<<endl;
00694 errh.Print();
00695
00696 return;
00697 }
00698
00699 void StiStEventFiller::fillEventPrimaries()
00700 {
00701
00702 mGloPri=1;
00703 if (!mTrkNodeMap.size())
00704 {
00705 cout <<"StiStEventFiller::fillEventPrimaries(). ERROR:\t"
00706 << "Mapping between the StTrackNodes and the StiKalmanTracks is empty. Exit." << endl;
00707 return;
00708 }
00709
00710 StPrimaryVertex* vertex = 0;
00711 StSPtrVecTrackDetectorInfo& detInfoVec = mEvent->trackDetectorInfo();
00712 cout << "StiStEventFiller::fillEventPrimaries() -I- Tracks in container:" << mTrackStore->size() << endl;
00713 int mTrackN=0,mVertN=0;
00714 int noPipe=0;
00715 int ifcOK=0;
00716 int fillTrackCount1=0;
00717 int fillTrackCount2=0;
00718 int fillTrackCountG=0;
00719 StErrorHelper errh;
00720 int nTracks = mTrackStore->size();
00721 StiKalmanTrack *kTrack = 0;
00722 StPrimaryTrack *pTrack = 0;
00723 StGlobalTrack *gTrack = 0;
00724 StTrackNode *nTRack = 0;
00725 mTrackNumber=0;
00726 for (mTrackN=0; mTrackN<nTracks;++mTrackN) {
00727 kTrack = (StiKalmanTrack*)(*mTrackStore)[mTrackN];
00728 if (!accept(kTrack)) continue;
00729 map<StiKalmanTrack*, StTrackNode*>::iterator itKtrack = mTrkNodeMap.find(kTrack);
00730 if (itKtrack == mTrkNodeMap.end()) continue;
00731 mTrackNumber++;
00732
00733 nTRack = (*itKtrack).second;
00734 assert(nTRack->entries()<=10);
00735 assert(nTRack->entries(global));
00736
00737
00738
00739
00740
00741
00742 gTrack = static_cast<StGlobalTrack*>(nTRack->track(global));
00743 assert(gTrack->key()==kTrack->getId());
00744 float minDca = 1e10;
00745
00746 pTrack = 0;
00747 for (mVertN=0; (vertex = mEvent->primaryVertex(mVertN));mVertN++) {
00748 StThreeVectorD vertexPosition = vertex->position();
00749 double zPrim = vertexPosition.z();
00750
00751 float globalDca = impactParameter(gTrack,vertexPosition);
00752 if (fabs(minDca) > fabs(globalDca)) minDca = globalDca;
00753
00754 if (!kTrack->isPrimary()) continue;
00755 StiKalmanTrackNode *lastNode = kTrack->getLastNode();
00756 StiHit *pHit = lastNode->getHit();
00757 assert (pHit);
00758 if (fabs(pHit->z_g()-zPrim)>0.1) continue;
00759
00760 fillTrackCount1++;
00761
00762 StTrackDetectorInfo* detInfo = new StTrackDetectorInfo;
00763 fillDetectorInfo(detInfo,kTrack,false);
00764 fillPulls(kTrack,1);
00765 StPrimaryTrack* pTrack = new StPrimaryTrack;
00766 pTrack->setKey( gTrack->key());
00767 fillTrack(pTrack,kTrack, detInfo);
00768
00769 detInfoVec.push_back(detInfo);
00770
00771 nTRack->addTrack(pTrack);
00772 vertex->addDaughter(pTrack);
00773 fillTrackCount2++;
00774 int ibad = pTrack->bad();
00775 errh.Add(ibad);
00776 if (ibad) {
00777
00778
00779 }
00780 if (pTrack->numberOfPossiblePoints()<10) break;
00781 if (pTrack->geometry()->momentum().mag()<0.1) break;
00782 fillTrackCountG++;
00783 break;
00784 }
00785 kTrack->setDca(minDca);
00786 gTrack->setImpactParameter(minDca);
00787 if (pTrack) pTrack->setImpactParameter(minDca);
00788
00789 }
00790 for (mVertN=0; (vertex = mEvent->primaryVertex(mVertN));mVertN++) {vertex->setTrackNumbers();}
00791
00792 mTrkNodeMap.clear();
00793 cout <<"StiStEventFiller::fillEventPrimaries() -I- Primaries (1):"<< fillTrackCount1<< " (2):"<< fillTrackCount2<< " no pipe node:"<<noPipe<<" with IFC:"<< ifcOK<<endl;
00794 cout <<"StiStEventFiller::fillEventPrimaries() -I- GOOD:"<< fillTrackCountG <<endl;
00795 errh.Print();
00796 return;
00797 }
00798
00803
00804 void StiStEventFiller::fillDetectorInfo(StTrackDetectorInfo* detInfo, StiKalmanTrack* track, bool refCountIncr)
00805 {
00806
00807 int dets[kMaxDetectorId][3];
00808 track->getAllPointCount(dets,kMaxDetectorId-1);
00809 for (int i=1;i<kMaxDetectorId;i++) {
00810 if (!dets[i][1]) continue;
00811 detInfo->setNumberOfPoints(dets[i][1],static_cast<StDetectorId>(i));
00812 }
00813 StiKTNIterator tNode = track->rbegin();
00814 StiKTNIterator eNode = track->rend();
00815 StiKalmanTrackNode *lastNode=0,*fistNode=0;
00816 for (;tNode!=eNode;++tNode)
00817 {
00818 StiKalmanTrackNode *node = &(*tNode);
00819 if(!node->isValid()) continue;
00820
00821 StiHit *stiHit = node->getHit();
00822 if (!stiHit) continue;
00823
00824 if (node->getChi2()>1000) continue;
00825 if (!node->isFitted()) continue;
00826
00827 const StiDetector *detector = node->getDetector();
00828 assert(detector == stiHit->detector());
00829 assert(!detector || stiHit->timesUsed());
00830 if (!fistNode) fistNode = node;
00831 lastNode = node;
00832 StHit *hh = (StHit*)stiHit->stHit();
00833 if (!detector) continue;
00834 if (!hh) continue;
00835 assert(detector->getGroupId()==hh->detector());
00836
00837 FillStHitErr(hh,node);
00838
00839 detInfo->addHit(hh,refCountIncr);
00840 if (!refCountIncr) continue;
00841 hh->setFitFlag(1);
00842
00843 fillResHack(hh,stiHit,node);
00844 }
00845 assert(lastNode && fistNode && (lastNode != fistNode));
00846
00847 StThreeVectorF posL(lastNode->x_g(),lastNode->y_g(),lastNode->z_g());
00848 detInfo->setLastPoint (posL);
00849 StThreeVectorF posF(fistNode->x_g(),fistNode->y_g(),fistNode->z_g());
00850 detInfo->setFirstPoint(posF);
00851
00852
00853
00854 }
00855
00856
00857 void StiStEventFiller::fillGeometry(StTrack* gTrack, StiKalmanTrack* track, bool outer)
00858 {
00859
00860 assert(gTrack);
00861 assert(track) ;
00862
00863 StiKalmanTrackNode* node = track->getInnOutMostNode(outer,3);
00864 StiHit *ihit = node->getHit();
00865 StThreeVectorF origin(node->x_g(),node->y_g(),node->z_g());
00866 StThreeVectorF hitpos(ihit->x_g(),ihit->y_g(),ihit->z_g());
00867
00868 double dif = (hitpos-origin).mag();
00869
00870 if (dif>3.) {
00871 dif = node->z_g()-ihit->z_g();
00872 double nowChi2 = node->evaluateChi2(ihit);
00873 printf("***Track(%d) DIFF TOO BIG %g chi2 = %g %g\n",track->getId(),dif,node->getChi2(),nowChi2);
00874 printf("H=%g %g %g N =%g %g %g\n",ihit->x() ,ihit->y() ,ihit->z()
00875 ,node->getX(),node->getY(),node->getZ());
00876 const StMeasuredPoint *mp = ihit->stHit();
00877 printf("H=%g %g %g N =%g %g %g\n",mp->position().x(),mp->position().y(),mp->position().z()
00878 ,origin.x(),origin.y(),origin.z());
00879
00880 assert(fabs(dif)<50.);
00881 }
00882
00883
00884
00885
00886 int ibad = origin.bad();
00887 if (ibad) {
00888 cout << "StiStEventFiller::fillGeometry() Encountered non-finite numbers!!!! Bail out completely!!! " << endl;
00889 cout << "StThreeVectorF::bad() = " << ibad << endl;
00890 cout << "Last node had:" << endl;
00891 cout << "Ref Position " << node->getRefPosition() << endl;
00892 cout << "node->getY() " << node->getY() << endl;
00893 cout << "node->getZ() " << node->getZ() << endl;
00894 cout << "Ref Angle " << node->getAlpha() << endl;
00895 cout << "origin " << origin << endl;
00896 cout << "curvature " << node->getCurvature() << endl;
00897 abort();
00898 }
00899 StTrackGeometry* geometry =new StHelixModel(short(track->getCharge()),
00900 node->getPsi(),
00901 fabs(node->getCurvature()),
00902 node->getDipAngle(),
00903 origin,
00904 node->getGlobalMomentumF(),
00905 node->getHelicity());
00906
00907 if (outer)
00908 gTrack->setOuterGeometry(geometry);
00909 else
00910 gTrack->setGeometry(geometry);
00911
00912
00913 return;
00914 }
00915
00916
00917
00918
00919
00920
00921
00922
00923
00924
00925
00926
00927
00928
00929 void StiStEventFiller::fillFitTraits(StTrack* gTrack, StiKalmanTrack* track){
00930
00931
00932 unsigned short geantIdPidHyp = 9999;
00933
00934 geantIdPidHyp = 9;
00935
00936
00937 StiKalmanTrackNode* node = track->getInnerMostHitNode(3);
00938 float x[6],covMFloat[15];
00939 node->getGlobalTpt(x,covMFloat);
00940 float chi2[2];
00941
00942 chi2[0] = track->getChi2();
00943 chi2[1] = -999;
00944
00945
00946 if (gTrack->type()==primary) {
00947 assert(node->getDetector()==0);
00948 chi2[1]=node->getChi2();
00949 }
00950
00951
00952
00953
00954 StTrackFitTraits fitTraits(geantIdPidHyp,0,chi2,covMFloat);
00955
00956
00957
00958 int dets[kMaxDetectorId][3];
00959 track->getAllPointCount(dets,kMaxDetectorId-1);
00960
00961 for (int i=1;i<kMaxDetectorId;i++) {
00962 if (!dets[i][2]) continue;
00963 fitTraits.setNumberOfFitPoints((unsigned char)dets[i][2],(StDetectorId)i);
00964 }
00965 if (gTrack->type()==primary) {
00966 fitTraits.setPrimaryVertexUsedInFit(true);
00967 }
00968 gTrack->setFitTraits(fitTraits);
00969 return;
00970 }
00971
01003
01004 void StiStEventFiller::fillFlags(StTrack* gTrack) {
01005 Int_t flag = 0;
01006 if (gTrack->type()==global) {
01007 flag = 101;
01008 }
01009 else if (gTrack->type()==primary) {
01010 flag = 301;
01011 }
01012 StTrackFitTraits& fitTrait = gTrack->fitTraits();
01013
01014 int svtFitPoints = fitTrait.numberOfFitPoints(kSvtId);
01015 int ssdFitPoints = fitTrait.numberOfFitPoints(kSsdId);
01016 int pxlFitPoints = fitTrait.numberOfFitPoints(kPxlId);
01017 int istFitPoints = fitTrait.numberOfFitPoints(kIstId);
01018
01022
01023
01024
01025
01026
01027
01028 if (svtFitPoints+ssdFitPoints+pxlFitPoints+istFitPoints>0) {
01029 if (gTrack->type()==global) {
01030 flag = 501;
01031 }
01032 else if (gTrack->type()==primary) {
01033 flag = 601;
01034 }
01035 }
01036 const StTrackDetectorInfo *dinfo = gTrack->detectorInfo();
01037 if (dinfo) {
01038 Int_t NoTpcFitPoints = dinfo->numberOfPoints(kTpcId);
01039 Int_t NoFtpcWestId = dinfo->numberOfPoints(kFtpcWestId);
01040 Int_t NoFtpcEastId = dinfo->numberOfPoints(kFtpcEastId);
01041
01042
01043 if (NoTpcFitPoints >= 11) {
01044 const StPtrVecHit& hits = dinfo->hits(kTpcId);
01045 Int_t Nhits = hits.size();
01046 Int_t NoWrongSignZ = 0;
01047 Int_t NoPositiveSignZ = 0;
01048 Int_t NoNegativeSignZ = 0;
01049 Int_t NoPromptHits = 0;
01050 for (Int_t i = 0; i < Nhits; i++) {
01051 const StTpcHit *hit = (StTpcHit *) hits[i];
01052 if ((hit->position().z() < -1.0 && hit->sector() <= 12) ||
01053 (hit->position().z() > 1.0 && hit->sector() > 12)) NoWrongSignZ++;
01054 else {
01055 if (hit->position().z() < -1.0) NoNegativeSignZ++;
01056 if (hit->position().z() > 1.0) NoPositiveSignZ++;
01057 }
01058 if (TMath::Abs(209.4 - TMath::Abs(hit->position().z())) < 3.0) NoPromptHits++;
01059 }
01060 if (NoWrongSignZ >= 2) gTrack->setPostCrossingTrack();
01061 else {
01062 if (NoPromptHits == 1) gTrack->setPromptTrack();
01063 if (NoPositiveSignZ >= 2 && NoNegativeSignZ >=2) gTrack->setMembraneCrossingTrack();
01064 }
01065 }
01066 if (NoTpcFitPoints < 11 && NoFtpcWestId < 5 && NoFtpcEastId < 5) {
01067
01068
01069 gTrack->setRejected();
01070 flag = - ((flag/100)*100 + 2);
01071 if (gTrack->geometry()) {
01072 const StThreeVectorF &momentum = gTrack->geometry()->momentum();
01073 if (momentum.pseudoRapidity() > 0.5) {
01074 const StTrackDetectorInfo *dinfo = gTrack->detectorInfo();
01075 const StPtrVecHit& hits = dinfo->hits();
01076 Int_t Nhits = hits.size();
01077 Bool_t ShortTrack2EMC = kFALSE;
01078 for (Int_t i = 0; i < Nhits; i++) {
01079 const StHit *hit = hits[i];
01080 if (hit->position().z() > 150.0) {
01081 ShortTrack2EMC = kTRUE;
01082 break;
01083 }
01084 }
01085 if (ShortTrack2EMC) {
01086 gTrack->setShortTrack2EMC();
01087 flag = (TMath::Abs(flag)/100)*100+11; ;
01088 }
01089 }
01090 }
01091 }
01092 }
01093
01094 gTrack->setFlag( flag);
01095
01096 StPhysicalHelixD hlx = gTrack->outerGeometry()->helix();
01097 TrackData t;
01098 mFastDetectorMatcher->matchTrack2FastDetectors(&hlx,&t);
01099 if (t.btofBin > 0) {
01100 if (t.mBtof > 0) gTrack->setToFMatched();
01101 else gTrack->setToFNotMatched();
01102 }
01103 if (t.ctbBin > 0) {
01104 if (t.mCtb > 0) gTrack->setCtbMatched();
01105 else gTrack->setCtbNotMatched();
01106 }
01107 if (t.bemcBin > 0) {
01108 if (t.mBemc > 0) gTrack->setBemcMatched();
01109 else gTrack->setBemcNotMatched();
01110 }
01111 if (t.eemcBin > 0) {
01112 if (t.mEemc > 0) gTrack->setEemcMatched();
01113 else gTrack->setEemcNotMatched();
01114 }
01115 }
01116
01117 void StiStEventFiller::fillTrack(StTrack* gTrack, StiKalmanTrack* track,StTrackDetectorInfo* detInfo )
01118 {
01119
01120
01121
01122
01123
01124
01125
01126
01127
01128
01129
01130
01131
01132
01133
01134
01135 gTrack->setEncodedMethod(mStiEncoded);
01136 double tlen = track->getTrackLength();
01137 assert(tlen >0.0 && tlen<1000.);
01138 gTrack->setLength(tlen);
01139 gTrack->setSeedQuality(track->getSeedHitCount());
01140
01141
01142
01143 int dets[kMaxDetectorId][3];
01144 track->getAllPointCount(dets,kMaxDetectorId-1);
01145 for (int i=1;i<kMaxDetectorId;i++) {
01146 if(!dets[i][0]) continue;
01147 gTrack->setNumberOfPossiblePoints((unsigned char)dets[i][0],(StDetectorId)i);
01148 }
01149 fillGeometry(gTrack, track, false);
01150 fillGeometry(gTrack, track, true );
01151 fillFitTraits(gTrack, track);
01152 gTrack->setDetectorInfo(detInfo);
01153 StuFixTopoMap(gTrack);
01154 fillFlags(gTrack);
01155 if (!track->isPrimary()) fillDca(gTrack,track);
01156 return;
01157 }
01158
01159 bool StiStEventFiller::accept(StiKalmanTrack* track)
01160 {
01161
01162
01163 int nFittedPoints = track->getFitPointCount(0);
01164 if (nFittedPoints < 5 ) return 0;
01165 #if 0
01166 if (nFittedPoints < 10 && nFittedPoints*2 < nPossiblePoints)return 0;
01167 if(track->getPt()<=0.1) return 0;
01168 #endif
01169 if(track->getTrackLength()<=0) return 0;
01170
01171
01172
01173 return 1;
01174 }
01175
01176 double StiStEventFiller::impactParameter(StiKalmanTrack* track
01177 ,StThreeVectorD &vertexPosition)
01178 {
01179 StiKalmanTrackNode* node;
01180
01181 node = track->getInnerMostNode(2);
01182
01183
01184 originD->setX(node->x_g());
01185 originD->setY(node->y_g());
01186 originD->setZ(node->z_g());
01187
01188
01189 physicalHelix->setParameters(fabs(node->getCurvature()),
01190 node->getDipAngle(),
01191 node->getPhase(),
01192 *originD,
01193 node->getHelicity());
01194
01195
01196
01197 return physicalHelix->distance(vertexPosition);
01198 }
01199
01200 double StiStEventFiller::impactParameter(StTrack* track, StThreeVectorD &vertex)
01201 {
01202 StPhysicalHelixD helix = track->geometry()->helix();
01203
01204
01205 return helix.distance(vertex);
01206 }
01207
01208 void StiStEventFiller::fillResHack(StHit *hh,const StiHit *stiHit, const StiKalmanTrackNode *node)
01209 {
01210
01211 if (!mAux) return;
01212 StiAux_t aux;
01213
01214 aux.xnl[0] = node->getX();
01215 aux.xnl[1] = node->getY();
01216 aux.xnl[2] = node->getZ();
01217
01218 aux.xhl[0] = stiHit->x();
01219 aux.xhl[1] = stiHit->y();
01220 aux.xhl[2] = stiHit->z();
01221
01222 aux.ca = node->getEta();
01223 aux.rho = node->getCurvature();
01224 aux.nYY = node->getCyy();
01225 aux.nZZ = node->getCzz();
01226 aux.hYY = node->getEyy();
01227 aux.hZZ = node->getEzz();
01228
01229 aux.unl[0] = node->getX();
01230 aux.unl[1] = node->unTouched().mPar[0];
01231 aux.unl[2] = node->unTouched().mPar[1];
01232 aux.uYY = sqrt(node->unTouched().mErr[0]);
01233 aux.uZZ = sqrt(node->unTouched().mErr[2]);
01234
01235
01236
01237 aux.xng[0] = node->x_g();
01238 aux.xng[1] = node->y_g();
01239 aux.xng[2] = node->z_g();
01240
01241 aux.xhg[0] = stiHit->x_g();
01242 aux.xhg[1] = stiHit->y_g();
01243 aux.xhg[2] = stiHit->z_g();
01244 aux.psi = node->getPsi();
01245 aux.dip = node->getDipAngle();
01246
01247 double chi2 = node->getChi2();if (chi2>1000) chi2=1000;
01248 aux.chi2 = chi2;
01249 int id = mAux->AddIt(&aux);
01250 assert(id);
01251 hh->setId(id);
01252 assert(hh->id());
01253
01254
01255
01256 }
01257
01258 void StiStEventFiller::fillDca(StTrack* stTrack, StiKalmanTrack* track)
01259 {
01260 StGlobalTrack *gTrack = dynamic_cast<StGlobalTrack*>(stTrack);
01261 assert(gTrack);
01262
01263 StiKalmanTrackNode *tNode = track->getInnerMostNode();
01264 if (!tNode->isDca()) return;
01265 const StiNodePars &pars = tNode->fitPars();
01266 const StiNodeErrs &errs = tNode->fitErrs();
01267 float alfa = tNode->getAlpha();
01268 Float_t setp[7] = {pars.y(), pars.z(), pars.eta(), pars.ptin(), pars.tanl(), pars.curv(), pars.hz()};
01269 setp[2]+= alfa;
01270 Float_t sete[15];
01271 for (int i=1,li=1,jj=0;i< kNPars;li+=++i) {
01272 for (int j=1;j<=i;j++) {sete[jj++]=errs.A[li+j];}}
01273 StDcaGeometry *dca = new StDcaGeometry;
01274 gTrack->setDcaGeometry(dca);
01275 dca->set(setp,sete);
01276
01277 }
01278
01279 void StiStEventFiller::FillStHitErr(StHit *hh,const StiKalmanTrackNode *node)
01280 {
01281 double stiErr[6],stErr[6];
01282 memcpy(stiErr,node->hitErrs(),sizeof(stiErr));
01283 double alfa = node->getAlpha();
01284 double c = cos(alfa);
01285 double s = sin(alfa);
01286 double T[3][3]={{c,-s, 0}
01287 ,{s, c, 0}
01288 ,{0, 0, 1}};
01289
01290 TCL::trasat(T[0],stiErr,stErr,3,3);
01291 StThreeVectorF f3(sqrt(stErr[0]),sqrt(stErr[2]),sqrt(stErr[5]));
01292 hh->setPositionError(f3);
01293 }
01294
01295 void StiStEventFiller::fillPulls(StiKalmanTrack* track, int gloPri)
01296 {
01297
01298 if (!mPullEvent) return;
01299 if (gloPri && track->isPrimary()!=1) return;
01300 int dets[kMaxDetectorId][3];
01301 track->getAllPointCount(dets,kMaxDetectorId-1);
01302 StiPullTrk aux;
01303 aux.mVertex = (unsigned char)track->isPrimary();
01304 aux.mTrackNumber=track->getId();
01305 aux.nAllHits = dets[0][2];
01306 aux.nTpcHits = dets[kTpcId][2];
01307 aux.nSvtHits = dets[kSvtId][2];
01308 aux.nSsdHits = dets[kSsdId][2];
01309 aux.nPxlHits = dets[kPxlId][2];
01310 aux.nIstHits = dets[kIstId][2];
01311 aux.mL = (unsigned char)track->getTrackLength();
01312 aux.mChi2 = track->getChi2();
01313 aux.mCurv = track->getCurvature();
01314 aux.mPt = track->getPt();
01315 aux.mPsi = track->getPhi();
01316 aux.mDip = atan(track->getTanL());
01317 StThreeVectorD v3 = track->getPoint();
01318 aux.mRxy = v3.perp();
01319 aux.mPhi = v3.phi();
01320 aux.mZ = v3.z();
01321 mPullEvent->Add(aux,gloPri);
01322
01323
01324 StiKTNIterator tNode = track->rbegin();
01325 StiKTNIterator eNode = track->rend();
01326 for (;tNode!=eNode;++tNode)
01327 {
01328 StiKalmanTrackNode *node = &(*tNode);
01329 if(!node->isValid()) continue;
01330
01331 StiHit *stiHit = node->getHit();
01332 if (!stiHit) continue;
01333
01334 if (node->getChi2()>1000) continue;
01335 if (!node->isFitted()) continue;
01336
01337 const StiDetector *detector = node->getDetector();
01338 assert(detector == stiHit->detector());
01339 assert(!detector || stiHit->timesUsed());
01340 StHit *hh = (StHit*)stiHit->stHit();
01341 fillPulls(hh,stiHit,node,track,dets,gloPri);
01342 if (gloPri) continue;
01343 fillPulls(hh,stiHit,node,track,dets,2);
01344 }
01345 }
01346
01347 void StiStEventFiller::fillPulls(StHit *stHit,const StiHit *stiHit
01348 ,const StiKalmanTrackNode *node
01349 ,const StiKalmanTrack *track
01350 ,int dets[1][3],int gloPriRnd)
01351 {
01352 double x,y,z,r,xp,yp,zp,rp;
01353 float untErrs[3];
01354
01355 const StiNodeInf *inf = 0;
01356 if (gloPriRnd==2) {
01357 inf = node->getInfo();
01358 if (!inf) return;
01359 }
01360 double timeFlight = node->getTime();
01361 const StiNodeErrs &mFE = (inf)? inf->mPE : node->fitErrs();
01362 const StiNodePars &mFP = (inf)? inf->mPP : node->fitPars();
01363 StiHitErrs mHrr;
01364 memcpy(mHrr.A, (inf)? inf->mHrr.A : node->hitErrs(),sizeof(StiHitErrs));
01365
01366 StiPullHit aux;
01367
01368
01369 aux.mVertex = (unsigned char)track->isPrimary();
01370 aux.nHitCand = node->getHitCand();
01371 aux.iHitCand = node->getIHitCand();
01372 if (!aux.nHitCand) aux.nHitCand=1;
01373 aux.lXHit = stiHit->x();
01374 aux.lYHit = stiHit->y(timeFlight);
01375 aux.lZHit = stiHit->z(timeFlight);
01376 aux.lYHitErr = sqrt(mHrr.hYY);
01377 aux.lZHitErr = sqrt(mHrr.hZZ);
01378 aux.lHitEmx[0] = mHrr.hYY;
01379 aux.lHitEmx[1] = mHrr.hZY;
01380 aux.lHitEmx[2] = mHrr.hZZ;
01381
01382
01383 aux.lXFit = mFP.x();
01384 aux.lYFit = mFP.y();
01385 aux.lZFit = mFP.z();
01386 aux.lYFitErr = sqrt(mFE._cYY);
01387 aux.lZFitErr = sqrt(mFE._cZZ);
01388 aux.lFitEmx[0] = mFE._cYY;
01389 aux.lFitEmx[1] = mFE._cZY;
01390 aux.lFitEmx[2] = mFE._cZZ;
01391
01392
01393
01394 xp = aux.lXHit;
01395 yp = (inf)? mFP.y(): (double)node->unTouched().mPar[0];
01396 zp = (inf)? mFP.z(): (double)node->unTouched().mPar[1];
01397 aux.lYPul = aux.lYHit-yp;
01398 aux.lZPul = aux.lZHit-zp;
01399 if (fabs(aux.lYPul)>10) StiDebug::Break(-1);
01400 if (fabs(aux.lZPul)>10) StiDebug::Break(-1);
01401 if (!inf) {TCL::ucopy(node->unTouched().mErr,untErrs,3);}
01402 else {TCL::ucopy(aux.lFitEmx ,untErrs,3);}
01403
01404 assert(untErrs[0]>=0);
01405 assert(untErrs[2]>=0);
01406
01407 TCL::vadd(untErrs,aux.lHitEmx,aux.lPulEmx,3);
01408 aux.lYPulErr = sqrt(aux.lPulEmx[0]);
01409 aux.lZPulErr = sqrt(aux.lPulEmx[2]);
01410
01411 aux.lPsi = mFP.eta();
01412 aux.lDip = atan(mFP.tanl());
01413
01414
01415 double alfa = node->getAlpha();
01416 float F[2][2];
01417
01418
01419 x = stiHit->x(); y = stiHit->y(timeFlight); z = stiHit->z(timeFlight);
01420 r = sqrt(x*x+y*y);
01421
01422 aux.gRHit = r;
01423 aux.gPHit = atan2(stiHit->y_g(),stiHit->x_g());
01424 aux.gZHit = stiHit->z_g();
01425 memset(F[0],0,sizeof(F));
01426 F[0][0]= x/(r*r);
01427 F[1][1]= 1;
01428 TCL::trasat(F[0],aux.lHitEmx,aux.gHitEmx,2,2);
01429 aux.gPHitErr = sqrt(aux.gHitEmx[0]);
01430 aux.gZHitErr = sqrt(aux.gHitEmx[2]);
01431
01432
01433
01434 x = mFP.x(); y = mFP.y();z = mFP.z();
01435 r = sqrt(x*x+y*y);
01436 aux.gRFit = r;
01437 aux.gPFit = NICE(atan2(y,x)+alfa);
01438 aux.gZFit = z;
01439
01440 memset(F[0],0,sizeof(F));
01441 F[0][0]= x/(r*r);
01442 F[1][1]= 1;
01443 TCL::trasat(F[0],aux.lFitEmx,aux.gFitEmx,2,2);
01444 aux.gPFitErr = sqrt(aux.gFitEmx[0]);
01445 aux.gZFitErr = sqrt(aux.gFitEmx[2]);
01446
01447
01448 rp = sqrt(xp*xp+yp*yp);
01449 aux.gPPul = ((aux.gPHit-alfa)-atan2(yp,xp))*rp;
01450 aux.gZPul = aux.lZHit-zp;
01451 memset(F[0],0,sizeof(F));
01452 F[0][0]= xp/(rp*rp);
01453 F[1][1]= 1;
01454 TCL::trasat(F[0],untErrs,aux.gPulEmx,2,2);
01455 TCL::vadd(aux.gHitEmx,aux.gPulEmx,aux.gPulEmx,3);
01456
01457 aux.gPulEmx[0]*= rp*rp;
01458 aux.gPulEmx[1]*= rp;
01459 aux.gPPulErr = sqrt(aux.gPulEmx[0]);
01460 aux.gZPulErr = sqrt(aux.gPulEmx[2]);
01461
01462 aux.gPsi = node->getPsi();
01463 aux.gDip = node->getDipAngle();
01464
01465
01466 aux.mCurv = mFP.curv();
01467 aux.mPt = fabs(1./mFP.ptin());
01468 aux.mCharge = stHit->charge();
01469 aux.mChi2 = node->getChi2();
01470 aux.mNormalRefAngle = alfa;
01471 aux.mHardwarePosition=0;
01472 aux.mDetector=0;
01473 aux.mTrackNumber=track->getId();
01474 aux.nAllHits = dets[0][2];
01475 aux.nTpcHits = dets[kTpcId][2];
01476 aux.nSvtHits = dets[kSvtId][2];
01477 aux.nSsdHits = dets[kSsdId][2];
01478 aux.nPxlHits = dets[kPxlId][2];
01479 aux.nIstHits = dets[kIstId][2];
01480 const StiDetector *stiDet = stiHit->detector();
01481 if (stiDet) {
01482 aux.mHardwarePosition=stHit->hardwarePosition();
01483 aux.mDetector=stHit->detector();
01484 const StiPlacement *place = stiDet->getPlacement();
01485 aux.mNormalRadius = place->getNormalRadius();
01486 aux.mNormalYOffset = place->getNormalYoffset();
01487 aux.mZCenter = 0;
01488 }
01489 mPullEvent->Add(aux,gloPriRnd);
01490
01491 }