BmnRoot
Loading...
Searching...
No Matches
SrcVertexFinder.cxx
Go to the documentation of this file.
1#include "SrcVertexFinder.h"
2
3#include "BmnMath.h"
4#include "TCanvas.h"
5#include "TFitResult.h"
6#include "TGraph.h"
7#include "vector"
8
9using namespace std;
10using namespace TMath;
11
12SrcVertexFinder::SrcVertexFinder(Int_t period, Bool_t isExp) {
13 fPeriodId = period;
14 fEventNo = 0;
15 fGlobalTracksArray = nullptr;
16 fGemHitsArray = nullptr;
17 fMwpcTracksArray = nullptr;
18 fTof400HitsArray = nullptr;
19 fVertexArray = nullptr;
20 fArmTracksArray = nullptr;
21 fTime = 0.0;
22 fisExp = isExp;
23}
24
26
28 if (fVerbose > 1)
29 cout << "=========================== Vertex finder init started ====================" << endl;
30
31 fKalman = new BmnKalmanFilter();
32
33 // Get ROOT Manager
34 FairRootManager *ioman = FairRootManager::Instance();
35 if (nullptr == ioman)
36 Fatal("Init", "FairRootManager is not instantiated");
37
38 fGlobalTracksArray = (TClonesArray *)ioman->GetObject("BmnGlobalTrack"); // in
39 if (!fGlobalTracksArray) {
40 cout << "SrcVertexFinder::Init(): branch BmnGlobalTrack not found! Task will be deactivated" << endl;
41 SetActive(kFALSE);
42 return kERROR;
43 }
44 fGemHitsArray = (TClonesArray *)ioman->GetObject("BmnGemStripHit"); // in
45 if (!fGemHitsArray) {
46 cout << "SrcVertexFinder::Init(): branch BmnGemStripHit not found! Task will be deactivated" << endl;
47 SetActive(kFALSE);
48 return kERROR;
49 }
50 fTof400HitsArray = (TClonesArray *)ioman->GetObject("BmnTof400Hit"); // in
51 if (!fTof400HitsArray) {
52 cout << "SrcVertexFinder::Init(): branch BmnTof400Hit not found! Task will be deactivated" << endl;
53 SetActive(kFALSE);
54 return kERROR;
55 }
56 fMwpcTracksArray = (TClonesArray *)ioman->GetObject("BmnMwpcTrack"); // in
57 if (!fMwpcTracksArray) {
58 cout << "SrcVertexFinder::Init(): branch BmnMwpcTrack not found! Task will be deactivated" << endl;
59 SetActive(kFALSE);
60 return kERROR;
61 }
62
63 fVertexArray = new TClonesArray("BmnVertex", 1); // out
64 ioman->Register("SrcVertex", "VERTEX", fVertexArray, kTRUE);
65 fArmTracksArray = new TClonesArray("BmnTrack", 1); // out
66 ioman->Register("SrcArmTrack", "ARM", fArmTracksArray, kTRUE);
67
68 if (fVerbose > 1) cout << "=========================== Vertex finder init finished ===================" << endl;
69
70 return kSUCCESS;
71}
72
73void SrcVertexFinder::Exec(Option_t *opt) {
74 TStopwatch sw;
75 sw.Start();
76
77 if (!IsActive())
78 return;
79
80 if (fVerbose > 1) {
81 cout << "======================== Vertex finder exec started ======================" << endl;
82 cout << "Event number: " << fEventNo++ << endl;
83 }
84
85 fVertexArray->Delete();
86 fArmTracksArray->Delete();
87
88 vector<BmnTrack> lTracks;
89 vector<BmnTrack> rTracks;
90 CreateArmCandidates(lTracks, rTracks);
91
92 Int_t nTrWithUpstream = 0;
93 for (Int_t iTrack = 0; iTrack < fGlobalTracksArray->GetEntriesFast(); iTrack++) {
94 BmnGlobalTrack *track = (BmnGlobalTrack *)fGlobalTracksArray->UncheckedAt(iTrack);
95 if (track->GetUpstreamTrackIndex() == -1)
96 track->SetFlag(13);
97 else
98 nTrWithUpstream++;
99 }
100
101 if (lTracks.size() > 100 || rTracks.size() > 100 || nTrWithUpstream > 10) { //Check this condition!
102 new ((*fVertexArray)[fVertexArray->GetEntriesFast()]) BmnVertex();
103 if (fVerbose > 0) cout << "SrcVertexFinder: Vertex NOT found" << endl;
104 } else {
105 //FindVertexByVirtualPlanes(lTracks, rTracks);
106 FindVertexAnalitically(lTracks, rTracks);
107 if (fVerbose > 0) {
108 BmnVertex *vert = (BmnVertex *)fVertexArray->At(0);
109 vert->Print();
110 }
111 }
112 lTracks.clear();
113 rTracks.clear();
114
115 if (fVerbose > 1)
116 cout << "\n======================== Vertex finder exec finished ======================" << endl;
117
118 sw.Stop();
119 fTime += sw.RealTime();
120}
121
122void SrcVertexFinder::FindVertexAnalitically(vector<BmnTrack> &lTracks, vector<BmnTrack> &rTracks) {
123 Double_t minVZ = 1000;
124 Double_t vy = 1000;
125 Double_t vx = 1000;
126 Double_t minDist = 1000;
127 vector<BmnTrack> trackCombination;
128 vector<BmnTrack> bestCombination;
129 if (lTracks.size() != 0 && rTracks.size() != 0) {
130 minDist = 1000;
131 minVZ = 1000;
132 for (auto lTr : lTracks) {
133 for (auto rTr : rTracks) {
134 trackCombination.push_back(lTr);
135 trackCombination.push_back(rTr);
136 Double_t dist;
137 Double_t vz = GetVzByVectorStraightTracks(trackCombination, dist);
138 if (dist < minDist) {
139 minDist = dist;
140 minVZ = vz;
141 bestCombination = trackCombination;
142 }
143 trackCombination.clear();
144 }
145 }
146 } else if (lTracks.size() != 0 && rTracks.size() == 0) {
147 for (auto lTr : lTracks) {
148 for (Int_t iGl = 0; iGl < fGlobalTracksArray->GetEntriesFast(); ++iGl) {
149 BmnGlobalTrack *gl = (BmnGlobalTrack *)fGlobalTracksArray->At(iGl);
150 if (gl->GetUpstreamTrackIndex() == -1) continue;
151 trackCombination.push_back(*((BmnTrack *)fGlobalTracksArray->At(iGl)));
152 }
153 for (Int_t iM = 0; iM < fMwpcTracksArray->GetEntriesFast(); ++iM) {
154 BmnTrack *mTr = (BmnTrack *)fMwpcTracksArray->At(iM);
155 if (mTr->GetParamFirst()->GetZ() < -700) trackCombination.push_back(*mTr);
156 }
157 trackCombination.push_back(lTr);
158 Double_t dist;
159 Double_t vz = GetVzByVectorStraightTracks(trackCombination, dist);
160 if (dist < minDist) {
161 minDist = dist;
162 minVZ = vz;
163 bestCombination = trackCombination;
164 }
165 trackCombination.clear();
166 }
167 } else if (rTracks.size() != 0 && lTracks.size() == 0) {
168 for (auto rTr : rTracks) {
169 for (Int_t iGl = 0; iGl < fGlobalTracksArray->GetEntriesFast(); ++iGl) {
170 BmnGlobalTrack *gl = (BmnGlobalTrack *)fGlobalTracksArray->At(iGl);
171 if (gl->GetUpstreamTrackIndex() == -1) continue;
172 trackCombination.push_back(*((BmnTrack *)fGlobalTracksArray->At(iGl)));
173 }
174 for (Int_t iM = 0; iM < fMwpcTracksArray->GetEntriesFast(); ++iM) {
175 BmnTrack *mTr = (BmnTrack *)fMwpcTracksArray->At(iM);
176 if (mTr->GetParamFirst()->GetZ() < -700) trackCombination.push_back(*mTr);
177 }
178 trackCombination.push_back(rTr);
179 Double_t dist;
180 Double_t vz = GetVzByVectorStraightTracks(trackCombination, dist);
181 if (dist < minDist) {
182 minDist = dist;
183 minVZ = vz;
184 bestCombination = trackCombination;
185 }
186 trackCombination.clear();
187 }
188 } else if (rTracks.size() == 0 && lTracks.size() == 0) {
189 for (Int_t iGl = 0; iGl < fGlobalTracksArray->GetEntriesFast(); ++iGl) {
190 BmnGlobalTrack *gl = (BmnGlobalTrack *)fGlobalTracksArray->At(iGl);
191 if (gl->GetUpstreamTrackIndex() == -1) continue;
192 trackCombination.push_back(*((BmnTrack *)fGlobalTracksArray->At(iGl)));
193 }
194 for (Int_t iM = 0; iM < fMwpcTracksArray->GetEntriesFast(); ++iM) {
195 BmnTrack *mTr = (BmnTrack *)fMwpcTracksArray->At(iM);
196 if (mTr->GetParamFirst()->GetZ() < -700) trackCombination.push_back(*mTr);
197 }
198 Double_t dist;
199 Double_t vz = GetVzByVectorStraightTracks(trackCombination, dist);
200 if (dist < minDist) {
201 minDist = dist;
202 minVZ = vz;
203 bestCombination = trackCombination;
204 }
205 trackCombination.clear();
206 } else
207 new ((*fVertexArray)[fVertexArray->GetEntriesFast()]) BmnVertex();
208
209 if (minDist > 10.0) {
210 new ((*fVertexArray)[fVertexArray->GetEntriesFast()]) BmnVertex();
211 } else {
212 for (Int_t iM = 0; iM < fMwpcTracksArray->GetEntriesFast(); ++iM) {
213 BmnTrack *mTr = (BmnTrack *)fMwpcTracksArray->At(iM);
214 if (mTr->GetParamFirst()->GetZ() > -700) continue; //we need incoming beam
215 Double_t ax = mTr->GetParamFirst()->GetTx();
216 Double_t ay = mTr->GetParamFirst()->GetTy();
217 Double_t x = mTr->GetParamFirst()->GetX();
218 Double_t y = mTr->GetParamFirst()->GetY();
219 Double_t z = mTr->GetParamFirst()->GetZ();
220 Double_t bx = x - ax * z;
221 Double_t by = y - ay * z;
222 vx = ax * minVZ + bx;
223 vy = ay * minVZ + by;
224 }
225
226 Bool_t isLeft = kFALSE;
227 Bool_t isRight = kFALSE;
228 Bool_t isFragm = kFALSE;
229 Bool_t isBeam = kFALSE;
230 for (size_t iTr = 0; iTr < bestCombination.size(); ++iTr)
231 {
232 BmnTrack track = bestCombination[iTr];
233 if (track.GetFlag() == 13) continue;
234 if (track.GetParamFirst()->GetZ() > -700 && track.GetParamFirst()->GetZ() < -400 && track.GetParamFirst()->GetX() < 0) { //right arm track
235 track.SetNHits(2);
236 new ((*fArmTracksArray)[fArmTracksArray->GetEntriesFast()]) BmnTrack(track);
237 isRight = kTRUE;
238 }
239 if (track.GetParamFirst()->GetZ() > -700 && track.GetParamFirst()->GetZ() < -400 && track.GetParamFirst()->GetX() > 0) { //left arm track
240 track.SetNHits(2);
241 new ((*fArmTracksArray)[fArmTracksArray->GetEntriesFast()]) BmnTrack(track);
242 isLeft = kTRUE;
243 }
244 if (track.GetParamFirst()->GetZ() > -400) { //fragment track
245 isFragm = kTRUE;
246 }
247 if (track.GetParamFirst()->GetZ() < -700) { //fragment track
248 isBeam = kTRUE;
249 }
250 }
251
252 Int_t type = (isRight && isLeft) ? 10 : (isLeft && isFragm && isBeam) ? 11 : (isRight && isFragm && isBeam) ? 12 : (isBeam && isFragm) ? 13 : -1;
253
254 vector<Int_t> idx;
255 for (Int_t iGl = 0; iGl < fGlobalTracksArray->GetEntriesFast(); ++iGl) {
256 BmnTrack gl = *((BmnTrack *)fGlobalTracksArray->At(iGl));
257 if (gl.GetFlag() != 13) idx.push_back(iGl);
258 }
259
260 new ((*fVertexArray)[fVertexArray->GetEntriesFast()]) BmnVertex(vx, vy, minVZ, minDist, 0, bestCombination.size(), TMatrixFSym(3), type, idx);
261 }
262}
263
264void SrcVertexFinder::FindVertexByVirtualPlanes(vector<BmnTrack> &lTracks, vector<BmnTrack> &rTracks) {
265 Float_t minDist = DBL_MAX; //minimal distance beetween tracks in point of vertex
266 Float_t minVZ = DBL_MAX; // VZ for tracks with minimal distance
267 vector<BmnTrack> trackCombination;
268 vector<BmnTrack> bestCombination;
269 if (lTracks.size() > 0 && rTracks.size() > 0) {
270 for (auto lTr : lTracks) {
271 for (auto rTr : rTracks) {
272 for (Int_t iGl = 0; iGl < fGlobalTracksArray->GetEntriesFast(); ++iGl) {
273 BmnGlobalTrack gl = *((BmnGlobalTrack *)fGlobalTracksArray->At(iGl));
274 if (gl.GetUpstreamTrackIndex() != -1)
275 trackCombination.push_back(*((BmnTrack *)fGlobalTracksArray->At(iGl)));
276 }
277 trackCombination.push_back(lTr);
278 trackCombination.push_back(rTr);
279 Float_t dist;
280 Float_t vz = FindVZByVirtualPlanes(-647, 100, trackCombination, dist);
281 if (vz < -999) {
282 trackCombination.clear();
283 continue;
284 }
285 if (dist < minDist) {
286 minDist = dist;
287 minVZ = vz;
288 bestCombination = trackCombination;
289 }
290 trackCombination.clear();
291 }
292 }
293 } else if (lTracks.size() == 0 && rTracks.size() > 0) {
294 for (auto rTr : rTracks) {
295 for (Int_t iGl = 0; iGl < fGlobalTracksArray->GetEntriesFast(); ++iGl) {
296 BmnGlobalTrack gl = *((BmnGlobalTrack *)fGlobalTracksArray->At(iGl));
297 if (gl.GetUpstreamTrackIndex() != -1)
298 trackCombination.push_back(*((BmnTrack *)fGlobalTracksArray->At(iGl)));
299 }
300 trackCombination.push_back(rTr);
301 Float_t dist;
302 Float_t vz = FindVZByVirtualPlanes(-647, 100, trackCombination, dist);
303 if (vz < -999) {
304 trackCombination.clear();
305 continue;
306 }
307 if (dist < minDist) {
308 minDist = dist;
309 minVZ = vz;
310 bestCombination = trackCombination;
311 }
312 trackCombination.clear();
313 }
314 } else if (lTracks.size() > 0 && rTracks.size() == 0) {
315 for (auto lTr : lTracks) {
316 for (Int_t iGl = 0; iGl < fGlobalTracksArray->GetEntriesFast(); ++iGl) {
317 BmnGlobalTrack gl = *((BmnGlobalTrack *)fGlobalTracksArray->At(iGl));
318 if (gl.GetUpstreamTrackIndex() != -1)
319 trackCombination.push_back(*((BmnTrack *)fGlobalTracksArray->At(iGl)));
320 }
321 trackCombination.push_back(lTr);
322 Float_t dist;
323 Float_t vz = FindVZByVirtualPlanes(-647, 100, trackCombination, dist);
324 if (vz < -999) {
325 trackCombination.clear();
326 continue;
327 }
328 if (dist < minDist) {
329 minDist = dist;
330 minVZ = vz;
331 bestCombination = trackCombination;
332 }
333 trackCombination.clear();
334 }
335 }
336 if (minDist > 10.0) {
337 new ((*fVertexArray)[fVertexArray->GetEntriesFast()]) BmnVertex();
338 } else {
339 vector<Double_t> xHits;
340 vector<Double_t> yHits;
341 for (size_t iTr = 0; iTr < bestCombination.size(); ++iTr)
342 {
343 BmnTrack *track = &bestCombination[iTr];
344 FairTrackParam par0 = (track->GetParamFirst()->GetZ() < -400) ? (*(track->GetParamLast())) : (*(track->GetParamFirst()));
345 Double_t len = track->GetLength();
346 if (fKalman->TGeoTrackPropagate(&par0, minVZ, 2212, nullptr, &len) == kBMNERROR) {
347 track->SetFlag(13);
348 continue;
349 }
350 track->SetLength(len);
351 track->SetB(len / track->GetB() / (TMath::C() * 1e-7)); //for arm tracks B contains beta. It is stupid, but simple way to store it...
352 xHits.push_back(par0.GetX());
353 yHits.push_back(par0.GetY());
354 }
355 if (xHits.size() < 2) {
356 new ((*fVertexArray)[fVertexArray->GetEntriesFast()]) BmnVertex();
357 return;
358 }
359
360 Bool_t isLeft = kFALSE;
361 Bool_t isRight = kFALSE;
362 Bool_t isFragm = kFALSE;
363 for (size_t iTr = 0; iTr < bestCombination.size(); ++iTr)
364 {
365 BmnTrack track = bestCombination[iTr];
366 if (track.GetFlag() == 13) continue;
367 if (track.GetParamFirst()->GetZ() < -400 && track.GetParamFirst()->GetX() < 0) { //right arm track
368 track.SetNHits(2);
369 new ((*fArmTracksArray)[fArmTracksArray->GetEntriesFast()]) BmnTrack(track);
370 isRight = kTRUE;
371 }
372 if (track.GetParamFirst()->GetZ() < -400 && track.GetParamFirst()->GetX() > 0) { //left arm track
373 track.SetNHits(2);
374 new ((*fArmTracksArray)[fArmTracksArray->GetEntriesFast()]) BmnTrack(track);
375 isLeft = kTRUE;
376 }
377 if (track.GetParamFirst()->GetZ() > -400) { //fragment track
378
379 isFragm = kTRUE;
380 }
381 }
382
383 Int_t type = (isRight && isLeft && isFragm) ? 10 : (isRight && isLeft) ? 11 : (isLeft && isFragm) ? 12 : 13;
384
385 Double_t vz = minVZ;
386 Double_t vx = Mean(xHits.begin(), xHits.end());
387 Double_t vy = Mean(yHits.begin(), yHits.end());
388
389 vector<Int_t> idx;
390 for (Int_t iGl = 0; iGl < fGlobalTracksArray->GetEntriesFast(); ++iGl) {
391 BmnTrack gl = *((BmnTrack *)fGlobalTracksArray->At(iGl));
392 if (gl.GetFlag() != 13) idx.push_back(iGl);
393 }
394
395 new ((*fVertexArray)[fVertexArray->GetEntriesFast()]) BmnVertex(vx, vy, vz, minDist, 0, xHits.size(), TMatrixFSym(3), type, idx);
396 }
397}
398
400 ofstream outFile;
401 outFile.open("QA/timing.txt", ofstream::app);
402 outFile << "Vertex Finder Time: " << fTime;
403 if (fVerbose == 0)
404 cout << "Work time of the GEM vertex finder: " << fTime << endl;
405}
406
407Double_t SrcVertexFinder::CalcMeanDist(vector<Double_t> x, vector<Double_t> y) {
408 Double_t sumDist = 0.0;
409 Int_t nPairs = 0;
410 for (size_t i = 0; i < x.size(); ++i) {
411 for (size_t j = i + 1; j < x.size(); ++j) {
412 sumDist += Sqrt(Sq(x[i] - x[j]) + Sq(y[i] - y[j]));
413 nPairs++;
414 }
415 }
416 return sumDist / nPairs; // calc. ave. dist value
417}
418
419Float_t SrcVertexFinder::FindVZByVirtualPlanes(Float_t z_0, Float_t range, vector<BmnTrack> tracks, Float_t &minDist) {
420 const Int_t nPlanes = 5;
421 Float_t minZ = z_0;
422
423 while (range >= 0.01) {
424 Float_t zMax = minZ + range;
425 Float_t zMin = minZ - range;
426 Float_t zStep = (zMax - zMin) / (nPlanes - 1);
427
428 vector<Double_t> xHits[nPlanes];
429 vector<Double_t> yHits[nPlanes];
430 Float_t zPlane[nPlanes];
431 Float_t rRMS[nPlanes] = {0};
432
433 for (Int_t iPlane = 0; iPlane < nPlanes; ++iPlane)
434 zPlane[iPlane] = zMax - iPlane * zStep;
435
436 Int_t nOkTr = 0;
437
438 for (size_t iTr = 0; iTr < tracks.size(); ++iTr)
439 {
440 BmnTrack track = tracks[iTr];
441 FairTrackParam par0 = *(track.GetParamFirst());
442 Double_t xTr[nPlanes];
443 Double_t yTr[nPlanes];
444 Bool_t trOk = kTRUE;
445 for (Int_t iPlane = 0; iPlane < nPlanes; ++iPlane) {
446 if (fKalman->TGeoTrackPropagate(&par0, zPlane[iPlane], 2212, nullptr, nullptr) == kBMNERROR) {
447 trOk = kFALSE;
448 break;
449 }
450
451 //cout << " " << par0.GetX() << " " << par0.GetY() << " " << par0.GetZ() << endl;
452 xTr[iPlane] = par0.GetX();
453 yTr[iPlane] = par0.GetY();
454 //xHits[iPlane].push_back(par0.GetX());
455 //yHits[iPlane].push_back(par0.GetY());
456 }
457
458 if (trOk) {
459 nOkTr++;
460 for (Int_t iPlane = 0; iPlane < nPlanes; ++iPlane) {
461 xHits[iPlane].push_back(xTr[iPlane]);
462 yHits[iPlane].push_back(yTr[iPlane]);
463 }
464 }
465 }
466
467 if (nOkTr < 2) {
468 minDist = DBL_MAX;
469 return -1000.0;
470 }
471
472 //Calculation minZ as minimum of parabola
473 for (Int_t iPlane = 0; iPlane < nPlanes; ++iPlane) {
474 rRMS[iPlane] = CalcMeanDist(xHits[iPlane], yHits[iPlane]);
475 }
476 TGraph *vertex = new TGraph(nPlanes, zPlane, rRMS);
477 TFitResultPtr ptr = vertex->Fit("pol2", "QFS");
478 Float_t c = ptr->Parameter(0);
479 Float_t b = ptr->Parameter(1);
480 Float_t a = ptr->Parameter(2);
481 minZ = -b / (2 * a);
482 minDist = a * minZ * minZ + b * minZ + c;
483 delete vertex;
484
485 range /= 2;
486 }
487 return minZ;
488}
489
490void SrcVertexFinder::CreateArmCandidates(vector<BmnTrack> &lTracks, vector<BmnTrack> &rTracks) {
491 vector<Int_t> lTofHitIdx;
492 vector<Int_t> rTofHitIdx;
493 for (Int_t iTof = 0; iTof < fTof400HitsArray->GetEntriesFast(); ++iTof) {
494 BmnHit *tHit = (BmnHit *)fTof400HitsArray->At(iTof);
495 if (tHit->GetX() > 0)
496 lTofHitIdx.push_back(iTof);
497 else
498 rTofHitIdx.push_back(iTof);
499 }
500 //Maybe we don't need next condition? Check it!
501 //if (lTofHitIdx.size() == 0 || rTofHitIdx.size() == 0) return; //one track must be on the left arm, another on the right
502
503 vector<Int_t> lGemHitIdx;
504 vector<Int_t> rGemHitIdx;
505 for (Int_t iGem = 0; iGem < fGemHitsArray->GetEntriesFast(); ++iGem) {
506 BmnHit *gHit = (BmnHit *)fGemHitsArray->At(iGem);
507 if (gHit->GetStation() > 3) continue;
508 if (gHit->GetX() > 0)
509 lGemHitIdx.push_back(iGem);
510 else
511 rGemHitIdx.push_back(iGem);
512 }
513 //Maybe we don't need next condition? Check it!
514 //if (lGemHitIdx.size() == 0 || rGemHitIdx.size() == 0) return; //one track must be on the left arm, another on the right
515
516 BmnHit tHit, gHit;
517
518 //Coeficients for arms alignment
519 Double_t dxall, dyall;
520 Double_t dy_tl, dy_tr;
521 Double_t gx, gy, gz, tx, ty, tz;
522 if (fisExp) {
523 dxall = -2.113;
524 dyall = 2.504;
525 dy_tl = -1.72;
526 dy_tr = 1.3;
527 } else {
528 dxall = 0;
529 dyall = 0;
530 dy_tl = 0;
531 dy_tr = 0;
532 }
533
534 for (auto gIdx : lGemHitIdx) {
535 for (auto tIdx : lTofHitIdx) {
536 tHit = *((BmnHit *)fTof400HitsArray->At(tIdx));
537 gHit = *((BmnHit *)fGemHitsArray->At(gIdx));
538 BmnTrack lTr;
539 gx = gHit.GetX() + dxall;
540 gy = gHit.GetY() + dyall;
541 gz = gHit.GetZ();
542 tx = tHit.GetX() + dxall;
543 ty = tHit.GetY() + dyall + dy_tl;
544 tz = tHit.GetZ();
545
546 lTr.GetParamFirst()->SetX(gx);
547 lTr.GetParamFirst()->SetY(gy);
548 lTr.GetParamFirst()->SetZ(gz);
549 lTr.GetParamFirst()->SetTx((tx - gx) / (tz - gz));
550 lTr.GetParamFirst()->SetTy((ty - gy) / (tz - gz));
551 lTr.GetParamLast()->SetX(tx);
552 lTr.GetParamLast()->SetY(ty);
553 lTr.GetParamLast()->SetZ(tz);
554 lTr.GetParamLast()->SetTx((tx - gx) / (tz - gz));
555 lTr.GetParamLast()->SetTy((ty - gy) / (tz - gz));
556 lTr.SetB(tHit.GetTimeStamp()); //temporary
557 lTracks.push_back(lTr);
558 }
559 }
560 for (auto gIdx : rGemHitIdx) {
561 for (auto tIdx : rTofHitIdx) {
562 tHit = *((BmnHit *)fTof400HitsArray->At(tIdx));
563 gHit = *((BmnHit *)fGemHitsArray->At(gIdx));
564 BmnTrack rTr;
565 gx = gHit.GetX() + dxall;
566 gy = gHit.GetY() + dyall;
567 gz = gHit.GetZ();
568 tx = tHit.GetX() + dxall;
569 ty = tHit.GetY() + dyall + dy_tr;
570 tz = tHit.GetZ();
571
572 rTr.GetParamFirst()->SetX(gx);
573 rTr.GetParamFirst()->SetY(gy);
574 rTr.GetParamFirst()->SetZ(gz);
575 rTr.GetParamFirst()->SetTx((tx - gx) / (tz - gz));
576 rTr.GetParamFirst()->SetTy((ty - gy) / (tz - gz));
577 rTr.GetParamLast()->SetX(tx);
578 rTr.GetParamLast()->SetY(ty);
579 rTr.GetParamLast()->SetZ(tz);
580 rTr.GetParamLast()->SetTx((tx - gx) / (tz - gz));
581 rTr.GetParamLast()->SetTy((ty - gy) / (tz - gz));
582 rTr.SetB(tHit.GetTimeStamp()); //temporary
583 rTracks.push_back(rTr);
584 }
585 }
586}
vector< Double_t > dist(vector< Double_t > qp, Double_t mu)
Definition BmnMath.cxx:869
Double_t GetVzByVectorStraightTracks(vector< BmnTrack > tr, Double_t &dist)
Definition BmnMath.cxx:999
int i
Definition P4_F32vec4.h:22
@ kBMNERROR
Definition BmnEnums.h:26
Short_t GetStation() const
Definition BmnHit.h:45
BmnStatus TGeoTrackPropagate(FairTrackParam *par, Double_t zOut, Int_t pdg, vector< Double_t > *F, Double_t *length, Bool_t isField)
Float_t GetLength() const
Definition BmnTrack.h:68
void SetFlag(Int_t flag)
Definition BmnTrack.h:93
FairTrackParam * GetParamFirst()
Definition BmnTrack.h:72
Float_t GetB() const
Definition BmnTrack.h:64
void SetNHits(Int_t n)
Definition BmnTrack.h:105
void SetB(Double_t b)
Definition BmnTrack.h:109
FairTrackParam * GetParamLast()
Definition BmnTrack.h:76
void SetLength(Double_t length)
Definition BmnTrack.h:113
Int_t GetFlag() const
Definition BmnTrack.h:52
void Print()
Definition BmnVertex.cxx:56
void FindVertexAnalitically(vector< BmnTrack > &lTracks, vector< BmnTrack > &rTracks)
virtual ~SrcVertexFinder()
void FindVertexByVirtualPlanes(vector< BmnTrack > &lTracks, vector< BmnTrack > &rTracks)
void CreateArmCandidates(vector< BmnTrack > &lTracks, vector< BmnTrack > &rTracks)
virtual void Finish()
Float_t FindVZByVirtualPlanes(Float_t z_0, Float_t range, vector< BmnTrack > tracks, Float_t &minDist)
virtual InitStatus Init()
virtual void Exec(Option_t *opt)
STL namespace.