00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 #include "StV0FinderMaker.h"
00013 #include "StMessMgr.h"
00014 #include "StEvent/StEventTypes.h"
00015 #include "TMath.h"
00016 #include "TVector2.h"
00017 #include "tables/St_V0FinderParameters_Table.h"
00018
00021 #include "StEvent/StTrack.h"
00022 #include "StMuDSTMaker/COMMON/StMuTypes.hh"
00024
00025 #include "math_constants.h"
00026 #include "phys_constants.h"
00027 #include "SystemOfUnits.h"
00028
00029
00030
00031
00032
00033
00034 static const int BLOCK=1024;
00035
00036
00037 StV0FinderMaker* StV0FinderMaker::mInstance = 0;
00038
00039 ClassImp(StV0FinderMaker)
00040
00041
00042 StV0FinderMaker::StV0FinderMaker(const char *name):StMaker(name),
00043 v0pars(0),pars(0),pars2(0),event(0),v0Vertex(0),
00044 prepared(kFALSE),useExistingV0s(kFALSE),dontZapV0s(kFALSE),
00045 useTracker(kTrackerUseBOTH),useSVT(kNoSVT),useEventModel(kUseStEvent),
00046 useV0Language(kV0LanguageUseCpp),useXiLanguage(kXiLanguageUseCppOnCppV0),
00047 useLanguage(kLanguageUseRun),useLikesign(kLikesignUseStandard),
00048 useRotating(kRotatingUseStandard)
00049 {
00050
00051 trks = 0;
00052 det_id_v0 = 0;
00053 ITTFflag = 0;
00054 TPTflag = 0;
00055 mainv.setX(0.);
00056 mainv.setY(0.);
00057 mainv.setZ(0.);
00058
00059 trkcnt = 0;
00060 trkmax = BLOCK;
00061 trkNodeRatio = 1.;
00062 trkNodeRatioCnt = 0.;
00063
00064
00065 if (mInstance != 0)
00066 gMessMgr->Warning() << "(" << name <<
00067 ") : MORE THAN ONE INSTANCE!" << endm;
00068 else mInstance = this;
00069
00070 Bfield = 1.e-10;
00071
00072
00073 }
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085 StV0FinderMaker::~StV0FinderMaker() {
00086 }
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100 void StV0FinderMaker::GetPars()
00101 {
00102 TDataSet* dbDataSet = GetDataBase("Calibrations/tracker");
00103 if (!dbDataSet) {
00104 gMessMgr->Error(
00105 "Init() : could not find Calibrations/tracker database.");
00106 return;
00107 }
00108 v0pars = (St_V0FinderParameters*) (dbDataSet->FindObject("V0FinderParameters"));
00109 if (!v0pars) {
00110 gMessMgr->Error(
00111 "Init() : could not find V0FinderParameters in database.");
00112 return;
00113 }
00115 }
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130 Int_t StV0FinderMaker::Init()
00131 {bool a,b,c;
00132
00133 if ((useTracker!=kTrackerUseTPT) && (useTracker!=kTrackerUseITTF) && (useTracker!=kTrackerUseBOTH))
00134 {gMessMgr->Error("Init() : wrong TrackerUsage parameter set.");
00135 return kStErr;
00136 }
00137 if ((useSVT!=kNoSVT) && (useSVT!=kUseSVT))
00138 {gMessMgr->Error("Init() : wrong SVTUsage parameter set.");
00139 return kStErr;
00140 }
00141 if ((useEventModel!=kUseStEvent) && (useEventModel!=kUseMuDst))
00142 {gMessMgr->Error("Init() : wrong EventModelUsage parameter set.");
00143 return kStErr;
00144 }
00145 if ((useLikesign!=kLikesignUseStandard) && (useLikesign!=kLikesignUseLikesign))
00146 {gMessMgr->Error("Init() : wrong LikesignUsage parameter set.");
00147 return kStErr;
00148 }
00149 if ((useRotating!=kRotatingUseStandard) && (useRotating!=kRotatingUseRotating) && (useRotating!=kRotatingUseSymmetry) && (useRotating!=kRotatingUseRotatingAndSymmetry))
00150 {gMessMgr->Error("Init() : wrong RotatingUsage parameter set.");
00151 return kStErr;
00152 }
00153
00154 if (useTracker == kTrackerUseTPT) gMessMgr->Info()<<"use TPT tracks."<<endm;
00155 if (useTracker == kTrackerUseITTF) gMessMgr->Info()<<"use ITTF tracks."<<endm;
00156 if (useTracker == kTrackerUseBOTH) gMessMgr->Info()<<"use TPT *and* ITTF tracks."<<endm;
00157 if (useSVT == kUseSVT) gMessMgr->Info()<<"use SVT points if possible."<<endm;
00158 if (useSVT == kNoSVT) gMessMgr->Info()<<"do not use SVT points."<<endm;
00159 if (useEventModel == kUseStEvent) gMessMgr->Info()<<"expect StEvent files in input."<<endm;
00160 if (useEventModel == kUseMuDst) gMessMgr->Info()<<"expect MuDst files in input."<<endm;
00161 if (useLikesign == kLikesignUseLikesign) gMessMgr->Info()<<"does like-sign finding."<<endm;
00162 if (useRotating == kRotatingUseRotating) gMessMgr->Info()<<"does rotating finding."<<endm;
00163 if (useRotating == kRotatingUseSymmetry) gMessMgr->Info()<<"does symmetry finding."<<endm;
00164 if (useRotating == kRotatingUseRotatingAndSymmetry) gMessMgr->Info()<<"does rotating + symmetry finding."<<endm;
00165
00166 if (useLanguage != kLanguageUseSpecial)
00167 {a=(bool)(1&(useLanguage>>2));
00168 b=(bool)(1&(useLanguage>>1));
00169 c=(bool)(1&useLanguage);
00170 useV0Language=2*(!(a^c))+(a|c);
00171 useXiLanguage=4*(b&(!(a^c)))+2*(a&b&(!c))+(a|c);
00172 }
00173 switch (useLanguage)
00174 {case kLanguageUseOldRun : gMessMgr->Info()<<"Fortran run."<<endm;
00175 break;
00176 case kLanguageUseRun : gMessMgr->Info()<<"C++ run."<<endm;
00177 gMessMgr->Info()<<"BE CAREFUL : you are NOT running the XiFinder !"<<endm;
00178 break;
00179 case kLanguageUseTestV0Finder : gMessMgr->Info()<<"Test V0Finder."<<endm;
00180 break;
00181 case kLanguageUseTestXiFinder : gMessMgr->Info()<<"Test XiFinder."<<endm;
00182 gMessMgr->Info()<<"BE CAREFUL : you are NOT running the XiFinder !"<<endm;
00183 break;
00184 case kLanguageUseTestBothFinders : gMessMgr->Info()<<"Test V0Finder and XiFinder."<<endm;
00185 gMessMgr->Info()<<"BE CAREFUL : you are NOT running the XiFinder !"<<endm;
00186 break;
00187 case kLanguageUseSpecial : break;
00188 default : gMessMgr->Error("Init() : wrong LanguageUsage parameter set.");
00189 return kStErr;
00190 }
00191 if ((useV0Language!=kV0LanguageUseFortran) && (useV0Language!=kV0LanguageUseCpp) && (useV0Language!=kV0LanguageUseBoth))
00192 {gMessMgr->Error("Init() : wrong V0LanguageUsage parameter set.");
00193 return kStErr;
00194 }
00195 if (1&useV0Language) gMessMgr->Info()<<" Will store Fortran V0."<<endm;
00196 if (2&useV0Language) gMessMgr->Info()<<" Will store C++ V0."<<endm;
00197 if (1&useXiLanguage) gMessMgr->Info()<<" Will store Fortran Xi."<<endm;
00198 if (2&useXiLanguage) gMessMgr->Info()<<" BE CAREFUL : will NOT store C++ Xi, although asked."<<endm;
00199 if (4&useXiLanguage) gMessMgr->Info()<<" BE CAREFUL : will NOT store C++ Xi, although asked."<<endm;
00200
00201 if (useEventModel) {
00202 mMuDstMaker = (StMuDstMaker*)GetMaker("myMuDstMaker");
00203 if(!mMuDstMaker) gMessMgr->Warning("Init can't find a valid MuDst");
00204 }
00205
00206 return StMaker::Init();
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 Int_t StV0FinderMaker::Prepare() {
00232
00233 if (prepared) return kStOk;
00234
00235 unsigned short i,j,nNodes,nTrksEstimate;
00236 StThreeVectorD p;
00237
00238
00239 GetPars();
00240 ITTFflag=kITKalmanFitId;
00241 TPTflag=kHelix3DIdentifier;
00242
00243
00244
00246 if(GetEventUsage()==kUseStEvent){
00247 event = (StEvent*) GetInputDS("StEvent");
00248 }
00250 else if(GetEventUsage()==kUseMuDst){
00251 StMuDst* mu = mMuDstMaker->muDst();
00252 if(mu) cout<<"V0Finder :: found a MuDst"<<endl;
00253 if(mu->event())cout<<"see a muEvent ... "<<endl;
00254 Int_t nV0s = mu->GetNV0();
00255 if(mu) event = mu->createStEvent();
00256
00257 if(event) AddData(event);
00258 if(event)cout<<"see a recreated StEvent!"<<endl;
00259 cout<<"the number of existing v0's is "<<nV0s<<endl;
00260 }
00262
00263 if (!event)
00264 {gMessMgr->Warning("no StEvent ; skipping event.");
00265 return kStWarn;
00266 }
00267
00268
00269 StPrimaryVertex* pvert = event->primaryVertex();
00270 if (!pvert)
00271 {gMessMgr->Warning("no primary vertex ; skipping event.");
00272 return kStWarn;
00273 }
00274 mainv = pvert->position();
00275
00276 StSPtrVecTrackNode& theNodes = event->trackNodes();
00277 nNodes = theNodes.size();
00278
00279
00280 nTrksEstimate = (unsigned int) (nNodes*trkNodeRatio);
00281 if (nTrksEstimate > trk.size()) ExpandVectors(nTrksEstimate);
00282
00283
00284 trks=0;
00285 double BfieldRunning = 0;
00286 double nBfieldRunning = 0;
00287 for (i=0; i<nNodes; i++) {
00288 int nj = theNodes[i]->entries(global);
00289 for (j=0; j<nj; j++) {
00290
00291 StTrack* tri = theNodes[i]->track(global,j);
00293
00294 if(useSVT){
00295 StTrack* svtTrack = theNodes[i]->track(estGlobal,j);
00296 if (svtTrack){
00297 tri=svtTrack;
00298 }
00299 }
00301
00302
00303
00304
00305 if (
00306
00307 (tri->fittingMethod() != ITTFflag && (GetTrackerUsage() == kTrackerUseITTF)) ||
00308 (tri->fittingMethod() == ITTFflag && (GetTrackerUsage() == kTrackerUseTPT))) continue;
00309
00310
00311 if (tri->flag() <= 0) continue;
00312
00313
00314 if (trks >= trk.size()) ExpandVectors(trks+1);
00315
00316
00317 const StTrackTopologyMap& map = tri->topologyMap();
00318 Bool_t tpcHit = map.hasHitInDetector(kTpcId);
00319 Bool_t silHit = map.hasHitInDetector(kSvtId) ||
00320 map.hasHitInDetector(kSsdId);
00321 if (tpcHit) {
00322 if (silHit)
00323 detId[trks] = 3;
00324 else
00325 detId[trks] = 1;
00326 } else if (silHit)
00327 detId[trks] = 2;
00328 else
00329
00330 continue;
00331
00332 trk[trks] = tri;
00333
00334 StTrackGeometry* triGeom = tri->geometry();
00335 if(!triGeom) continue;
00336
00337 heli[trks] = triGeom->helix();
00338
00339 p = triGeom->momentum();
00340
00341 pt[trks] = p.perp();
00342 ptot[trks] = p.mag();
00343 trkID[trks]=tri->key();
00344
00345
00346 hits[trks] = map.numberOfHits(kTpcId) +
00347 map.numberOfHits(kSvtId) +
00348 map.numberOfHits(kSsdId);
00349
00350 pars2 = v0pars->GetTable(detId[trks]-1);
00351 if (hits[trks] < pars2->n_point) continue;
00352
00353
00354 if (nBfieldRunning<1e2) {
00355 StThreeVectorD p1 = triGeom->momentum();
00356 StThreeVectorD p2 = heli[trks].momentum(Bfield);
00357
00358 if (fabs(p2.x()) > fabs(p2.y())) Bfield *= p1.x()/p2.x();
00359 else if (fabs(p2.y()) > 1.e-20) Bfield *= p1.y()/p2.y();
00360 else continue;
00361
00362
00363
00364
00365
00366
00367 if (fabs(Bfield)<1.e-20) return kStWarn;
00368 Bfield = TMath::Sign(Bfield,-1.0*triGeom->charge()*triGeom->helicity());
00369
00370 BfieldRunning += Bfield;
00371 nBfieldRunning += 1;
00372 }
00373 if (triGeom->charge() > 0) ptrk.push_back(trks);
00374 else if (triGeom->charge() < 0) ntrk.push_back(trks);
00375 trks++;
00376 }
00377 }
00378 if (nBfieldRunning>0) Bfield = BfieldRunning/nBfieldRunning;
00379
00380
00381
00382 if (trks < (trk.size()/2)) {
00383 if (trks > trkmax) trkmax = trks;
00384 trkcnt++;
00385 if (trkcnt >= 10) {
00386 ExpandVectors(trkmax);
00387 trkcnt = 0;
00388 trkmax = BLOCK;
00389 }
00390 } else {
00391 trkcnt = 0;
00392 trkmax = BLOCK;
00393 }
00394 trkNodeRatio = ((trkNodeRatio*trkNodeRatioCnt)+((float) trks)/((float) nNodes)) /
00395 (trkNodeRatioCnt + 1.);
00396 trkNodeRatioCnt++;
00397
00398 gMessMgr->Info() << "No. of nodes is : "
00399 << nNodes << endm;
00400 gMessMgr->Info() << "No. of tracks is : "
00401 << trks << endm;
00402
00403 prepared = kTRUE;
00404 return kStOk;
00405 }
00406
00407
00408
00409
00410
00411
00412
00413
00414
00415
00416
00417
00418 Int_t StV0FinderMaker::Make() {
00419
00420
00421 StThreeVectorD xi,xj,pi,pj,xpp,pp,impact;
00422 StThreeVectorD xi1,xj1,pi1,pj1,tmp3V;
00423 TVector2 ri,rj,xci,xcj,tmp2V;
00424 double rad_i,rad_j,separation,solution,dxc;
00425 double dca_ij,dca_ij1,rmin,dlen,pperpsq,ppmag;
00426 double alpha,ptArm_sq,pPosAlongV0,pNegAlongV0;
00427 double cosij,sin2ij,t1,t2;
00428 unsigned short i,j,ii,jj;
00429 pairD paths,paths1,path2;
00430 double temp;
00431 Bool_t doSecond, isPrimaryV0, usedV0;
00432 Int_t iRes;
00433 long keepTrack;
00434
00435 if (! (2&useV0Language)) return kStOk;
00436
00437 gMessMgr->Info("Make() : Starting...");
00438
00439
00440 iRes = Prepare();
00441 if (iRes != kStOk) return iRes;
00442
00443
00444 StSPtrVecV0Vertex& v0Vertices = event->v0Vertices();
00445 gMessMgr->Info()<<"coming in I have "<<v0Vertices.size()<<" V0s."<<endm;
00446
00447 if (!(1&useV0Language)) {
00448
00449
00450 gMessMgr->Info()<<"pre-existing V0s and Xis deleted."<<endm;
00451 StSPtrVecV0Vertex v0Vertices2;
00452 v0Vertices = v0Vertices2;
00453 StSPtrVecXiVertex& xiVertices = event->xiVertices();
00454 StSPtrVecXiVertex xiVertices2;
00455 xiVertices = xiVertices2;
00456 }
00457
00458
00459
00460
00461
00462 int nii = ptrk.size();
00463 for (ii=0; ii<nii; ii++) {
00464 i = ptrk[ii];
00465
00466 xci.Set(heli[i].xcenter(),heli[i].ycenter());
00467 ri.Set(heli[i].origin().x(),heli[i].origin().y());
00468 ri -= xci;
00469 rad_i = ri.Mod();
00470
00471
00472 int njj = ntrk.size();
00473 for (jj=0; jj<njj; jj++) {
00474 j = ntrk[jj];
00475
00476 if (GetTrackerUsage() == kTrackerUseBOTH)
00477 {
00478
00479 if ((trk[i]->fittingMethod() == ITTFflag) && (trk[j]->fittingMethod() != ITTFflag)) continue;
00480
00481 if ((trk[i]->fittingMethod() != ITTFflag) && (trk[j]->fittingMethod() == ITTFflag)) continue;
00482 }
00483
00484
00485 det_id_v0 = TMath::Max(detId[i],detId[j]);
00486
00487
00488 pars = v0pars->GetTable(det_id_v0+2);
00489
00490 pars2 = v0pars->GetTable(det_id_v0-1);
00491
00492
00493
00494
00495
00496
00497
00498
00499
00500 temp = pt[i] + pt[j];
00501 if ((temp < 0.99 * pars2->dcapn_pt) &&
00502 ((trk[i]->impactParameter() <= pars2->dcapnmin) ||
00503 (trk[j]->impactParameter() <= pars2->dcapnmin))) continue;
00504
00505 xcj.Set(heli[j].xcenter(),heli[j].ycenter());
00506 rj.Set(heli[j].origin().x(),heli[j].origin().y());
00507 rj -= xcj;
00508 rad_j = rj.Mod();
00509
00510 tmp2V = xci - xcj;
00511 dxc = tmp2V.Mod();
00512 separation = dxc - (rad_i + rad_j);
00513 doSecond = kFALSE;
00514 dca_ij1 = -9999;
00515
00516
00517
00518 if (separation < 0)
00519 {
00520 if (dxc < TMath::Abs(rad_i - rad_j)) continue;
00521
00522
00523
00524
00525 path2 = heli[i].pathLength(rad_j,xcj.X(),xcj.Y());
00526
00527 if (!isnan(path2.first))
00528 {solution = path2.first;
00529 if ((!isnan(path2.second)) && (path2.second != path2.first))
00530 {doSecond = kTRUE;
00531 }
00532 goto ProcessSolution;
00533 }
00534 else if (isnan(path2.second))
00535 {
00536 continue;
00537 }
00538
00539 SecondSolution:
00540 solution = path2.second;
00541 doSecond = kFALSE;
00542 ProcessSolution:
00543
00544
00545
00546 paths.first = solution;
00547 xi = heli[i].at(paths.first );
00548 paths.second = heli[j].pathLength(xi.x(),xi.y());
00549 xj = heli[j].at(paths.second);
00550 }
00551 else if (separation < pars2->dca)
00552 {
00553
00554 tmp2V = (xci + xcj) * 0.5;
00555 paths.first = heli[i].pathLength(tmp2V.X(),tmp2V.Y());
00556 paths.second = heli[j].pathLength(tmp2V.X(),tmp2V.Y());
00557 xi = heli[i].at(paths.first );
00558 xj = heli[j].at(paths.second);
00559 }
00560 else
00561 {
00562 continue;
00563 }
00564
00565 dca_ij = xi.z() - xj.z();
00566 if (doSecond) {
00567
00568 dca_ij1 = dca_ij;
00569 xi1=xi;
00570 xj1=xj;
00571 paths1=paths;
00572 goto SecondSolution;
00573 }
00574 if ((dca_ij1 != -9999) &&
00575 (TMath::Abs(dca_ij1) < TMath::Abs(dca_ij))) {
00576
00577 dca_ij = dca_ij1;
00578 xi=xi1;
00579 xj=xj1;
00580 paths=paths1;
00581 }
00582
00583
00584
00585 pi = heli[i].momentumAt(paths.first ,Bfield);
00586 pj = heli[j].momentumAt(paths.second,Bfield);
00587
00588
00589 if ((pi.dot(xi-mainv) < 0.0) ||
00590 (pj.dot(xj-mainv) < 0.0)) continue;
00591
00593
00594
00595
00596
00597
00599
00600
00601
00602
00603 pi1 = pi/ptot[i];
00604 pj1 = pj/ptot[j];
00605
00606 cosij = pi1.dot(pj1);
00607 sin2ij = 1.0 - cosij*cosij;
00608 if (sin2ij) {
00609 temp = dca_ij/sin2ij;
00610 t1 = (-pi1.z()+pj1.z()*cosij)*temp;
00611 t2 = ( pj1.z()-pi1.z()*cosij)*temp;
00612
00613 temp = rad_i*(ptot[i]/pt[i]);
00614 temp *= sin(t1/temp);
00615 xi1 = xi + pi1.pseudoProduct(temp,temp,t1);
00616
00617 temp = rad_j*(ptot[j]/pt[j]);
00618 temp *= sin(t2/temp);
00619 xj1 = xj + pj1.pseudoProduct(temp,temp,t2);
00620
00621 dca_ij1 = (xi1 - xj1).mag2();
00622 dca_ij *= dca_ij;
00623
00624 if (dca_ij1 < dca_ij) {
00625 paths.first = heli[i].pathLength(xi1.x(),xi1.y());
00626 paths.second = heli[j].pathLength(xj1.x(),xj1.y());
00629 xi1 = heli[i].at(paths.first);
00630 xj1 = heli[j].at(paths.second);
00631 dca_ij1 = (xi1 - xj1).mag2();
00632 if (dca_ij1 < dca_ij) {
00633 xi = xi1;
00634 xj = xj1;
00635 pi = heli[i].momentumAt(paths.first ,Bfield);
00636 pj = heli[j].momentumAt(paths.second,Bfield);
00637 dca_ij = dca_ij1;
00638 }
00639 }
00640
00641
00642
00643
00644
00645
00646
00647
00648
00649
00650
00651
00652
00653
00654 }
00655
00656
00657
00658
00659
00660
00661
00662
00663 if (dca_ij >= (pars2->dca*pars2->dca)) continue;
00664
00665 pp = pi + pj;
00666 pperpsq = pp.perp2();
00667
00668
00669 if ((pperpsq < pars2->dcapn_pt * pars2->dcapn_pt) &&
00670 ((trk[i]->impactParameter() <= pars2->dcapnmin) ||
00671 (trk[j]->impactParameter() <= pars2->dcapnmin))) continue;
00672
00673 xpp = (xi + xj) * 0.5;
00674 impact = xpp - mainv;
00675 dlen = impact.mag2();
00676
00677
00678 if (dlen <= (pars2->dlen*pars2->dlen)) continue;
00679
00680
00681 if (pp.dot(impact) < 0.0) continue;
00682
00683 ppmag = pperpsq + pp.z()*pp.z();
00684 rmin = impact.cross(pp).mag2()/ppmag;
00685
00686
00687 if (rmin >= (pars2->dcav0*pars2->dcav0)) continue;
00688
00689 tmp3V = pp/::sqrt(ppmag);
00690 pPosAlongV0 = pi.dot(tmp3V);
00691 pNegAlongV0 = pj.dot(tmp3V);
00692 alpha = (pPosAlongV0-pNegAlongV0) /
00693 (pPosAlongV0+pNegAlongV0);
00694
00695
00696 if (TMath::Abs(alpha) > pars2->alpha_max) continue;
00697
00698 ptArm_sq = ptot[i]*ptot[i] - pPosAlongV0*pPosAlongV0;
00699
00700
00701 if (ptArm_sq > (pars2->ptarm_max*pars2->ptarm_max)) continue;
00702
00703 rmin = ::sqrt(rmin);
00704 dca_ij = ::sqrt(dca_ij);
00705 if (trk[i]->fittingMethod() == ITTFflag) dca_ij=-dca_ij;
00706
00707
00708 v0Vertex = new StV0Vertex();
00709 v0Vertex->setPosition(xpp);
00710 v0Vertex->addDaughter(trk[i]);
00711 v0Vertex->addDaughter(trk[j]);
00712 v0Vertex->setDcaDaughterToPrimaryVertex(positive,trk[i]->impactParameter());
00713 v0Vertex->setDcaDaughterToPrimaryVertex(negative,trk[j]->impactParameter());
00715 v0Vertex->setMomentumOfDaughter(positive,pi);
00716 v0Vertex->setMomentumOfDaughter(negative,pj);
00717 v0Vertex->setDcaDaughters(dca_ij);
00718 v0Vertex->setDcaParentToPrimaryVertex(rmin);
00719
00721
00722
00723
00724
00725
00726
00728 keepTrack=0;
00729 keepTrack |=((long)1 << 4);
00730 if (useSVT)
00731 {keepTrack |=((long)1 << 3);
00732 }
00733 if (detId[i]==2 || detId[i]==3){
00734 keepTrack |=((long)1 << 2);
00735 }
00736 if(detId[j]==2 || detId[j]==3){
00737 keepTrack |=((long)1 << 1);
00738 }
00739
00740 keepTrack *= -1;
00741 v0Vertex->setChiSquared((float)keepTrack);
00743
00744
00745 isPrimaryV0 =
00746 (rmin < pars->dcav0) &&
00747 ((pperpsq >= pars->dcapn_pt * pars->dcapn_pt) ||
00748 ((trk[i]->impactParameter() > pars->dcapnmin) &&
00749 (trk[j]->impactParameter() > pars->dcapnmin)));
00750
00751
00752 usedV0 = UseV0();
00753
00754
00755 if (usedV0 && !isPrimaryV0)
00756 v0Vertex->setDcaParentToPrimaryVertex(-rmin);
00757
00758
00759 if (usedV0 || isPrimaryV0) {
00760 v0Vertices.push_back(v0Vertex);
00761 } else {
00762 delete v0Vertex;
00763 v0Vertex = 0;
00764 }
00765
00766 }
00767 }
00768
00769 gMessMgr->Info()<<"now I have "<<v0Vertices.size()<<" V0s."<<endm;
00770 gMessMgr->Info()<<"using magnetic field : "<<Bfield/tesla<<" T."<<endm;
00771
00772
00773
00774 return kStOk;
00775 }
00776
00777
00778
00779
00780
00781
00782
00783
00784
00785
00786
00787
00788 void StV0FinderMaker::Clear(Option_t *option){
00789 prepared = kFALSE;
00790 ptrk.clear();
00791 ntrk.clear();
00792 if(useEventModel){
00793 if(mMuDstMaker){
00794 if(event){
00795 delete event;
00796 event=0;
00797 }
00798 }
00799 }
00800 }
00801
00802
00803
00804
00805
00806
00807
00808
00809
00810
00811
00812
00813
00814 void StV0FinderMaker::Trim() {
00815
00816
00817 gMessMgr->Info() << "Trim() : Starting..." << endm;
00818
00819 event = (StEvent*) GetInputDS("StEvent");
00820 pars = v0pars->GetTable(3);
00821 StSPtrVecV0Vertex& v0Vertices = event->v0Vertices();
00822 int iV0s = v0Vertices.size();
00823 for (int i=iV0s-1; i>=0; i--) {
00824
00825 v0Vertex = v0Vertices[i];
00826 if ((v0Vertex) &&
00827
00828 (v0Vertex->dcaParentToPrimaryVertex() >= 0) &&
00829
00830 ! ((v0Vertex->dcaParentToPrimaryVertex() < pars->dcav0) &&
00831 ((v0Vertex->momentum().perp2() >= pars->dcapn_pt * pars->dcapn_pt) ||
00832 ((v0Vertex->dcaDaughterToPrimaryVertex(positive) > pars->dcapnmin) &&
00833 (v0Vertex->dcaDaughterToPrimaryVertex(negative) > pars->dcapnmin))))) {
00834 v0Vertex->makeZombie();
00835 iV0s--;
00836 }
00837 }
00838
00839 gMessMgr->Info() << "Trim() : saving " << iV0s <<
00840 " V0 candidates" << endm;
00841 }
00842
00843 void StV0FinderMaker::ExpandVectors(unsigned short size) {
00844 unsigned int oldsize = trk.size();
00845 unsigned int newsize = oldsize;
00846 if (newsize > size) newsize = BLOCK;
00847 while (newsize <= size) newsize += BLOCK;
00848 if (newsize == oldsize) return;
00849 gMessMgr->Info() << IsA()->GetName() << "::ExpandVectors(" << newsize
00850 << ") for " << GetName() << endm;
00851 trk.resize(newsize);
00852 hits.resize(newsize);
00853 detId.resize(newsize);
00854 pt.resize(newsize);
00855 ptot.resize(newsize);
00856 heli.resize(newsize);
00857 trkID.resize(newsize);
00858 }
00859
00860
00861
00862
00863
00864
00865
00866
00867
00868
00869
00870
00871
00872
00873
00874
00875
00876
00877
00878
00879
00880
00881
00882
00883
00884
00885
00886
00887
00888
00889
00890
00891
00892
00893
00894
00895
00896
00897
00898
00899
00900
00901
00902
00903
00904
00905
00906
00907
00908
00909
00910
00911
00912
00913
00914
00915
00916
00917
00918
00919
00920
00921
00922
00923
00924
00925
00926
00927
00928
00929
00930
00931
00932
00933
00934
00935
00936
00937
00938
00939
00940
00941
00942
00943
00944
00945
00946
00947
00948
00949
00950
00951
00952
00953
00954
00955
00956
00957
00958