BmnRoot
Loading...
Searching...
No Matches
BmnVSPHitMaker.cxx
Go to the documentation of this file.
1#include "BmnVSPHitMaker.h"
2
3#include "BmnVSPModule.h"
4#include "FairRootManager.h"
5#include "TClonesArray.h"
6#include "TSystem.h"
7#include "UniDetectorParameter.h"
8
9#include <TStopwatch.h>
10
11static Double_t workTime = 0.0;
12
14 : fHitMatching(kTRUE)
15{
16
17 fInputPointsBranchName = "VSPPoint";
18 fInputDigitsBranchName = "BmnVSPDigit";
19 fInputDigitMatchesBranchName = "BmnVSPDigitMatch";
20
21 fOutputHitsBranchName = "BmnVSPHit";
22
23 fCurrentConfig = BmnVSPConfiguration::None;
24 StationSet = nullptr;
25 TransfSet = nullptr;
26}
27
28BmnVSPHitMaker::BmnVSPHitMaker(Int_t run_period, Int_t run_number, Bool_t isExp)
29 : fHitMatching(kTRUE)
30{
31 fIsExp = isExp;
32
33 fRunPeriod = run_period;
34 fRunNumber = run_number;
35
36 fInputPointsBranchName = "VSPPoint";
37 fInputDigitsBranchName = (!isExp) ? "BmnVSPDigit" : "VSP";
38 fInputDigitMatchesBranchName = "BmnVSPDigitMatch";
39
40 fOutputHitsBranchName = "BmnVSPHit";
41
42 fCurrentConfig = BmnVSPConfiguration::None;
43 StationSet = nullptr;
44 TransfSet = nullptr;
45
46 fSignalLow = 0.;
47 fSignalUp = DBL_MAX;
48
49 switch (fRunPeriod) {
50 case 9: // BM@N RUN-9
51 fCurrentConfig = BmnVSPConfiguration::Run9;
52 break;
53
54 default:
55 fCurrentConfig = BmnVSPConfiguration::None;
56 }
57}
58
59void BmnVSPHitMaker::createVSPDetector()
60{
61 TString gPathVSPConfig = gSystem->Getenv("VMCWORKDIR");
62 gPathVSPConfig += "/parameters/vsp/XMLConfigs/";
63
64 // Create VSP detector --------------------------------------------------
65 switch (fCurrentConfig) {
67 XMLConfigFile = gPathVSPConfig + "VSP_Run9.xml";
68 if (fVerbose)
69 cout << " Current VSP Configuration : Run9"
70 << "\n";
71 break;
72
73 default:;
74 }
75
76 if (!gSystem->AccessPathName(XMLConfigFile)) { // returns FALSE if the file exists
77 if (fVerbose)
78 cout << " XMLConfigFile : " << XMLConfigFile << "\n";
79 StationSet = new BmnVSPStationSet(XMLConfigFile);
80 TransfSet = new BmnVSPTransform();
81 TransfSet->LoadFromXMLFile(XMLConfigFile);
82 } else {
83 Fatal("BmnVSPHitMaker::createVSPDetector()", " !!! Current configuration is not set !!! ");
84 }
85
86 if (fIsExp) {
87 const Int_t nStat = StationSet->GetNStations();
88 UniDetectorParameter* coeffAlignCorrs =
89 UniDetectorParameter::GetDetectorParameter("VSP", "alignment_shift", fRunPeriod, fRunNumber);
90 vector<UniValue*> algnShifts;
91 if (coeffAlignCorrs)
92 coeffAlignCorrs->GetValue(algnShifts);
93 fAlignCor = new Double_t**[nStat];
94 for (Int_t iStat = 0; iStat < nStat; iStat++) {
95 const Int_t nParams = 3;
96 Int_t nModul = StationSet->GetVSPStation(iStat)->GetNModules();
97 fAlignCor[iStat] = new Double_t*[nModul];
98 for (Int_t iMod = 0; iMod < nModul; iMod++) {
99 fAlignCor[iStat][iMod] = new Double_t[nParams];
100 for (Int_t iPar = 0; iPar < nParams; iPar++) {
101 fAlignCor[iStat][iMod][iPar] = 0.0;
102 }
103 }
104 }
105
106 for (size_t i = 0; i < algnShifts.size(); ++i) {
107 // cout << ((AlignmentValue*)algnShifts[i])->value[0] << endl;
108 // cout << ((AlignmentValue*)algnShifts[i])->value[1] << endl;
109 // cout << ((AlignmentValue*)algnShifts[i])->value[2] << endl;
110 Int_t st = ((AlignmentValue*)algnShifts[i])->station;
111 Int_t mod = ((AlignmentValue*)algnShifts[i])->module;
112 fAlignCor[st][mod][0] = ((AlignmentValue*)algnShifts[i])->value[0];
113 fAlignCor[st][mod][1] = ((AlignmentValue*)algnShifts[i])->value[1];
114 fAlignCor[st][mod][2] = ((AlignmentValue*)algnShifts[i])->value[2];
115 }
116 }
117}
118
120{
121
122 if (fIsExp && fRunPeriod >= 9) {
123 for (Int_t iStat = 0; iStat < StationSet->GetNStations(); iStat++) {
124 for (Int_t iMod = 0; iMod < StationSet->GetVSPStation(iStat)->GetNModules(); iMod++)
125 delete[] fAlignCor[iStat][iMod];
126 delete[] fAlignCor[iStat];
127 }
128 delete[] fAlignCor;
129 }
130
131 if (StationSet) {
132 delete StationSet;
133 }
134
135 if (TransfSet) {
136 delete TransfSet;
137 }
138}
139
141{
142
143 if (fVerbose > 1)
144 cout << "=================== BmnVSPHitMaker::Init() started ================" << endl;
145
146 if (fRunPeriod < 9) {
147 SetActive(kFALSE);
148 return kERROR;
149 }
150
151 createVSPDetector();
152
153 FairRootManager* ioman = FairRootManager::Instance();
154
155 fBmnVSPDigitsArray = (TClonesArray*)ioman->GetObject(fInputDigitsBranchName);
156 if (!fBmnVSPDigitsArray) {
157 cout << "BmnVSPHitMaker::Init(): branch " << fInputDigitsBranchName << " not found! Task will be deactivated"
158 << endl;
159 SetActive(kFALSE);
160 return kERROR;
161 }
162 fBmnVSPDigitMatchesArray = (TClonesArray*)ioman->GetObject(fInputDigitMatchesBranchName);
163
164 if (fVerbose > 1) {
165 if (fBmnVSPDigitMatchesArray)
166 cout << " Strip matching information exists!\n";
167 else
168 cout << " Strip matching information doesn`t exist!\n";
169 }
170
171 fBmnVSPHitsArray = new TClonesArray(fOutputHitsBranchName);
172 ioman->Register(fOutputHitsBranchName, "VSP", fBmnVSPHitsArray, kTRUE);
173 fBmnVSPUpperClustersArray = new TClonesArray("StripCluster");
174 ioman->Register("BmnVSPUpperCluster", "VSP", fBmnVSPUpperClustersArray, kTRUE);
175 fBmnVSPLowerClustersArray = new TClonesArray("StripCluster");
176 ioman->Register("BmnVSPLowerCluster", "VSP", fBmnVSPLowerClustersArray, kTRUE);
177
178 if (fVerbose > 1)
179 cout << "=================== BmnVSPHitMaker::Init() finished ===============" << endl;
180
181 return kSUCCESS;
182}
183
184void BmnVSPHitMaker::Exec(Option_t* opt)
185{
186
187 TStopwatch sw;
188 sw.Start();
189
190 if (!IsActive())
191 return;
192
193 fBmnVSPHitsArray->Delete();
194 fBmnVSPUpperClustersArray->Delete();
195 fBmnVSPLowerClustersArray->Delete();
196
199
200 if (fVerbose > 1)
201 cout << "=================== BmnVSPHitMaker::Exec() started ================" << endl;
202 if (fVerbose > 1)
203 cout << " BmnVSPHitMaker::Exec(), Number of BmnVSPDigits = " << fBmnVSPDigitsArray->GetEntriesFast() << "\n";
204
206
207 if (fVerbose > 1)
208 cout << "=================== BmnVSPHitMaker::Exec() finished ===============" << endl;
209
210 sw.Stop();
211 workTime += sw.RealTime();
212}
213
215{
216
217 // FairMCPoint* MCPoint;
218 BmnVSPDigit* digit;
219 BmnMatch* strip_match; // MC-information for a strip
220
221 BmnVSPStation* station;
222 BmnVSPModule* module;
223
224 // Loading digits ---------------------------------------------------------------
225 Int_t AddedDigits = 0;
226 Int_t AddedStripDigitMatches = 0;
227
228 for (Int_t idigit = 0; idigit < fBmnVSPDigitsArray->GetEntriesFast(); idigit++) {
229 digit = (BmnVSPDigit*)fBmnVSPDigitsArray->At(idigit);
230 if (!digit->IsGoodDigit())
231 continue;
232
233 if (digit->GetStripSignal() < fSignalLow || digit->GetStripSignal() > fSignalUp)
234 continue;
235
236 station = StationSet->GetVSPStation(digit->GetStation());
237 module = station->GetModule(digit->GetModule());
238 if (!module)
239 return;
240
241 if (module->SetStripSignalInLayerByZoneId(digit->GetStripLayer(), digit->GetStripNumber(),
242 digit->GetStripSignal()))
243 AddedDigits++;
244
245 // Add a MC-match to the current strip if this MC-match array exists
246 if (fBmnVSPDigitMatchesArray) {
247 strip_match = (BmnMatch*)fBmnVSPDigitMatchesArray->At(idigit);
248 if (module->SetStripMatchInLayerByZoneId(digit->GetStripLayer(), digit->GetStripNumber(), *strip_match))
249 AddedStripDigitMatches++;
250 }
251
252 // Add a digit number match to the current strip
253 BmnMatch stripDigitNumberMatch; // digit number information for the current strip
254 stripDigitNumberMatch.AddLink(1.0, idigit);
255 module->SetStripDigitNumberMatchInLayerByZoneId(digit->GetStripLayer(), digit->GetStripNumber(),
256 stripDigitNumberMatch);
257 }
258
259 if (fVerbose > 1)
260 cout << " Processed strip digits : " << AddedDigits << "\n";
261 if (fVerbose > 1 && fBmnVSPDigitMatchesArray)
262 cout << " Added strip digit MC-matches : " << AddedStripDigitMatches << "\n";
263 //------------------------------------------------------------------------------
264
265 // Processing digits
266 StationSet->ProcessPointsInDetector();
267
268 Int_t NCalculatedPoints = StationSet->CountNProcessedPointsInDetector();
269 if (fVerbose > 1)
270 cout << " Calculated points : " << NCalculatedPoints << "\n";
271 if (fVerbose == 1)
272 cout << "BmnVSPHitMaker: " << NCalculatedPoints << " hits\n";
273
274 Int_t clear_matched_points_cnt = 0; // points with the only one match-indexes
275
276 map<Int_t, StripCluster> UniqueUpperClusters;
277 map<Int_t, StripCluster> UniqueLowerClusters;
278
279 for (Int_t iStation = 0; iStation < StationSet->GetNStations(); ++iStation) {
280 station = StationSet->GetVSPStation(iStation);
281
282 for (Int_t iModule = 0; iModule < station->GetNModules(); ++iModule) {
283 module = station->GetModule(iModule);
284
285 for (Int_t iLayer = 0; iLayer < module->GetNStripLayers(); ++iLayer) {
286 BmnVSPLayer layer = module->GetStripLayer(iLayer);
287
288 for (size_t iCluster = 0; iCluster < layer.GetStripClusters().size(); ++iCluster) {
289 StripCluster cluster = layer.GetStripClusters()[iCluster];
290 cluster.SetModule(iModule);
291 cluster.SetStation(iStation);
292
293 // Saving all strip clusters into branches (each cluster have its own unique id in the current
294 // event)
295 if (cluster.GetType() == UpperStripLayer) {
296 new ((*fBmnVSPUpperClustersArray)[fBmnVSPUpperClustersArray->GetEntriesFast()])
297 StripCluster(cluster);
298 } else {
299 new ((*fBmnVSPLowerClustersArray)[fBmnVSPLowerClustersArray->GetEntriesFast()])
300 StripCluster(cluster);
301 }
302 }
303 }
304
305 Int_t NIntersectionPointsInModule = module->GetNIntersectionPoints();
306
307 // Int_t NPseudoIntersectionsInModule = module->GetNPseudoIntersections();
308 Int_t NPseudoIntersectionsInModule = 0; // !!! IGNORE PSEUDO INTERSECTIONS !!!
309
310 for (Int_t iPoint = 0; iPoint < NIntersectionPointsInModule + NPseudoIntersectionsInModule; ++iPoint) {
311
312 // Double_t threshold = 2000000;// ADC
313
314 // Double_t sigL;
315 // Double_t sigU;
316 Double_t x, y, z;
317 Double_t x_err, y_err, z_err;
318
319 if (iPoint < NIntersectionPointsInModule) {
320 // sigL = module->GetIntersectionPoint_LowerLayerSripTotalSignal(iPoint);
321 // sigU = module->GetIntersectionPoint_UpperLayerSripTotalSignal(iPoint);
322 x = module->GetIntersectionPointX(iPoint);
323 y = module->GetIntersectionPointY(iPoint);
324 z = module->GetZPositionRegistered();
325 x_err = module->GetIntersectionPointXError(iPoint);
326 y_err = module->GetIntersectionPointYError(iPoint);
327 z_err = 0.0;
328 } else {
329 // sigL =
330 // module->GetPseudoIntersection_LowerLayerSripTotalSignal(iPoint-NIntersectionPointsInModule); sigU
331 // = module->GetPseudoIntersection_UpperLayerSripTotalSignal(iPoint-NIntersectionPointsInModule);
332 x = module->GetPseudoIntersectionX(iPoint - NIntersectionPointsInModule);
333 y = module->GetPseudoIntersectionY(iPoint - NIntersectionPointsInModule);
334 z = module->GetZPositionRegistered();
335 x_err = module->GetPseudoIntersectionXError(iPoint - NIntersectionPointsInModule);
336 y_err = module->GetPseudoIntersectionYError(iPoint - NIntersectionPointsInModule);
337 z_err = 0.0;
338 }
339
340 z += fIsExp ? fAlignCor[iStation][iModule][2] : 0.; // alignment shift
341
342 // Transform hit coordinates from local coordinate system to global
343 if (TransfSet) {
344 Plane3D::Point glob_point = TransfSet->ApplyTransforms(Plane3D::Point(-x, y, z), iStation, iModule);
345 x = -glob_point.X();
346 y = glob_point.Y();
347 z = glob_point.Z();
348 }
349
350 Int_t RefMCIndex = -1;
351
352 // MC-matching for the current hit (define RefMCIndex)) ---------
353 BmnMatch mc_match_hit;
354 if (iPoint < NIntersectionPointsInModule) {
355 mc_match_hit = module->GetIntersectionPointMatch(iPoint);
356 } else {
357 mc_match_hit = module->GetPseudoIntersectionMatch(iPoint - NIntersectionPointsInModule);
358 }
359
360 Int_t most_probably_index = -1;
361 Double_t max_weight = 0;
362
363 Int_t n_links = mc_match_hit.GetNofLinks();
364 if (n_links == 1)
365 clear_matched_points_cnt++;
366 for (Int_t ilink = 0; ilink < n_links; ilink++) {
367 Int_t index = mc_match_hit.GetLink(ilink).GetIndex();
368 Double_t weight = mc_match_hit.GetLink(ilink).GetWeight();
369 if (weight > max_weight) {
370 max_weight = weight;
371 most_probably_index = index;
372 }
373 }
374
375 RefMCIndex = most_probably_index;
376 //--------------------------------------------------------------
377
378 // Add hit ------------------------------------------------------
379 x *= -1; // invert to global X
380 Double_t deltaX = fIsExp ? fAlignCor[iStation][iModule][0] : 0.;
381 Double_t deltaY = fIsExp ? fAlignCor[iStation][iModule][1] : 0.;
382
383 x += deltaX;
384 y += deltaY;
385
386 new ((*fBmnVSPHitsArray)[fBmnVSPHitsArray->GetEntriesFast()])
387 BmnVSPHit(0, TVector3(x, y, z), TVector3(x_err, y_err, z_err), RefMCIndex);
388
389 BmnVSPHit* hit = (BmnVSPHit*)fBmnVSPHitsArray->At(fBmnVSPHitsArray->GetEntriesFast() - 1);
390 hit->SetStation(iStation);
391 hit->SetModule(iModule);
392 hit->SetIndex(fBmnVSPHitsArray->GetEntriesFast() - 1);
393
394 if (iPoint < NIntersectionPointsInModule) {
396 module->GetIntersectionPointDigitNumberMatch(iPoint)); // digit number match for the hit
397 hit->SetUpperClusterIndex(module->GetUpperCluster(iPoint).GetUniqueID());
398 hit->SetLowerClusterIndex(module->GetLowerCluster(iPoint).GetUniqueID());
399 hit->SetPseudo(false);
400 } else {
401 hit->SetDigitNumberMatch(module->GetPseudoIntersectionDigitNumberMatch(
402 iPoint - NIntersectionPointsInModule)); // digit number match for the hit
404 module->GetUpperCluster_PseudoIntersections(iPoint - NIntersectionPointsInModule)
405 .GetUniqueID());
407 module->GetLowerCluster_PseudoIntersections(iPoint - NIntersectionPointsInModule)
408 .GetUniqueID());
409 hit->SetPseudo(true);
410 }
411 //--------------------------------------------------------------
412 /*
413 StripCluster ucls = module->GetUpperCluster(iPoint);
414 StripCluster lcls = module->GetLowerCluster(iPoint);
415 ucls.SetModule(iModule);
416 lcls.SetModule(iModule);
417 ucls.SetStation(iStation);
418 lcls.SetStation(iStation);
419
420 Int_t upperClusterIndex = ucls.GetUniqueID();
421 Int_t lowerClusterIndex = lcls.GetUniqueID();
422
423 //Filtration of duplicate clusters in the current module for the same type of strips
424 Bool_t duplicateUpperCluster = false;
425 for(auto it : UniqueUpperClusters) {
426 StripCluster cluster = it.second;
427
428 if(cluster.GetClusterSize() != ucls.GetClusterSize()) continue;
429
430 Int_t identic_cnt = 0;
431 for(Int_t istrip = 0; istrip < cluster.GetClusterSize(); ++istrip) {
432 Int_t strip_num = cluster.Strips.at(istrip);
433 Double_t strip_signal = cluster.Signals.at(istrip);
434
435 if( (strip_num == ucls.Strips.at(istrip)) && (strip_signal ==
436 ucls.Signals.at(istrip)) ) identic_cnt++; else break;
437 }
438 if(identic_cnt == ucls.GetClusterSize()) {
439 upperClusterIndex = cluster.GetUniqueID();
440 duplicateUpperCluster = true;
441 break;
442 }
443 }
444
445 Bool_t duplicateLowerCluster = false;
446 for(auto it : UniqueLowerClusters) {
447 StripCluster cluster = it.second;
448
449 if(cluster.GetClusterSize() != lcls.GetClusterSize()) continue;
450
451 Int_t identic_cnt = 0;
452 for(Int_t istrip = 0; istrip < cluster.GetClusterSize(); ++istrip) {
453 Int_t strip_num = cluster.Strips.at(istrip);
454 Double_t strip_signal = cluster.Signals.at(istrip);
455
456 if( (strip_num == lcls.Strips.at(istrip)) && (strip_signal ==
457 lcls.Signals.at(istrip)) ) identic_cnt++; else break;
458 }
459 if(identic_cnt == lcls.GetClusterSize()) {
460 lowerClusterIndex = cluster.GetUniqueID();
461 duplicateLowerCluster = true;
462 break;
463 }
464 }
465 //--------------------------------------------------------------
466
467 if(!duplicateUpperCluster) UniqueUpperClusters[upperClusterIndex] = ucls;
468 if(!duplicateLowerCluster) UniqueLowerClusters[lowerClusterIndex] = lcls;
469 hit->SetUpperClusterIndex(upperClusterIndex);
470 hit->SetLowerClusterIndex(lowerClusterIndex);
471 */
472 if (fHitMatching) {
473 BmnMatch digiMatch;
474 Int_t idx0, idx1;
475 if (iPoint < NIntersectionPointsInModule) {
476 digiMatch = module->GetIntersectionPointDigitNumberMatch(iPoint);
477 idx0 = digiMatch.GetLink(0).GetIndex();
478 idx1 = digiMatch.GetLink(1).GetIndex();
479 } else {
480 digiMatch = module->GetPseudoIntersectionDigitNumberMatch(iPoint - NIntersectionPointsInModule);
481 idx0 = digiMatch.GetLink(0).GetIndex();
482 idx1 = digiMatch.GetLink(0).GetIndex();
483 }
484 BmnMatch* digiMatch0 = (BmnMatch*)fBmnVSPDigitMatchesArray->At(idx0);
485 BmnMatch* digiMatch1 = (BmnMatch*)fBmnVSPDigitMatchesArray->At(idx1);
486
487 Bool_t hitOk = kFALSE;
488 for (Int_t ilink = 0; ilink < digiMatch0->GetNofLinks(); ilink++) {
489 Int_t iindex = digiMatch0->GetLink(ilink).GetIndex();
490 for (Int_t jlink = 0; jlink < digiMatch1->GetNofLinks(); jlink++) {
491 Int_t jindex = digiMatch1->GetLink(jlink).GetIndex();
492 if (iindex == jindex) {
493 hitOk = kTRUE;
494 break;
495 }
496 }
497 if (hitOk)
498 break;
499 }
500
501 hit->SetType(hitOk);
502 if (!hitOk)
503 hit->SetRefIndex(-1);
504
505 // hit MC-matching ----------------------------------------------
506 FairRootManager::Instance()->SetUseFairLinks(kTRUE);
507 BmnMatch hitMatch;
508 if (iPoint < NIntersectionPointsInModule) {
509 hitMatch = module->GetIntersectionPointMatch(iPoint);
510 } else {
511 hitMatch = module->GetPseudoIntersectionMatch(iPoint - NIntersectionPointsInModule);
512 }
513
514 for (BmnLink lnk : hitMatch.GetLinks())
515 hit->AddLink(FairLink(-1, lnk.GetIndex(), lnk.GetWeight()));
516 FairRootManager::Instance()->SetUseFairLinks(kFALSE);
517 }
518 //--------------------------------------------------------------
519 }
520 }
521 }
522
523 /*
524 for (auto it : UniqueUpperClusters) {
525 for (Int_t i = 0; i < fBmnVSPHitsArray->GetEntriesFast(); i++) {
526 BmnVSPHit* hit = (BmnVSPHit*)fBmnVSPHitsArray->At(i);
527 if (hit->GetUpperClusterIndex() != it.first) continue;
528 hit->SetUpperClusterIndex(fBmnVSPUpperClustersArray->GetEntriesFast());
529 }
530 it.second.SetUniqueID(fBmnVSPUpperClustersArray->GetEntriesFast());
531 new ((*fBmnVSPUpperClustersArray)[fBmnVSPUpperClustersArray->GetEntriesFast()]) StripCluster(it.second);
532 }
533 for (auto it : UniqueLowerClusters) {
534 for (Int_t i = 0; i < fBmnVSPHitsArray->GetEntriesFast(); i++) {
535 BmnVSPHit* hit = (BmnVSPHit*)fBmnVSPHitsArray->At(i);
536 if (hit->GetLowerClusterIndex() != it.first) continue;
537 hit->SetLowerClusterIndex(fBmnVSPLowerClustersArray->GetEntriesFast());
538 }
539 it.second.SetUniqueID(fBmnVSPLowerClustersArray->GetEntriesFast());
540 new ((*fBmnVSPLowerClustersArray)[fBmnVSPLowerClustersArray->GetEntriesFast()]) StripCluster(it.second);
541 }
542 */
543
544 if (fVerbose > 1)
545 cout << " N clear matches with MC-points = " << clear_matched_points_cnt << "\n";
546 //------------------------------------------------------------------------------
547 StationSet->Reset();
548}
549
550void BmnVSPHitMaker::OnlineWrite(const std::unique_ptr<TTree>& resultTree)
551{
552 if (!IsActive())
553 return;
554
555 resultTree->Branch(fOutputHitsBranchName, &fBmnVSPHitsArray);
556 resultTree->Branch("BmnVSPUpperCluster", &fBmnVSPUpperClustersArray);
557 resultTree->Branch("BmnVSPLowerCluster", &fBmnVSPLowerClustersArray);
558 resultTree->Fill();
559}
560
562{
563 printf("Work time of BmnVSPHitMaker: %4.2f sec.\n", workTime);
564}
int i
Definition P4_F32vec4.h:22
@ UpperStripLayer
void SetModule(Int_t mod)
Definition BmnHit.h:73
void SetType(Int_t type)
Definition BmnHit.h:81
void SetUpperClusterIndex(Int_t idx)
Definition BmnHit.h:133
void SetDigitNumberMatch(BmnMatch match)
Definition BmnHit.h:121
void SetIndex(Int_t id)
Definition BmnHit.h:57
void SetPseudo(Bool_t pseudo)
Definition BmnHit.h:141
void SetStation(Short_t st)
Definition BmnHit.h:69
void SetLowerClusterIndex(Int_t idx)
Definition BmnHit.h:137
const BmnLink & GetLink(Int_t i) const
Definition BmnMatch.h:37
void AddLink(const BmnMatch &match)
Definition BmnMatch.cxx:43
const vector< BmnLink > & GetLinks() const
Definition BmnMatch.h:38
Int_t GetNofLinks() const
Definition BmnMatch.h:40
Bool_t IsGoodDigit()
Int_t GetStripNumber()
Int_t GetStripLayer()
Int_t GetStation()
Double_t GetStripSignal()
virtual void OnlineWrite(const std::unique_ptr< TTree > &resultTree)
Write task resul to tree.
virtual void Finish()
virtual InitStatus Init()
virtual void Exec(Option_t *opt)
virtual ~BmnVSPHitMaker()
static void SetLowerUniqueID(Int_t id)
static void SetUpperUniqueID(Int_t id)
vector< StripCluster > GetStripClusters()
Int_t CountNProcessedPointsInDetector()
BmnVSPStation * GetVSPStation(Int_t station_num)
Int_t GetNModules()
Bool_t LoadFromXMLFile(TString xml_config_file)
Plane3D::Point ApplyTransforms(Plane3D::Point point, Int_t station, Int_t module)
void SetStation(Int_t st)
Int_t GetType()
void SetModule(Int_t m)
int GetValue(vector< UniValue * > &parameter_value)
get value of detector parameter presented by an array
static UniDetectorParameter * GetDetectorParameter(int value_id)
get detector parameter from the database