StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StMinuitVertexFinder.cxx
1 /***************************************************************************
2  *
3  * $Id: StMinuitVertexFinder.cxx,v 1.59 2017/05/09 12:29:40 smirnovd Exp $
4  *
5  * Author: Thomas Ullrich, Feb 2002
6  ***************************************************************************
7  *
8  * Description:
9  *
10  **************************************************************************/
11 #include <vector>
12 
13 #include "StEmcUtil/geometry/StEmcGeom.h"
14 #include "StEvent/StCtbTriggerDetector.h"
15 #include "StEvent/StDcaGeometry.h"
16 #include "StEvent/StEmcCollection.h"
17 #include "StEvent/StEmcDetector.h"
18 #include "StEvent/StEmcModule.h"
19 #include "StEvent/StEmcRawHit.h"
20 #include "StEvent/StEvent.h"
21 #include "StEvent/StGlobalTrack.h"
22 #include "StEvent/StTrackDetectorInfo.h"
23 #include "StEvent/StTrackNode.h"
24 #include "StEvent/StTriggerDetectorCollection.h"
25 #include "StGenericVertexMaker/Minuit/StMinuitVertexFinder.h"
26 #include "StGenericVertexMaker/Minuit/St_VertexCutsC.h"
27 #include "StGenericVertexMaker/StCtbMatcher.h"
28 #include "St_base/StMessMgr.h"
29 
30 
31 std::vector<StPhysicalHelixD> StMinuitVertexFinder::mHelices;
32 std::vector<UShort_t> StMinuitVertexFinder::mHelixFlags;
33 std::vector<Double_t > StMinuitVertexFinder::mZImpact;
34 Bool_t StMinuitVertexFinder::requireCTB;
35 Int_t StMinuitVertexFinder::nCTBHits;
36 //==========================================================
37 //==========================================================
38 void StMinuitVertexFinder::Clear(){
39  StGenericVertexFinder::Clear();
40  mStatusMin = 0;
41  mNSeed = 0;
42 }
43 
44 void
45 StMinuitVertexFinder::setExternalSeed(const StThreeVectorD& s)
46 {
47  mExternalSeedPresent = kTRUE;
48  mExternalSeed = s;
49 }
50 
51 
52 StMinuitVertexFinder::StMinuitVertexFinder(VertexFit_t fitMode) :
53  StGenericVertexFinder(SeedFinder_t::MinuitVF, fitMode)
54 {
55  mExternalSeedPresent = kFALSE;
56  mRequireCTB = kFALSE;
57  requireCTB = kFALSE;
58  mUseITTF = kFALSE;
59  mUseOldBEMCRank = kFALSE;
60  mLowerSplitVtxRank = kFALSE;
61  mVertexOrderMethod = orderByRanking; // change ordering by ranking
62  mMinTrack = -1;
63  mCTBSum = 0;
64 }
65 
66 
67 StMinuitVertexFinder::~StMinuitVertexFinder()
68 {
69  LOG_WARN << "Skipping delete Minuit in StMinuitVertexFinder::~StMinuitVertexFinder()" << endm;
70  mHelices.clear();
71  mHelixFlags.clear();
72  mZImpact.clear();
73 }
74 //________________________________________________________________________________
75 void StMinuitVertexFinder::InitRun(int run_number, const St_db_Maker* db_maker)
76 {
77  StGenericVertexFinder::InitRun(run_number, db_maker);
78 
79  St_VertexCutsC *cuts = St_VertexCutsC::instance();
80  mMinNumberOfFitPointsOnTrack = cuts->MinNumberOfFitPointsOnTrack();
81  mDcaZMax = cuts->DcaZMax(); // Note: best to use integer numbers
82  mMinTrack = (mMinTrack<0 ? cuts->MinTrack() : mMinTrack);
83  mRImpactMax = cuts->RImpactMax();
84  mZMin = cuts->ZMin(); // note: best to use integer numbers
85  mZMax = cuts->ZMax(); // note: best to use integer numbers
86  if (mZMin == mZMax) {
87  // historical defaults
88  mZMin = -200.0;
89  mZMax = 200.0;
90  }
91  LOG_INFO << "Set cuts: MinNumberOfFitPointsOnTrack = " << mMinNumberOfFitPointsOnTrack
92  << " DcaZMax = " << mDcaZMax
93  << " MinTrack = " << mMinTrack
94  << " RImpactMax = " << mRImpactMax
95  << " ZMin = " << mZMin
96  << " ZMax = " << mZMax
97  << endm;
98 }
99 //________________________________________________________________________________
100 
101 
102 void
103 StMinuitVertexFinder::setFlagBase(){
104  if(mUseITTF){
105  mFlagBase = 8000;
106  } else {
107  mFlagBase = 1000;
108  }
109 }
110 
111 Int_t StMinuitVertexFinder::findSeeds() {
112  mNSeed = 0;
113 
114  int zIdxMax = (int) (mZMax - mZMin);
115  Int_t zImpactArr[zIdxMax]; // simple array to 'histogram' zImpacts
116  for (Int_t i=0; i < zIdxMax; i++)
117  zImpactArr[i]=0;
118 
119  Int_t nTrk = mZImpact.size();
120  for (Int_t iTrk=0; iTrk < nTrk; iTrk++) {
121  if ((mZImpact[iTrk] > mZMin) &&
122  (mZImpact[iTrk] < mZMax))
123  zImpactArr[int(mZImpact[iTrk]-mZMin)]++;
124  }
125 
126  // Search for maxima using sliding 3-bin window
127  Int_t nOldBin = 0;
128  Int_t slope = 0;
129  Int_t nBinZ = 3;
130  for (Int_t iBin=0; iBin < zIdxMax - nBinZ; iBin++) {
131  Int_t nTrkBin = 0;
132  for (Int_t iBin2=0; iBin2 < nBinZ; iBin2++) {
133  nTrkBin += zImpactArr[iBin + iBin2];
134  }
135  if (nTrkBin > nOldBin)
136  slope = 1;
137  else if (nTrkBin < nOldBin) {
138  if (slope == 1) {
139  if (mNSeed < maxSeed) {
140  Float_t seed_z = mZMin + iBin + (Float_t)nBinZ / 2 - 1;
141  Double_t meanZ = 0;
142  Int_t nTrkZ = 0;
143  for (Int_t iTrk = 0; iTrk < nTrk; iTrk ++ ) {
144  if ( fabs(mZImpact[iTrk] - seed_z ) < mDcaZMax ) {
145  meanZ += mZImpact[iTrk];
146  nTrkZ++;
147  }
148  }
149  if (nTrkZ >= mMinTrack) {
150  if (mDebugLevel) {
151  LOG_INFO << "Seed " << mNSeed << ", z " << seed_z << " nTrk " << nTrkZ << " meanZ/nTrkZ " << meanZ/nTrkZ << endm;
152  }
153  seed_z = meanZ/nTrkZ;
154  mSeedZ[mNSeed] = seed_z;
155  mNSeed++;
156  }
157  }
158  else {
159  LOG_WARN << "Too many vertex seeds, limit=" << maxSeed << endm;
160  }
161  }
162  slope = -1;
163  }
164  if (mDebugLevel > 1) {
165  LOG_INFO << "iBin " << iBin << " nTrkBin " << nTrkBin << " nOldBin " << nOldBin << ", slope " << slope << " mNSeed " << mNSeed << endm;
166  }
167  nOldBin = nTrkBin;
168  }
169 
170  LOG_INFO << "Found " << mNSeed << " seeds" << endm;
171  return mNSeed;
172 }
173 
174 void StMinuitVertexFinder::fillBemcHits(StEvent *event){
175  static Int_t nMod = 120;
176  static Int_t nEta = 20;
177  static Int_t nSub = 2;
178  static Float_t mEmcThres = 0.15;
179  for (Int_t m=0; m < nMod; m++)
180  for (Int_t e=0; e < nEta; e++)
181  for (Int_t s=0; s < nSub; s++)
182  mBemcHit[m][e][s]=0;
183 
184  Int_t n_emc_hit=0;
185  if (event->emcCollection() && event->emcCollection()->detector(kBarrelEmcTowerId)) {
186  StEmcDetector* bemcDet = event->emcCollection()->detector(kBarrelEmcTowerId);
187  for (Int_t iMod=0; iMod < nMod; iMod++) {
188  if (!bemcDet->module(iMod))
189  continue;
190  StEmcModule *mod = bemcDet->module(iMod);
191  StSPtrVecEmcRawHit& hits = mod->hits();
192  for (StEmcRawHitIterator hitIter = hits.begin(); hitIter != hits.end(); hitIter++) {
193  StEmcRawHit *hit = *hitIter;
194  if (hit->energy() > mEmcThres) {
195  mBemcHit[hit->module()-1][hit->eta()-1][hit->sub()-1]=1;
196  n_emc_hit++;
197  }
198  }
199  }
200  }
201  if (mDebugLevel) {
202  LOG_INFO << "Found " << n_emc_hit << " emc hits" << endm;
203  }
204 }
205 
206 Int_t
207 StMinuitVertexFinder::matchTrack2BEMC(const StTrack *track){
208  static const Double_t rBemc = 242; // middle of tower
209  static StEmcGeom *bemcGeom = StEmcGeom::getEmcGeom("bemc");
210  //static Double_t rBemc = bemcGeom->Radius(); // front face??
211  //LOG_INFO << "rBemc: " << rBemc << endm;
212 
213  if (track->outerGeometry()==0) {
214  if (mDebugLevel) { // Happens only rarely
215  LOG_INFO << "No outer track geom" << endm;
216  }
217  return 0;
218  }
219 
220  StPhysicalHelixD helix = track->outerGeometry()->helix();
221 
222  if (!helix.valid()) {
223  if (mDebugLevel) { // Happens only rarely
224  LOG_INFO << "Invalid helix" << endm;
225  }
226  return 0;
227  }
228 
229  pairD d2;
230  d2 = helix.pathLength(rBemc);
231  if (d2.first > 99999999 && d2.second > 99999999)
232  return 0;
233  Double_t path=d2.second;
234  if (d2.first >= 0 && d2.first < d2.second)
235  path = d2.first;
236  else if(d2.first>=0 || d2.second<=0) {
237  LOG_WARN << "Unexpected solution for track crossing Cyl:"
238  << " track mom = "
239  << track->geometry()->momentum().mag()
240  << ", d2.first=" << d2.first
241  << ", second=" << d2.second <<", trying first" << endm;
242  path=d2.first;
243  }
244  StThreeVectorD posCyl = helix.at(path);
245 
246  Float_t phi=posCyl.phi();
247  Float_t eta=posCyl.pseudoRapidity();
248 
249  if (fabs(eta) < 1) {
250  Int_t mod, ieta, sub;
251  if (bemcGeom->getBin(phi, eta, mod, ieta, sub) == 0) {
252  // There is some edge effect leading to sub=-1.
253  // Strange, but leave in for now// HOW CAN IT be: sub == -1????
254  if (sub > 0 && mBemcHit[mod-1][ieta-1][sub-1]) {
255  return 1;
256  }
257  }
258  }
259  return 0;
260 }
261 
262 Int_t StMinuitVertexFinder::checkCrossMembrane(const StTrack *track) {
263  Int_t n_pnt_neg=0, n_pnt_pos=0;
264 
265  StPtrVecHit hits = track->detectorInfo()->hits(kTpcId);
266  for (StHitIterator hitIter = hits.begin(); hitIter != hits.end(); hitIter++) {
267  if ((*hitIter)->position().z() < 0)
268  n_pnt_neg++;
269  if ((*hitIter)->position().z() > 0)
270  n_pnt_pos++;
271  }
272  return (n_pnt_pos > 5 && n_pnt_neg > 5);
273 }
274 
275 void StMinuitVertexFinder::calculateRanks() {
276 
277  // Calculation of veretx ranks to select 'best' (i.e. triggered)
278  // vertex
279  // Three ranks are used, based on average dip, number of BEMC matches
280  // and tarcks crossing central membrane
281  // Each rank is normalised to be 0 for 'average' valid vertices and
282  // has a sigma of 1.
283  //
284  // Values are limited to [-5,1] for each rank
285  //
286  // Note that the average dip angle ranking is naturally peaked to 1,
287  // while the others are peaked at 1 (intentionally) due to rounding.
288 
289  // A fancier way would be to calculate something more like
290  // a likelihood based on the expected distributions.
291  // That's left as an excercise to the reader.
292 
293  // Added by Ant: split vertices have 3 rank units deducted near the
294  // end of function (if mLowerSplitVtxRank). See hypernews post:
295  // http://www.star.bnl.gov/HyperNews-star/get/vertex/127.html
296 
297  Int_t nBemcMatchTot = 0;
298  Int_t nVtxTrackTot = 0;
299  for (Int_t iVertex=0; iVertex < size(); iVertex++) {
300  StPrimaryVertex *primV = getVertex(iVertex);
301  nVtxTrackTot += primV->numTracksUsedInFinder();
302  nBemcMatchTot += primV->numMatchesWithBEMC();
303  }
304 
305  mBestRank = -999;
306  mBestVtx = 0;
307  for (Int_t iVertex=0; iVertex < size(); iVertex++) {
308  StPrimaryVertex *primV = getVertex(iVertex);
309  // expected values based on Cu+Cu data
310  Float_t avg_dip_expected = -0.0033*primV->position().z();
311  Float_t n_bemc_expected = 0;
312  if (nVtxTrackTot) {
313  if (mUseOldBEMCRank)
314  n_bemc_expected = (1-0.25*(1-(float)primV->numTracksUsedInFinder()/nVtxTrackTot))*nBemcMatchTot;
315  else
316  n_bemc_expected = (float)primV->numTracksUsedInFinder()/nVtxTrackTot*nBemcMatchTot;
317  }
318 
319  Float_t n_cross_expected = fabs(primV->position().z())*0.0020*primV->numTracksUsedInFinder(); // old coeff 0.0016 with dca 3 and 10 points on track
320 
321  if (mDebugLevel) {
322  LOG_INFO << "vertex z " << primV->position().z() << " dip expected " << avg_dip_expected << " bemc " << n_bemc_expected << " cross " << n_cross_expected << endm;
323  }
324 
325  Float_t rank_avg_dip = 1 - fabs(primV->meanDip() - avg_dip_expected)*sqrt((float)primV->numTracksUsedInFinder())/0.67; // Sigma was 0.8 for old cuts
326 
327  Float_t rank_bemc = 0;
328  if (n_bemc_expected >= 1) {
329  //Float_t sigma = 0.12*n_bemc_match_tot;
330  Float_t sigma = 0.5*sqrt(n_bemc_expected);
331 
332  // limit sigma to avoid large weights at small multiplicity
333  sigma = ( sigma < 0.75 ? 0.75 : sigma);
334 
335  rank_bemc = (primV->numMatchesWithBEMC() - n_bemc_expected)/sigma;
336  if (mUseOldBEMCRank)
337  rank_bemc += 0.5; // distribution is asymmetric; add 0.5
338  }
339 
340  Float_t rank_cross = 0;
341  if ( n_cross_expected >= 1 ) {
342  Float_t sigma=1.1*sqrt(n_cross_expected);
343  rank_cross = (primV->numTracksCrossingCentralMembrane() - n_cross_expected)/sigma;
344  }
345 
346  // Handle possible overflows
347  rank_avg_dip = ( rank_avg_dip < -5 ? -5 : rank_avg_dip );
348  rank_bemc = ( rank_bemc < -5 ? -5 : (rank_bemc > 1 ? 1 : rank_bemc) );
349  rank_cross = ( rank_cross < -5 ? -5 : (rank_cross > 1 ? 1 : rank_cross) );
350 
351  if (mDebugLevel) {
352  LOG_INFO << "rankings: " << rank_avg_dip << " " << rank_bemc << " " << rank_cross << endm;
353  }
354 
355  //Give split vertices a lower rank...
356  if (mLowerSplitVtxRank && mCTBSum > 6000.0*primV->numTracksUsedInFinder()/80. + 2000)
357  primV->setRanking(rank_cross+rank_bemc+rank_avg_dip-3);
358  else
359  primV->setRanking(rank_cross+rank_bemc+rank_avg_dip);
360 
361  if (primV->ranking() > mBestRank) {
362  mBestRank = primV->ranking();
363  mBestVtx = primV;
364  }
365  }
366 }
367 
368 
369 int StMinuitVertexFinder::fit(StEvent* event)
370 {
371  setFlagBase();
372 
373  // get CTB info
374  StCtbTriggerDetector* ctbDet = 0;
375  std::vector<ctbHit> ctbHits;
376 
377  StTriggerDetectorCollection* trigCol = event->triggerDetectorCollection();
378  mCTBSum = 0;
379  if(trigCol){
380  ctbDet = &(trigCol->ctb());
381 
382  for (UInt_t slat = 0; slat < ctbDet->numberOfSlats(); slat++) {
383  for (UInt_t tray = 0; tray < ctbDet->numberOfTrays(); tray++) {
384  ctbHit curHit;
385  curHit.adc = ctbDet->mips(tray,slat,0);
386  if(curHit.adc > 0){
387  mCTBSum += curHit.adc;
388  ctb_get_slat_from_data(slat, tray, curHit.phi, curHit.eta);
389  ctbHits.push_back(curHit);
390  }
391  }
392  }
393  }
394 
395  // Loop over all global tracks (TPC) and store the refering helices and
396  // their estimated DCA resolution in vectors. Quality cuts are applied (see
397  // StMinuitVertexFinder::accept()). The helices and the sigma are used in
398  // fcn to calculate the fit potential which gets minimized by Minuit.
399  mDCAs.clear();
400  mHelices.clear();
401  mHelixFlags.clear();
402  mZImpact.clear();
403 
404  fillBemcHits(event);
405 
406  Int_t n_ctb_match_tot = 0;
407  Int_t n_bemc_match_tot = 0;
408  Int_t n_cross_tot = 0;
409 
410  for (const StTrackNode* stTrack : event->trackNodes())
411  {
412  StGlobalTrack* g = ( StGlobalTrack*) stTrack->track(global);
413  if (!accept(g)) continue;
414  StDcaGeometry* gDCA = g->dcaGeometry();
415  if (! gDCA) continue;
416  if (TMath::Abs(gDCA->impact()) > mRImpactMax) continue;
417  mDCAs.push_back(gDCA);
418  // StPhysicalHelixD helix = gDCA->helix();
419  // mHelices.push_back(helix);
420  mHelices.push_back(g->geometry()->helix());
421  mHelixFlags.push_back(1);
422  Double_t z_lin = gDCA->z();
423  mZImpact.push_back(z_lin);
424 
425  Bool_t shouldHitCTB = kFALSE;
426  Double_t etaInCTBFrame = -999;
427  bool ctb_match = EtaAndPhiToOrriginAtCTB(g,&ctbHits,shouldHitCTB,etaInCTBFrame);
428  if (ctb_match) {
429  mHelixFlags[mHelixFlags.size()-1] |= kFlagCTBMatch;
430  n_ctb_match_tot++;
431  }
432 
433  if (matchTrack2BEMC(g)) {
434  mHelixFlags[mHelixFlags.size()-1] |= kFlagBEMCMatch;
435  n_bemc_match_tot++;
436  }
437 
438  if (checkCrossMembrane(g)) {
439  mHelixFlags[mHelixFlags.size()-1] |= kFlagCrossMembrane;
440  n_cross_tot++;
441  }
442  }
443 
444  if (mDebugLevel) {
445  LOG_INFO << "Found " << n_ctb_match_tot << " ctb matches, " << n_bemc_match_tot << " bemc matches, " << n_cross_tot << " tracks crossing central membrane" << endm;
446  }
447 
448  // In case there are no tracks left we better quit
449  if (mHelices.empty()) {
450  LOG_WARN << "StMinuitVertexFinder::fit: no tracks to fit." << endm;
451  mStatusMin = -1;
452  return 0;
453  }
454 
455  LOG_INFO << "StMinuitVertexFinder::fit size of helix vector: " << mHelices.size() << endm;
456 
457  // Set some global pars
458  if (mRequireCTB) requireCTB = kTRUE;
459 
460  // Reset and clear Minuit parameters mStatusMin
461  mMinuit->mnexcm("CLEar", 0, 0, mStatusMin);
462 
463  // Make sure the global pointer points to valid object so Minuit uses correct data
465 
466  // Set parameters and start values. We do constrain the parameters since it
467  // harms the fit quality (see Minuit documentation).
468  //
469  // Initialize the seed with a z value which is not one of the discrete
470  // values which it can tend to, implies zero not allowed. Also need
471  // different initialization when vertex constraint.
472 
473  static Double_t step[3] = {0.03, 0.03, 0.03};
474 
475  // Scan z to find best seed for the actual fit.
476  // Skip this step if an external seed is given.
477  if (!mExternalSeedPresent) {
478  findSeeds();
479  }
480  else {
481  mNSeed = 1;
482  }
483 
484  Float_t old_vtx_z = -999;
485  Double_t seed_z = -999;
486  Double_t chisquare = 0;
487 
488  for (Int_t iSeed = 0; iSeed < mNSeed; iSeed++) {
489 
490  // Reset and clear Minuit parameters mStatusMin
491  mMinuit->mnexcm("CLEar", 0, 0, mStatusMin);
492 
493  seed_z= mSeedZ[iSeed];
494 
495  if (mExternalSeedPresent)
496  seed_z = mExternalSeed.z();
497  if (mDebugLevel) {
498  LOG_INFO << "Vertex seed = " << seed_z << endm;
499  }
500 
501  if (mVertexFitMode == VertexFit_t::Beamline1D){
502  mMinuit->mnparm(0, "z", seed_z, step[2], 0, 0, mStatusMin);
503  }
504  else {
505  mMinuit->mnparm(0, "x", 0, step[0], 0, 0, mStatusMin);
506  mMinuit->mnparm(1, "y", 0, step[1], 0, 0, mStatusMin);
507  mMinuit->mnparm(2, "z", seed_z, step[2], 0, 0, mStatusMin);
508  }
509 
510  Int_t done = 0;
511  Int_t iter = 0;
512 
513  Int_t n_trk_vtx = 0;
514  Int_t n_helix = mHelices.size();
515  do {
516  // For most vertices one pass is fine, but multiple passes can be done
517  n_trk_vtx = 0;
518  for (Int_t i=0; i < n_helix; i++) {
519  if (fabs(mZImpact[i]-seed_z) < mDcaZMax) {
520  mHelixFlags[i] |= kFlagDcaz;
521  n_trk_vtx++;
522  }
523  else
524  mHelixFlags[i] &= ~kFlagDcaz;
525  }
526 
527  if (mDebugLevel) {
528  LOG_INFO << n_trk_vtx << " tracks within dcaZ cut (iter " << iter <<" )" << endm;
529  }
530  if (n_trk_vtx < mMinTrack) {
531  if (mDebugLevel) {
532  LOG_INFO << "Less than mMinTrack (=" << mMinTrack << ") tracks, skipping vtx" << endm;
533  }
534  continue;
535  }
536  mMinuit->mnexcm("MINImize", 0, 0, mStatusMin);
537  done = 1;
538 
539  // Check fit result
540  if (mStatusMin) {
541  LOG_WARN << "StMinuitVertexFinder::fit: error in Minuit::mnexcm(), check status flag. ( iter=" << iter << ")" << endm;
542  done = 0; // refit
543  }
544 
545  Double_t fedm, errdef;
546  Int_t npari, nparx;
547 
548  mMinuit->mnstat(chisquare, fedm, errdef, npari, nparx, mStatusMin);
549 
550  if (mStatusMin != 3) {
551  LOG_INFO << "Warning: Minuit Status: " << mStatusMin << ", func val " << chisquare<< endm;
552  done = 0; // refit
553  }
554  mMinuit->mnhess();
555 
556  Double_t new_z, zerr;
557  if (mVertexFitMode == VertexFit_t::Beamline1D) {
558  mMinuit->GetParameter(0, new_z, zerr);
559  }
560  else {
561  mMinuit->GetParameter(2, new_z, zerr);
562  }
563 
564  if (fabs(new_z - seed_z) > 1) // refit if vertex shifted
565  done = 0;
566 
567  Int_t n_trk = 0;
568  for (Int_t i=0; i < n_helix; i++) {
569  if ( fabs(mZImpact[i] - new_z) < mDcaZMax ) {
570  n_trk++;
571  }
572  }
573  if ( 10 * abs(n_trk - n_trk_vtx) >= n_trk_vtx ) // refit if number of track changed by more than 10%
574  done = 0;
575 
576  iter++;
577  seed_z = new_z; // seed for next iteration
578  } while (!done && iter < 5 && n_trk_vtx >= mMinTrack);
579 
580  if (n_trk_vtx < mMinTrack)
581  continue;
582 
583  if (!done) {
584  LOG_WARN << "Vertex unstable: no convergence after " << iter << " iterations. Skipping vertex " << endm;
585  continue;
586  }
587 
588  if (!mExternalSeedPresent && fabs(seed_z-mSeedZ[iSeed]) > mDcaZMax) {
589  LOG_WARN << "Vertex walks during fits: seed was " << mSeedZ[iSeed] << ", vertex at " << seed_z << endm;
590  }
591 
592  if (fabs(seed_z - old_vtx_z) < mDcaZMax) {
593  if (mDebugLevel) {
594  LOG_INFO << "Vertices too close (<mDcaZMax). Skipping" << endm;
595  }
596  continue;
597  }
598 
599  // Store vertex
600  StThreeVectorD XVertex;
601  Float_t cov[6];
602  memset(cov,0,sizeof(cov));
603 
604  Double_t val, verr;
605  if (mVertexFitMode == VertexFit_t::Beamline1D) {
606  mMinuit->GetParameter(0, val, verr);
607  XVertex.setZ(val); cov[5]=verr*verr;
608 
609  // LSB Really error in x and y should come from error on constraint
610  // At least this way it is clear that those were fixed paramters
611  XVertex.setX(beamX(val)); cov[0]=0.1; // non-zero error values needed for Sti
612  XVertex.setY(beamY(val)); cov[2]=0.1;
613  }
614  else {
615  XVertex = StThreeVectorD(mMinuit->fU[0],mMinuit->fU[1],mMinuit->fU[2]);
616  Double_t emat[9];
617  /* 0 1 2
618  3 4 5
619  6 7 8 */
620  mMinuit->mnemat(emat,3);
621  cov[0] = emat[0];
622  cov[1] = emat[3];
623  cov[2] = emat[4];
624  cov[3] = emat[6];
625  cov[4] = emat[7];
626  cov[5] = emat[8];
627  }
628  StPrimaryVertex primV;
629  primV.setPosition(XVertex);
630  primV.setChiSquared(chisquare); // this is not really a chisquare, but anyways
631  primV.setCovariantMatrix(cov);
632  primV.setVertexFinderId(minuitVertexFinder);
633  primV.setFlag(1); // was not set earlier by this vertex finder ?? Jan
634  primV.setRanking(333);
635  primV.setNumTracksUsedInFinder(n_trk_vtx);
636 
637  Int_t n_ctb_match = 0;
638  Int_t n_bemc_match = 0;
639  Int_t n_cross = 0;
640  n_trk_vtx = 0;
641 
642  Double_t mean_dip = 0;
643  Double_t sum_pt = 0;
644  for (UInt_t i=0; i < mHelixFlags.size(); i++) {
645  if (!(mHelixFlags[i] & kFlagDcaz))
646  continue;
647  n_trk_vtx++;
648  if (mHelixFlags[i] & kFlagCTBMatch)
649  n_ctb_match++;
650  if (mHelixFlags[i] & kFlagBEMCMatch)
651  n_bemc_match++;
652  if (mHelixFlags[i] & kFlagCrossMembrane)
653  n_cross++;
654  mean_dip += mHelices[i].dipAngle();
655  sum_pt += mHelices[i].momentum(0).perp();
656  }
657 
658  mean_dip /= n_trk_vtx;
659 
660  if (mDebugLevel) {
661  LOG_INFO << "check n_trk_vtx " << n_trk_vtx << ", found "
662  << n_ctb_match << " ctb matches, "
663  << n_bemc_match << " bemc matches, "
664  << n_cross << " tracks crossing central membrane\n"
665  << "mean dip " << mean_dip << endm;
666  }
667  primV.setNumMatchesWithCTB(n_ctb_match);
668  primV.setNumMatchesWithBEMC(n_bemc_match);
669  primV.setNumTracksCrossingCentralMembrane(n_cross);
670  primV.setMeanDip(mean_dip);
671  primV.setSumOfTrackPt(sum_pt);
672 
673  //..... add vertex to the list
674  addVertex(primV);
675 
676  old_vtx_z = XVertex.z();
677  }
678 
679  calculateRanks();
680 
681  // Reset the flag which tells us about external
682  // seeds. This needs to be provided for every event.
683  mExternalSeedPresent = kFALSE;
684 
685  requireCTB = kFALSE;
686 
687  return 1;
688 }
689 
690 
691 //________________________________________________________________________________
692 double StMinuitVertexFinder::CalcChi2DCAs(const StThreeVectorD &vtx) {
693  Double_t f = 0;
694  Double_t e;
695  nCTBHits = 0;
696  if (fabs(vtx.x())> 10) return 1e6;
697  if (fabs(vtx.y())> 10) return 1e6;
698  if (fabs(vtx.z())>300) return 1e6;
699  for (UInt_t i=0; i<mDCAs.size(); i++) {
700 
701  if ( !(mHelixFlags[i] & kFlagDcaz) || (requireCTB && !(mHelixFlags[i] & kFlagCTBMatch)) )
702  continue;
703 
704  const StDcaGeometry* gDCA = mDCAs[i];
705  if (! gDCA) continue;
706  const StPhysicalHelixD helix = gDCA->helix();
707  e = helix.distance(vtx, kFALSE); // false: don't do multiple loops
708  //VP version
709  //VP Double_t chi2 = e*e/(errMatrix[0] + errMatrix[2]);
710  Double_t err2;
711  Double_t chi2 = gDCA->thelix().Dca(&(vtx.x()),&err2);
712  chi2*=chi2/err2;
713  //EndVP
714  static double scale = 100;
715  f += scale*(1. - TMath::Exp(-chi2/scale)); // robust potential
716  // f -= scale*TMath::Exp(-chi2/scale); // robust potential
717  if((mHelixFlags[i] & kFlagCTBMatch) && e<3.0) nCTBHits++;
718  }
719  return f;
720 }
721 
722 
726 bool StMinuitVertexFinder::accept(StTrack* track) const
727 {
728  return (track &&
729  track->flag() >= 0 &&
730  track->fitTraits().numberOfFitPoints() >= mMinNumberOfFitPointsOnTrack &&
731  !track->topologyMap().trackFtpc() &&
732  finite(track->length()) && //LSB another temporary check
733  track->geometry()->helix().valid());
734 }
735 
736 
739 {
740  mMinuit->SetPrintLevel(level);
741 }
742 
743 
744 void StMinuitVertexFinder::printInfo(ostream& os) const
745 {
746  os << "StMinuitVertexFinder - Statistics:" << endl;
747  os << "Number of vertices found ......." << size() << endl;
748  os << "Rank of best vertex ............" << mBestRank << endl;
749  if(mBestVtx) {
750  os << "Properties of best vertex:" << endl;
751  os << "position ..................... " << mBestVtx->position() << endl;
752  os << "position errors .............. " << mBestVtx->positionError()<< endl;
753  os << "# of used tracks ............. " << mBestVtx->numTracksUsedInFinder() << endl;
754  os << "Chisquare .................... " << mBestVtx->chiSquared() << endl;
755  }
756  os << "min # of fit points for tracks . " << mMinNumberOfFitPointsOnTrack << endl;
757 }
758 
759 
760 Int_t StMinuitVertexFinder::NCtbMatches() {
761  return nCTBHits;
762 }
bool valid(double world=1.e+5) const
checks for valid parametrization
Definition: StHelix.hh:128
double beamX(double z) const
double beamY(double z) const
double distance(const StThreeVector< double > &p, bool scanPeriods=true) const
minimal distance between point and helix
Definition: StHelix.cc:240
pair< double, double > pathLength(double r) const
path length at given r (cylindrical r)
Definition: StHelix.cc:351
void setPrintLevel(Int_t=0)
Use mMinuit print level.
VertexFit_t mVertexFitMode
The type of vertex fit to use in derived concrete implementation.
static StGenericVertexFinder * sSelf
By default point to invalid object.
Int_t getBin(const Float_t phi, const Float_t eta, Int_t &m, Int_t &e, Int_t &s) const
Definition: StEmcGeom.h:321
double Dca(const double point[3], double *dcaErr=0) const
DCA to given space point (with error matrix)