FairRoot/PandaRoot
PndGemFindHits.cxx
Go to the documentation of this file.
1 //* $Id: */
2 
3 // -------------------------------------------------------------------------
4 // ----- PndGemFindHits source file -----
5 // ----- Created 15/02/2009 by R. Karabowicz -----
6 // -------------------------------------------------------------------------
7 
17 #include "PndGemFindHits.h"
18 
19 #include "PndGemDigi.h"
20 #include "PndGemCluster.h"
21 #include "PndGemDigiPar.h"
22 #include "PndGemHit.h"
23 #include "PndGemSensor.h"
24 #include "PndGemStation.h"
25 #include "PndGemMonitor.h"
26 
27 #include "PndDetectorList.h"
28 
29 #include "FairRootManager.h"
30 #include "FairRunAna.h"
31 #include "FairRuntimeDb.h"
32 
33 #include "TClonesArray.h"
34 #include "TMath.h"
35 #include "TH2F.h"
36 
37 #include <iomanip>
38 
39 using std::cerr;
40 using std::cout;
41 using std::endl;
42 using std::fixed;
43 using std::flush;
44 using std::left;
45 using std::map;
46 using std::right;
47 using std::set;
48 using std::setprecision;
49 using std::setw;
50 
51 // ----- Default constructor ------------------------------------------
53  PndPersistencyTask("GEM Hit Finder", 1),
54  fMonitor(NULL),
55  fDigiPar(NULL),
56  fDigis (NULL),
57  fHits (NULL),
58  fHitsTemp (NULL),
59  fMCPointBranchId(-1),
60  fUseClusters(kFALSE),
61  fTimeOrderedDigi(kFALSE),
62  fPrepTime(0.),
63  fSortTime(0.),
64  fCreateTime(0.),
65  fConfirmTime(0.),
66  fActivateTime(0.),
67  fAllTime(0.),
68  fHitWindow(1.5),
69  fTNofEvents(0),
70  fTNofDigis (0),
71  fTNofHits (0),
72  fTNofHitsTemp (0)
73 {
74  SetPersistency(kTRUE);
75 }
76 // -------------------------------------------------------------------------
77 
78 
79 
80 // ----- Standard constructor ------------------------------------------
82  : PndPersistencyTask("GEMFindHits", iVerbose),
83  fMonitor(NULL),
84  fDigiPar(NULL),
85  fDigis (NULL),
86  fHits (NULL),
87  fHitsTemp (NULL),
88  fMCPointBranchId(-1),
89  fUseClusters(kFALSE),
90  fTimeOrderedDigi(kFALSE),
91  fPrepTime(0.),
92  fSortTime(0.),
93  fCreateTime(0.),
94  fConfirmTime(0.),
95  fActivateTime(0.),
96  fAllTime(0.),
97  fHitWindow(1.5),
98  fTNofEvents(0),
99  fTNofDigis (0),
100  fTNofHits (0),
101  fTNofHitsTemp (0)
102 {
103  SetPersistency(kTRUE);
104 }
105 // -------------------------------------------------------------------------
106 
107 
108 
109 // ----- Constructor with name -----------------------------------------
111  : PndPersistencyTask(name, iVerbose),
112  fMonitor(NULL),
113  fDigiPar(NULL),
114  fDigis (NULL),
115  fHits (NULL),
116  fHitsTemp (NULL),
117  fMCPointBranchId(-1),
118  fUseClusters(kFALSE),
119  fTimeOrderedDigi(kFALSE),
120  fPrepTime(0.),
121  fSortTime(0.),
122  fCreateTime(0.),
123  fConfirmTime(0.),
124  fActivateTime(0.),
125  fAllTime(0.),
126  fHitWindow(1.5),
127  fTNofEvents(0),
128  fTNofDigis (0),
129  fTNofHits (0),
130  fTNofHitsTemp (0)
131 {
132  SetPersistency(kTRUE);
133 }
134 // -------------------------------------------------------------------------
135 
136 
137 
138 // ----- Destructor ----------------------------------------------------
140  if ( fHits ) {
141  fHits->Delete();
142  delete fHits;
143  }
144  if ( fHitsTemp ) {
145  fHitsTemp->Delete();
146  delete fHitsTemp;
147  }
148 }
149 // -------------------------------------------------------------------------
150 
151 
152 
153 // ----- Public method Exec --------------------------------------------
154 void PndGemFindHits::Exec(Option_t*) {
155  if ( fVerbose > 0 ) cout << "-I- PndGemFindHits::Exec()" << endl;
156  // cout << "======== PndGemFindHits::Exec(Event = " << fTNofEvents << " ) ====================" << endl;
157 
158  fTimer.Start();
159 
160  fTNofEvents++;
161 
162  Bool_t warn = kFALSE;
163 
164  // Clear output array
165  fHits->Delete();
166  fHitsTemp->Delete();
167 
168  fPrepTime+=fTimer.RealTime(); // preparation time
169  fTimer.Continue();
170 
171  // Sort GEM digis with respect to sectors
172  SortDigis();
173 
174  fSortTime+=fTimer.RealTime();
175  fTimer.Continue();
176 
177  // Find hits in sectors
178  Int_t nDigisF = 0;
179  Int_t nDigisB = 0;
180  Int_t nHits = 0;
181  Int_t nHitsTemp = 0;
182 
183  // cout << "GEM found hits: " << endl;
184 
185  Int_t nStations = fDigiPar->GetNStations();
186  for (Int_t iStation=0; iStation<nStations; iStation++) {
187  PndGemStation* station = (PndGemStation*)fDigiPar->GetStation(iStation);
188  Int_t nDigisFInStation = 0;
189  Int_t nDigisBInStation = 0;
190  Int_t nHitsInStation = 0;
191  Int_t nSensors = station->GetNSensors();
192  for (Int_t iSensor=0; iSensor<nSensors; iSensor++) {
193  PndGemSensor* sensor = (PndGemSensor*)station->GetSensor(iSensor);
194 
195  set <Int_t> fSet, bSet;
196  if ( fDigiMapF.find(sensor) == fDigiMapF.end() ) {
197  cout << "-E- " << fName << "::Exec: sensor "
198  << sensor->GetSensorNr() << " of station "
199  << station->GetStationNr() << "not found in front map!"
200  << endl;
201  warn = kTRUE;
202  continue;
203  }
204  fSet = fDigiMapF[sensor];
205  if ( fDigiMapB.find(sensor) == fDigiMapB.end() ) {
206  cout << "-E- " << fName << "::Exec: sensor "
207  << sensor->GetSensorNr() << " of station "
208  << station->GetStationNr() << "not found in back map!"
209  << endl;
210  warn = kTRUE;
211  continue;
212  }
213  bSet = fDigiMapB[sensor];
214 
215  Int_t nDigisFInSensor = fSet.size();
216  Int_t nDigisBInSensor = bSet.size();
217 
218 
219  Int_t nHitsInSensor;
220  //Int_t nHitsInSensor = FindHits(sensor, fSet, bSet);
221 
222  // TimeBased use basically same method..
223  if ( fTimeOrderedDigi ) {
224  nHitsInSensor = FindHits(sensor, fSet, bSet);//return non-0 now.. why it was 0??
225  }
226  else {
227  nHitsInSensor = FindHits2(sensor, fSet, bSet);//test version ok
228  }
229 
230  if ( fVerbose > 2 )
231  cout << "Sensor " << sensor->GetSensorNr()
232  << ", Digis front " << nDigisFInSensor
233  << ", Digis Back " << nDigisBInSensor
234  << ", Hit combinations " << nHitsInSensor << endl;
235  nHitsInStation += nHitsInSensor;
236  nDigisFInStation += nDigisFInSensor;
237  nDigisBInStation += nDigisBInSensor;
238  } // Sector loop
239 
240  if ( fVerbose > 1 ) cout << "Total for station "
241  << station->GetStationNr() << ": Digis front "
242  << nDigisFInStation << ", digis back "
243  << nDigisBInStation << ", hit combinations "
244  << nHitsInStation << endl;
245  nDigisB += nDigisBInStation;
246  nDigisF += nDigisFInStation;
247  nHitsTemp += nHitsInStation;
248 
249  } // Station loop
250 
251  if ( fVerbose > 1 ) cout << "Total for GEM "
252  << ": Digis front "<< nDigisF
253  << ", digis back "<< nDigisB
254  << ", hit combinations "<< nHitsTemp << endl;
255 
256  fCreateTime+=fTimer.RealTime();
257  fTimer.Continue();
258 
259  // can only confirm hits using clusters fromt the same event
260  if ( fTimeOrderedDigi ) {
261  if ( fUseClusters ) {
262  ConfirmHits();
263  fConfirmTime+=fTimer.RealTime();
264  fTimer.Continue();
265  // there should be another funciton, that will confirm the hits using fired strips in previous events, that could have obscured the cluster/hit finding in this event because of the small time difference
266  // activate digis from the clusters in the PndGemMonitor
267  ActivateDigis();
268  fActivateTime+=fTimer.RealTime();
269  fTimer.Continue();
270  // check the digis in the previous events...
271  }
272  }
273  else {
274  ConfirmHits2();//test version ok
275  }
276 
277  if ( fVerbose > 1 ) {
278  // cout << "-I- PndGemFindHits::Exec() Created " << fHits->GetEntriesFast() << " hits from "
279  //<< nDigisF << " front and " << nDigisB << " back digis" << endl;
280  cout << endl;
281  cout << "-I- " << fName << ":Event summary" << endl;
282  cout << " Active channels front side: " << nDigisF << endl;
283  cout << " Active channels back side : " << nDigisB << endl;
284  cout << " Hits created : " << fHits->GetEntriesFast() << endl;
285  cout << " Real time : " << fTimer.RealTime()
286  << endl;
287  }
288  if ( fVerbose == 1 ) {
289  if ( warn ) cout << "- ";
290  else cout << "+ ";
291  cout << setw(15) << left << fName << ": " << setprecision(4) << setw(8)
292  << fixed << right << fTimer.RealTime()
293  << " s, digis " << nDigisF << " / " << nDigisB << ", hits: "
294  << nHits << endl;
295  }
296 
297  // cout << " PndGemFindHits::Exec. Created " << fHits->GetEntries() << " hits out of " << fDigis->GetEntries() << " digis/clusters." << endl;
298 
299  fAllTime+=fTimer.RealTime();
300 
301  fTimer.Stop();
302 }
303 // -------------------------------------------------------------------------
304 
305 
306 
307 
308 // ----- Private method SetParContainers -------------------------------
310 
311  // Get run and runtime database
312  FairRun* run = FairRun::Instance();
313  if ( ! run ) Fatal("SetParContainers", "No analysis run");
314 
315  FairRuntimeDb* db = run->GetRuntimeDb();
316  if ( ! db ) Fatal("SetParContainers", "No runtime database");
317 
318  // Get GEM digitisation parameter container
319  fDigiPar = (PndGemDigiPar*) db->getContainer("PndGemDetectors");
320 }
321 // -------------------------------------------------------------------------
322 
323 
324 
325 // ----- Private method Init -------------------------------------------
326 InitStatus PndGemFindHits::Init() {
327 
328  // Get input array
329  FairRootManager* ioman = FairRootManager::Instance();
330  if ( ! ioman ) Fatal("Init", "No FairRootManager");
331 
332  if ( fUseClusters )
333  fDigis = (TClonesArray*) ioman->GetObject("GEMCluster");
334  else
335  fDigis = (TClonesArray*) ioman->GetObject("GEMDigi");
336 
337  fMCPointBranchId = ioman->GetBranchId("GEMPoint");
338 
339  // Register output array
340  fHits = new TClonesArray("PndGemHit", 1000);
341  ioman->Register("GEMHit", "Hit in GEM", fHits, GetPersistency());
342 
343  // Test Register output array
344  fHitsTemp = new TClonesArray("PndGemHit", 1000);
345  ioman->Register("GEMHitTemp", "temp Hit in GEM", fHitsTemp, kFALSE);
346 
347  // Create sectorwise digi sets
348  MakeSets();
349 
350  cout << "-I- " << fName.Data() << "::Init(). There are " << fDigiPar->GetNStations() << " GEM stations." << endl;
351  cout << "-I- " << fName.Data() << "::Init(). Initialization succesfull." << endl;
352 
354  PndGemSensor* sensor = (PndGemSensor*)station->GetSensor(0);
355  cout << "sensor out rad is " << sensor->GetOuterRadius() << endl;
356 
357  // fChargeCorrHist = new TH2F("fChargeCorrHist","fChargeCorrHist",1000,0.,100.,1000,0.,100.);
358  // fChargeDiffHist = new TH2F("fChargeDiffHist","fChargeDiffHist",1000,0.,100.,1000,-5.,5.);
359 
360  // Open GEM Monitor
362 
363  for ( Int_t istat = 0 ; istat < fDigiPar->GetNStations() ; istat++ ) {
364  station = (PndGemStation*)fDigiPar->GetStation(istat);
365  for ( Int_t isens = 0 ; isens < station->GetNSensors() ; isens++ ) {
366  sensor = (PndGemSensor*)station->GetSensor(isens);
367  fMonitor->CreateSensorMonitor(*sensor);
368  }
369  }
370 
371  return kSUCCESS;
372 }
373 // -------------------------------------------------------------------------
374 
375 
376 
377 
378 // ----- Private method ReInit -----------------------------------------
380 
381  // Create sectorwise digi sets
382  MakeSets();
383 
384  return kSUCCESS;
385 }
386 // -------------------------------------------------------------------------
387 
388 
389 
390 
391 // ----- Private method MakeSets ---------------------------------------
393 
394  fDigiMapF.clear();
395  fDigiMapB.clear();
396  Int_t nStations = fDigiPar->GetNStations();
397  for (Int_t iStation=0; iStation<nStations; iStation++) {
398  PndGemStation* station = (PndGemStation*)fDigiPar->GetStation(iStation);
399  Int_t nSensors = station->GetNSensors();
400  for (Int_t iSensor=0; iSensor<nSensors; iSensor++) {
401  PndGemSensor* sensor = (PndGemSensor*)station->GetSensor(iSensor);
402  set<Int_t> a;
403  fDigiMapF[sensor] = a;
404  set<Int_t> b;
405  fDigiMapB[sensor] = b;
406  }
407  }
408 }
409 // -------------------------------------------------------------------------
410 
411 // ----- Private method ActivateDigis --------------------------------------
413  PndGemCluster* cluster;
414  for ( Int_t iclus = 0 ; iclus < fDigis->GetEntries() ; iclus++ ) {
415  cluster = (PndGemCluster*)fDigis->At(iclus);
416 
417  fMonitor->EnableCluster(fTNofEvents,iclus,cluster);
418  }
419 }
420 // -------------------------------------------------------------------------
421 
422 // ----- Private method ConfirmHits --------------------
424  if ( fVerbose > 0 ) cout << "-I- PndGemFindHits::ConfirmHits()" << endl;
425  Bool_t printInfo = kFALSE;
426  if( fVerbose > 1 ) printInfo = kTRUE;
427 
428  PndGemSensor* sensor = NULL;
429  PndGemHit* hit = NULL;
430  PndGemCluster* cluster = NULL;
431 
432  for ( Int_t ihit = 0 ; ihit < fHits->GetEntriesFast() ; ihit++ ) {
433  hit = (PndGemHit*)fHits->At(ihit);
434 
435  sensor = (PndGemSensor*)fDigiPar->GetSensor(hit->GetStationNr(),
436  (3-hit->GetSensorNr ()));
437  Double_t hitProjX = hit->GetX()*(sensor->GetZ0()/hit->GetZ());
438  Double_t hitProjY = hit->GetY()*(sensor->GetZ0()/hit->GetZ());
439  // if ( abs(hit->GetX()+3.8)<0.1 &&
440  // abs(hit->GetY()-2.0)<0.1 &&
441  // hit->GetStationNr()==2 )
442  // printInfo = kTRUE;
443  // else
444  // printInfo = kFALSE;
445 
446  if ( printInfo ) {
447  cout << "COMPARING HIT " << ihit << " at "
448  << hit->GetX() << ","
449  << hit->GetY() << ","
450  << hit->GetZ() << " @ "
451  << hit->GetTimeStamp() << endl;
452  }
453 
454  Bool_t hitConfirmed[2] = {kFALSE,kFALSE};
455  Int_t channelNr [2] = {sensor->GetChannel(hitProjX,hitProjY,0),
456  sensor->GetChannel(hitProjX,hitProjY,1)};
457 
458  // check in previous events....
459  for ( Int_t iside = 0 ; iside < 2 ; iside++ ) {
460  Double_t lastChannelActivateTime = fMonitor->ChannelLastActiveAt(hit->GetStationNr(),
461  (3-hit->GetSensorNr()),
462  iside,
463  channelNr[iside]);
464  if ( printInfo ) {
465  cout << " on " << hit->GetStationNr() << "." << 3-hit->GetSensorNr() << "." << iside << " channel " << channelNr[iside] << " was last active at " << lastChannelActivateTime << endl;
466  // fMonitor->Print();
467  }
468  if ( lastChannelActivateTime < 0 ) continue;
469  if ( hit->GetTimeStamp()-lastChannelActivateTime < 100 &&
470  hit->GetTimeStamp()-lastChannelActivateTime > -10 ) {
471  if (printInfo) cout << "HIT CONFIRMED ON SIDE " << iside << " IN SOME PREVIOUS EVENT (" << hit->GetTimeStamp()-lastChannelActivateTime << " ns ago)" << endl;
472  hitConfirmed[iside] = kTRUE;
473  }
474  }
475 
476  if ( hitConfirmed[0] == kFALSE || hitConfirmed[1] == kFALSE ) {
477  for ( Int_t iclus = 0 ; iclus < fDigis->GetEntries() ; iclus++ ) {
478  cluster = (PndGemCluster*)fDigis->At(iclus);
479 
480  if ( cluster->GetTimeStamp() < hit->GetTimeStamp() - 111. ||
481  cluster->GetTimeStamp() > hit->GetTimeStamp() + 11. )
482  continue;
483 
484  if ( cluster->GetStationNr() != hit->GetStationNr() ) continue;
485  if ( cluster->GetSensorNr () == hit->GetSensorNr () ) continue;
486 
487  // should have some function in sensor for that:
488  if ( channelNr[cluster->GetSide()] <= cluster->GetClusterEnd() &&
489  channelNr[cluster->GetSide()] >= cluster->GetClusterBeg() ) {
490  hitConfirmed[cluster->GetSide()] = kTRUE;
491  if ( printInfo ) {
492  cout << "HIT CONFIRMED ON SIDE " << cluster->GetSide() << " BY CLUSTER " << iclus << " at "
493  << (cluster->GetSide()==0?"F":"B") << "side on "
494  << cluster->GetStationNr() << "."
495  << cluster->GetSensorNr() << " from "
496  << (cluster->GetClusterBeg()) << " to "
497  << (cluster->GetClusterEnd()) << " **** "
498  << sensor->GetChannel(hitProjX,hitProjY,cluster->GetSide())
499  << ", time_diff = " << hit->GetTimeStamp() << " - " << cluster->GetTimeStamp() << "ns"
500  << endl;
501  }
502  }
503  }
504  }
505 
506  if ( hitConfirmed[0] == kTRUE && hitConfirmed[1] == kTRUE ) {
507  if ( printInfo ) cout << "GOOD HIT" << endl;
508  ((PndGemHit*)fHits->At(ihit))->SetCharge(10.);
509  }
510  else {
511  if ( printInfo ) cout << " NO CLUSTER MATCHES TO THIS HIT" << endl;
512 
513 
514 
515  }
516  }
517 }
518 // -------------------------------------------------------------------------
519 
520 // ----- Private method ConfirmHits --------------------
521 //Check corresponding sensor and store hit when both sensors have same location of hits
522 //Calculate acceptable width of hit matching based on Dx Dy
523 //So far only for EB.. Channel activation is not considered...
525  if ( fVerbose > 0 ) cout << "-I- PndGemFindHits::ConfirmHits2()" << endl;
526 
527  PndGemHit* hitTemp = NULL;//hit in all combination of clusters to read
528  PndGemHit* hitTemp2 = NULL;//hit in all combination of clusters to compare
529  PndGemSensor* sensor = NULL;
530 
531  Double_t xHit[2],yHit[2],zHit[2];
532  Double_t dx[2], dy[2], dz[2];
533  TVector3 pos, dpos;
534 
535  Int_t nHits = fHits->GetEntriesFast();
536  Int_t sensorDetId, hitDetId;
537  Bool_t hitMatch = kFALSE;
538  Bool_t hitMatch2 = kFALSE;
539 
540  // window to find corresponding hit on other sensor
541 // Double_t xHitW=1.5; //3rd attempt, 1 seems to be too narrow (7% missed), 2 seems too wide (10% more Hits)
542 // Double_t yHitW=1.5;
543  Double_t xHitW = fHitWindow;
544  Double_t yHitW = fHitWindow;
545 
546  if (! fUseClusters ) {
547  if ( fVerbose > 1 ) cout << " Hit finding window is constant: xHitW=" << xHitW << " yHitW=" << yHitW << endl;
548  }
549 
550  for ( Int_t ihitTemp = 0 ; ihitTemp < fHitsTemp->GetEntriesFast() ; ihitTemp++ ) {
551  hitTemp = (PndGemHit*)fHitsTemp->At(ihitTemp);
552 
553  Int_t hitTempDetId = hitTemp->GetDetectorID();
554  Int_t testGemHit = hitTempDetId & kGemHit << 21;
555  if ( fVerbose > 1 ) {
556  cout << "Check kGemHit bit of hitTempDetId=" << hitTempDetId
557  << " kGemHit=" << kGemHit << " masked value testGemHit=" << testGemHit <<endl;
558  }
559  if( testGemHit != 0 ) {
560  if ( fVerbose > 1 ) {
561  cout << "ihitTemp " << ihitTemp << " : already used... skip" <<endl;
562  }
563  continue;
564  }
565 
566  Int_t staNr = hitTemp->GetStationNr();
567  Int_t senNr = hitTemp->GetSensorNr();
568 
569  xHit[0] = hitTemp->GetX();
570  yHit[0] = hitTemp->GetY();
571  zHit[0] = hitTemp->GetZ();
572  dx[0] = hitTemp->GetDx();
573  dy[0] = hitTemp->GetDy();
574  dz[0] = hitTemp->GetDz();
575 
576  sensor = (PndGemSensor*)fDigiPar->GetSensor(hitTemp->GetStationNr(),
577  hitTemp->GetSensorNr ());
578  sensorDetId = sensor->GetDetectorId();
579  Double_t senD = sensor->GetD();// sensor thickness 1cm (?)
580 
581  if ( fVerbose > 1 ) {
582  cout << "COMPARING HitTemp ihitTemp " << ihitTemp
583  << " of staNr " << staNr
584  << " senNr " << senNr
585  << " at (" << xHit[0] << "," << yHit[0] << "," << zHit[0]
586  << ") time " << hitTemp->GetTimeStamp()
587  <<" dx " << dx[0] << " dy " << dy[0]
588  <<" sensorDetId " << sensorDetId << " senThick " << senD
589  << endl;
590  }
591 
592  hitMatch = kFALSE;
593  for ( Int_t ihitTemp2 = ihitTemp+1 ; ihitTemp2 < fHitsTemp->GetEntriesFast() ; ihitTemp2++ ) {
594  hitTemp2 = (PndGemHit*)fHitsTemp->At(ihitTemp2);
595  if ( fVerbose > 1 ) cout << "WITH ihitTemp2 " << ihitTemp2;
596 
597  hitMatch2 = kFALSE;
598  Int_t hitTemp2DetId = hitTemp2->GetDetectorID();
599  Int_t test2GemHit = hitTemp2DetId & kGemHit << 21;
600  if ( fVerbose > 2 ) {
601  cout << " hitTemp2DetId=" << hitTempDetId
602  << " kGemHit=" << kGemHit << " masked value test2GemHit=" << test2GemHit <<endl;
603  }
604 // if( test2GemHit != 0 ) {
605 // if ( fVerbose > 1 ) cout << " ihitTemp2 " << ihitTemp2 << " : already used... skip" <<endl;
606 // continue;
607 // }
608 
609  Int_t staNr2 = hitTemp2->GetStationNr();
610  if ( fVerbose > 1 ) cout << " ihitTemp2 " << ihitTemp2 << " staNr " << staNr2;
611  if( staNr != staNr2 ) {
612  if ( fVerbose > 1 ) cout << " : diff Station .. skip " << endl;
613  continue;
614  }
615  Int_t senNr2 = hitTemp2->GetSensorNr();
616  if ( fVerbose > 1 ) cout << " of senNr " << senNr2 ;
617  if( senNr == senNr2 ) {
618  if ( fVerbose > 1 ) cout << " : same sensor.. skip " << endl;
619  continue;
620  }
621  xHit[1] = hitTemp2->GetX();
622  yHit[1] = hitTemp2->GetY();
623  zHit[1] = hitTemp2->GetZ();
624  dx[1] = hitTemp2->GetDx();
625  dy[1] = hitTemp2->GetDy();
626  dz[1] = hitTemp2->GetDz();
627  if ( fVerbose > 1 ) {
628  cout << " at (" << xHit[1] << "," << yHit[1] << "," << zHit[1]
629  << ") time " << hitTemp2->GetTimeStamp() <<" dx " << dx[1] << " dy " << dy[1] << endl;
630  }
631 
632  //calc window to search hit on other sensor based on Dx Dy of Hit
633  //assuming straight path (NOT assuming path from target)
634  if ( fUseClusters ) {
635  xHitW = TMath::Sqrt(12.) * dx[0] * TMath::Abs( zHit[0] - zHit[1] ) / senD + TMath::Sqrt(12.) * dx[0] * 0.5;
636  yHitW = TMath::Sqrt(12.) * dy[0] * TMath::Abs( zHit[0] - zHit[1] ) / senD + TMath::Sqrt(12.) * dx[0] * 0.5;
637  if ( fVerbose > 1 ) {
638  cout << " Hit finding window based on dx dy: xHitW=" << xHitW << " yHitW=" << yHitW << endl;
639  }
640  }
641 
642  if( xHit[0] - xHitW < xHit[1] && xHit[1] < xHit[0] + xHitW ){
643  if( yHit[0] - yHitW < yHit[1] && yHit[1] < yHit[0] + yHitW ){
644 
645  hitTempDetId = hitTemp->GetDetectorID();
646  testGemHit = hitTempDetId & kGemHit << 21;
647  if( testGemHit != 0 ) {
648  if ( fVerbose > 1 ) {
649  cout << "ihitTemp " << ihitTemp << " : already stored... " <<endl;
650  }
651  }
652  else {
653  pos.SetXYZ(xHit[0], yHit[0], zHit[0]);
654  dpos.SetXYZ(dx[0], dy[0], dz[0]);
655  hitDetId = hitTemp->GetDetectorID() | kGemHit << 21;
656 
657  TString fromStr = "GEMDigi";
658  if ( fUseClusters ) fromStr = "GEMCluster";
659 
660  new ((*fHits)[nHits++]) PndGemHit(hitDetId,
661  pos, dpos,
662  hitTemp->GetCharge(),
663  hitTemp->GetTimeStamp(),
664  hitTemp->GetDigiNr(0),hitTemp->GetDigiNr(1),
665  hitTemp->GetDr(), hitTemp->GetDp(),
666  hitTemp->GetRefIndex(),fromStr);
667  if ( fVerbose > 1 ) {
668  cout << "ihitTemp "<< ihitTemp<<" xHit yHit zHit " << pos.X() << " / " << pos.Y() << " / " << pos.Z() << " stored as " << nHits << " from " << fromStr.Data() << "." << hitTemp->GetDigiNr(0) << " / " << hitTemp->GetDigiNr(1) << endl;
669  }
670  ((PndGemHit*)fHitsTemp->At(ihitTemp))->SetDetectorID(hitDetId);
671  fTNofHits++;
672  }
673 
674 
675  if( test2GemHit != 0 ) {
676  if ( fVerbose > 1 ) {
677  cout << "ihitTemp2 " << ihitTemp2 << " : already stored... " <<endl;
678  }
679  }
680  else {
681  pos.SetXYZ(xHit[1], yHit[1], zHit[1]);
682  dpos.SetXYZ(dx[1], dy[1], dz[1]);
683  hitDetId = hitTemp2->GetDetectorID() | kGemHit << 21;
684 
685  TString fromStr = "GEMDigi";
686  if ( fUseClusters ) fromStr = "GEMCluster";
687 
688  new ((*fHits)[nHits++]) PndGemHit(hitDetId,
689  pos, dpos,
690  hitTemp2->GetCharge(),
691  hitTemp2->GetTimeStamp(),
692  hitTemp2->GetDigiNr(0),hitTemp2->GetDigiNr(1),
693  hitTemp2->GetDr(), hitTemp2->GetDp(),
694  hitTemp2->GetRefIndex(),fromStr);
695 
696  if ( fVerbose > 1 ) {
697  cout << "ihitTemp2 "<< ihitTemp<<" xHit yHit zHit " << pos.X() << " / " << pos.Y() << " / " << pos.Z() << " stored as " << nHits << " from " << fromStr.Data() << "." << hitTemp2->GetDigiNr(0) << " / " << hitTemp2->GetDigiNr(1) << endl;
698  }
699  ((PndGemHit*)fHitsTemp->At(ihitTemp2))->SetDetectorID(hitDetId);
700  fTNofHits++;
701  }
702 
703  hitMatch = kTRUE;
704  hitMatch2 = kTRUE;
705  }
706  }
707  if ( hitMatch2 == kTRUE ) {
708  if ( fVerbose > 1 ) cout << "GOOD HIT" << endl;
709  //break;
710  }
711  else {
712  if ( fVerbose > 2 ) cout << "NO MATCHE TO THIS HIT ihitTemp2 = " << ihitTemp2<< endl;
713  }
714  }
715  if ( hitMatch == kFALSE )
716  if ( fVerbose > 1 ) cout << "NO MATCHES OF HIT ihitTemp = " << ihitTemp << endl;
717  }
718 }
719 // -------------------------------------------------------------------------
720 
721 
722 // // ----- Private method ConfirmHits --------------------------------------
723 // void PndGemFindHits::ConfirmHits() {
724 // PndGemSensor* sensor = NULL;
725 // PndGemHit* hit = NULL;
726 // PndGemCluster* cluster = NULL;
727 
728 // for ( Int_t ihit = 0 ; ihit < fHits->GetEntriesFast() ; ihit++ ) {
729 // hit = (PndGemHit*)fHits->At(ihit);
730 
731 // cout << "trying to confirm hit " << ihit << " at "
732 // << hit->GetX() << ","
733 // << hit->GetY() << ","
734 // << hit->GetZ() << " @ "
735 // << hit->GetTimeStamp() << endl;
736 
737 // sensor = (PndGemSensor*)fDigiPar->GetSensor(hit->GetStationNr(),
738 // (3-hit->GetSensorNr ()));
739 // Double_t hitProjX = hit->GetX()*(sensor->GetZ0()/hit->GetZ());
740 // Double_t hitProjY = hit->GetY()*(sensor->GetZ0()/hit->GetZ());
741 
742 // cout << " with hits at sensor @ z = " << sensor->GetZ0() << endl;
743 // Int_t hitConfirmed = 0;
744 
745 // }
746 // if ( hitConfirmed < 2 )
747 // cout << " NO CLUSTER MATCHES TO THIS HIT" << endl;
748 // else
749 // ((PndGemHit*)fHits->At(ihit))->SetCharge(10.);
750 // }
751 // }
752 // -------------------------------------------------------------------------
753 
754 // ----- Private method SortDigis --------------------------------------
756 
757  // Check input array
758  if ( ! fDigis ) {
759  cout << "-E- " << fName << "::SortDigis: No input array!" << endl;
760  return;
761  }
762 
763  // Clear sensor digi sets
764  map<PndGemSensor*, set<Int_t> >::iterator mapIt;
765  for (mapIt=fDigiMapF.begin(); mapIt!=fDigiMapF.end(); mapIt++)
766  ((*mapIt).second).clear();
767  for (mapIt=fDigiMapB.begin(); mapIt!=fDigiMapB.end(); mapIt++)
768  ((*mapIt).second).clear();
769 
770  // Fill digis into sets
771  PndGemDigi* digi = NULL;
772  //PndGemCluster* clust = NULL; //[R.K.01/2017]unused variable?
773  PndGemSensor* sensor = NULL;
774  Int_t stationNr = -1;
775  Int_t sensorNr = -1;
776  Int_t iSide = -1;
777  Int_t nDigis = fDigis->GetEntriesFast();
778 
779  fTNofDigis += nDigis;
780 
781  for (Int_t iDigi=0; iDigi<nDigis; iDigi++) {
782  digi = (PndGemDigi*) fDigis->At(iDigi);
783  //clust = (PndGemCluster*) fDigis->At(iDigi); //[R.K.01/2017]unused variable?
784  stationNr = digi->GetStationNr();
785  sensorNr = digi->GetSensorNr();
786  iSide = digi->GetSide();
787 
788  //Strange...this case to do continue..??
789  //it doesn't exist in apr13, but exists already on trunk r23148.
790  //(this code is taken from r24527). Comment out ok?
791  //if ( !( sensorNr == 1 && iSide == 0 ) && clust->GetClusterEnd()-clust->GetClusterBeg() <= 1 ) {
792  // // cout << "THE DIGI " << iDigi << ", DET ID = " << digi->GetDetectorId() << " in STATION " << stationNr << " SENSOR " << sensorNr << " SIDE " << iSide << " -> " << clust->GetClusterBeg() << " : " << clust->GetClusterEnd() << endl;
793  // continue;
794  //}
795 
796  // cout << "THE DIGI " << iDigi << ", DET ID = " << digi->GetDetectorId() << " in STATION " << stationNr << " SENSOR " << sensorNr << " SIDE " << iSide << endl;
797  // cout << "digi #" << iDigi+1 << " in station " << stationNr << " sensor " << sensorNr << " channel " << digi->GetChannelNr() << endl;
798 
799  sensor = fDigiPar->GetSensor(stationNr, sensorNr);
800  // cout << "Looking for station " << stationNr << " sensor " << sensorNr << " side " << iSide << endl;
801  if (iSide == 0 ) {
802  if ( fDigiMapF.find(sensor) == fDigiMapF.end() ) {
803  cerr << "-E- " << fName << "::SortDigits:: sensor " << sensorNr
804  << " of station " << stationNr
805  << " not found in digi scheme (F)!" << endl;
806  continue;
807  }
808  fDigiMapF[sensor].insert(iDigi);
809  }
810  else if (iSide == 1 ) {
811  if ( fDigiMapB.find(sensor) == fDigiMapB.end() ) {
812  cerr << "-E- " << fName << "::SortDigits:: sensor " << sensorNr
813  << " of station " << stationNr
814  << " not found in digi scheme (B)!" << endl;
815  continue;
816  }
817  fDigiMapB[sensor].insert(iDigi);
818  }
819  }
820 
821 }
822 // -------------------------------------------------------------------------
823 
824 
825 
826 
827 // ----- Private method FindHits ---------------------------------------
828 //one should add case of UseCluster to carry information of cluster Beg and End (or?)
830  set<Int_t>& fSet, set<Int_t>& bSet) {
831 
832  //Int_t iType = sensor->GetType(); //[R.K. 01/2017] unused variable?
833 
834  Double_t sigmaX = 1., sigmaY = 2.;
835 
836  Int_t iDigiF = -1;
837  Int_t iDigiB = -1;
838  Double_t iChanF = -1;
839  Double_t iChanB = -1;
840  Int_t nHits = fHits->GetEntriesFast();
841  Int_t nHitsInSensor = 0;
842  Double_t xHit;
843  Double_t yHit;
844  Double_t zHit;
845  Double_t dx, dy, dr, dp;
846  TVector3 pos, dpos;
847  PndGemDigi* digiF = NULL;
848  PndGemDigi* digiB = NULL;
849  //PndGemCluster* clusterB = NULL; //[R.K.01/2017]unused variable?
850 
851  set<Int_t>::iterator it1;
852  set<Int_t>::iterator it2;
853 
854  for (it1=fSet.begin(); it1!=fSet.end(); it1++) {
855  iDigiF = (*it1);
856  digiF = (PndGemDigi*) fDigis->At(iDigiF);
857  if ( ! digiF ) {
858  cout << "-W- " << GetName() << "::FindHits: Invalid digi index "
859  << iDigiF << " in front set of sensor "
860  << sensor->GetDetectorName() << endl;
861  continue;
862  }
863  iChanF = digiF->GetChannelNr();
864  for (it2=bSet.begin(); it2!=bSet.end(); it2++) {
865  iDigiB = (*it2);
866  digiB = (PndGemDigi*) fDigis->At(iDigiB);
867  //clusterB = (PndGemCluster*) fDigis->At(iDigiB); //[R.K.01/2017]unused variable?
868  if ( ! digiB ) {
869  cout << "-W- " << GetName() << "::FindHits: Invalid digi index "
870  << iDigiB << " in front set of sensor "
871  << sensor->GetDetectorName() << endl;
872  continue;
873  }
874  iChanB = digiB->GetChannelNr();
875  Int_t sensorDetId = sensor->Intersect(iChanF,iChanB,xHit,yHit,zHit,dx,dy,dr,dp);
876  // cout << "intersecting channels " << iChanF << " and " << iChanB << " gave " << sensorDetId << endl;
877  // cout << "got the following position from sensor " << sensorDetId << " : ("
878  // << xHit << ", " << yHit << ", " << zHit << ")" << endl;
879  if ( sensorDetId == -1 ) continue;
880 
881  Double_t rad = TMath::Sqrt(xHit*xHit+yHit*yHit);
882  if ( rad < sensor->GetInnerRadius() || rad > sensor->GetOuterRadius() ) {
883  cout << " point " << xHit << "," << yHit << " (" << rad << ") is still ok??? at station "
884  << sensor->GetStationNr() << "."
885  << sensor->GetSensorNr() << endl;
886  }
887 
888  //Don't we use information of Cluster pos Min Max to set error of Hit??
889  sigmaX = dx;
890  sigmaY = dy;
891 
892 // if ( fUseClusters ) {
893 // zHit += sensor->GetD()/2.;
894 // if ( sensor->GetType() )
895 // zHit -= sensor->GetD();
896 // }
897  // cout << "what about sensorz0 = " << sensor->GetZ0() << " or zHit = " << zHit << endl;
898 
899  pos.SetXYZ(xHit, yHit, zHit);
900  dpos.SetXYZ(sigmaX, sigmaY, sensor->GetD());
901 
902  // cout << "SETTING ERROR TO " << sigmaX << " " << sigmaY << " " << sensor->GetD() << " at " << zHit << endl;
903 
904  // if ( TMath::Abs(xHit) < 2. ) continue;
905 
906  Int_t refIndex = -1;
907  for ( Int_t irf1 = digiF->GetNIndices()-1 ; irf1 >= 0 ; irf1-- ) {
908  if ( digiF->GetLink(irf1).GetType() != fMCPointBranchId ) continue;
909  for ( Int_t irf2 = digiB->GetNIndices()-1 ; irf2 >= 0 ; irf2-- ) {
910  if ( digiB->GetLink(irf2).GetType() != fMCPointBranchId ) continue;
911  if ( digiF->GetIndex(irf1) == digiB->GetIndex(irf2) )
912  refIndex = digiF->GetIndex(irf1);
913  }
914  }
915 
916  Int_t hitDetId = sensorDetId | kGemHit << 21;
917 
918  // cout << nHits << " : " << pos.X() << " " << pos.Y() << " " << pos.Z() << endl;
919 
920  // fChargeCorrHist->Fill(digiF->GetCharge(),digiB->GetCharge());
921  // fChargeDiffHist->Fill((digiF->GetCharge()+digiB->GetCharge())/2.,
922  // (digiF->GetCharge()-digiB->GetCharge()));
923 
924  // if ( TMath::Abs(digiF->GetCharge()-digiB->GetCharge()) > 0.02+(digiF->GetCharge()+digiB->GetCharge())/200. ||
925 
926  // How critical are they??? may I open more??
927  if ( TMath::Abs(digiF->GetCharge()-digiB->GetCharge()) > 5. ||
928  TMath::Abs(digiF->GetTimeStamp()-digiB->GetTimeStamp()) > 5. )
929  continue;
930 
931  TString fromStr = "GEMDigi";
932  if ( fUseClusters ) fromStr = "GEMCluster";
933  new ((*fHits)[nHits++]) PndGemHit(hitDetId, pos, dpos,
934  (digiF->GetCharge()+digiB->GetCharge())/2.,
935  (digiF->GetTimeStamp()+digiB->GetTimeStamp())/2.,
936  iDigiF, iDigiB,
937  dr, dp, refIndex,fromStr);
938 
939  fTNofHits++;
940  nHitsInSensor++;
941  }
942  }
943 
944  // return 0;
945  return nHitsInSensor;
946 }
947 // -------------------------------------------------------------------------
948 
949 
950 // ----- Private method FindHits ---------------------------------------
951 //create hits for all combination and store in fHitsTemp..
952 //carry information of cluster width to dpos..
954  set<Int_t>& fSet, set<Int_t>& bSet) {
955  if ( fVerbose > 0 ) cout << "-I- PndGemFindHits::FindHits2(). Station " << sensor->GetStationNr() << " sensor " << sensor->GetSensorNr() << endl;
956 
957  //Int_t iType = sensor->GetType(); //[R.K. 01/2017] unused variable?
958 
959  Double_t sigmaX = 1., sigmaY = 2., sigmaZ;
960  Double_t sigmaR, sigmaP;
961 
962  Int_t iDigiF = -1;
963  Int_t iDigiB = -1;
964  Double_t iChanF = -1;
965  Double_t iChanB = -1;
966  Double_t iChanFMn,iChanFMx;
967  Double_t iChanBMn,iChanBMx;
968  Int_t nHitsTemp = fHitsTemp->GetEntriesFast();
969  Int_t nHitsInSensor = 0;
970  Double_t xHit;
971  Double_t yHit;
972  Double_t zHit;
973  Double_t dx, dy, dr, dp;
974  TVector3 pos, dpos;
975  PndGemDigi* digiF = NULL;
976  PndGemDigi* digiB = NULL;
977  PndGemCluster* clusterF = NULL;
978  PndGemCluster* clusterB = NULL;
979 
980  Double_t xMnMn,xMxMx,xMnMx,xMxMn;
981  Double_t yMnMn,yMxMx,yMnMx,yMxMn;
982 
983  set<Int_t>::iterator it1;
984  set<Int_t>::iterator it2;
985 
986  for (it1=fSet.begin(); it1!=fSet.end(); it1++) {
987  iDigiF = (*it1);
988  digiF = (PndGemDigi*) fDigis->At(iDigiF);
989  clusterF = (PndGemCluster*) fDigis->At(iDigiF);
990  if ( ! digiF ) {
991  cout << "-W- " << GetName() << "::FindHits: Invalid digi index "
992  << iDigiF << " in front set of sensor "
993  << sensor->GetDetectorName() << endl;
994  continue;
995  }
996  iChanF = digiF->GetChannelNr();
997  if ( fUseClusters ) {
998  iChanFMn = clusterF->GetClusterBeg();
999  iChanFMx = clusterF->GetClusterEnd();
1000  }
1001  for (it2=bSet.begin(); it2!=bSet.end(); it2++) {
1002  iDigiB = (*it2);
1003  digiB = (PndGemDigi*) fDigis->At(iDigiB);
1004  clusterB = (PndGemCluster*) fDigis->At(iDigiB);
1005  if ( ! digiB ) {
1006  cout << "-W- " << GetName() << "::FindHits: Invalid digi index "
1007  << iDigiB << " in front set of sensor "
1008  << sensor->GetDetectorName() << endl;
1009  continue;
1010  }
1011  iChanB = digiB->GetChannelNr();
1012  if ( fUseClusters ) {
1013  iChanBMn = clusterB->GetClusterBeg();
1014  iChanBMx = clusterB->GetClusterEnd();
1015  }
1016 
1017  Int_t sensorDetId = sensor->Intersect(iChanF,iChanB,xHit,yHit,zHit,dx,dy,dr,dp);
1018  if ( fVerbose > 1 ) {
1019  cout << "Intersect chF " << iChanF << " chB " << iChanB << " detID " << sensorDetId << endl;
1020  cout << "at (" << xHit << ", " << yHit << ", " << zHit << ")"
1021  << "dx " << dx << " dy " << dy << " dr " << dr << " dp " << dp
1022  << endl;
1023  }
1024  if ( sensorDetId == -1 ){
1025  if ( fVerbose > 2 ) {
1026  cout << "no intersect.. skip " <<endl;
1027  }
1028  continue;
1029  }
1030 
1031  pos.SetXYZ(xHit, yHit, zHit);
1032  Double_t rad = TMath::Sqrt(xHit*xHit+yHit*yHit);
1033  if ( rad < sensor->GetInnerRadius() || rad > sensor->GetOuterRadius() ) {
1034  cout << "-W- !! point " << xHit << "," << yHit << " (" << rad << ") is still ok??? at station "
1035  << sensor->GetStationNr() << "." << sensor->GetSensorNr() << endl;
1036  }
1037  //store resolution of the hit pos
1038  sigmaX = dx;//basically reso of single strip
1039  sigmaY = dy;//basically reso of single strip
1040  sigmaZ = sensor->GetD();//
1041  sigmaR = dr;//basically reso of single strip
1042  sigmaP = dp;//basically reso of single strip
1043 
1044  //take channel Beg and End for dpos when Cluster is used
1045  //Diagonal lines of the overlaped part are used to extract dx dy and longer one is taken
1046  //no check on dr dp, no comparison with original resolution sigmaX sigmaY.. ok?
1047  if ( fUseClusters ) {
1048  Int_t testMnMn = sensor->Intersect(iChanFMn,iChanBMn,xMnMn,yMnMn,zHit,dx,dy,dr,dp);
1049  Int_t testMxMx = sensor->Intersect(iChanFMx,iChanBMx,xMxMx,yMxMx,zHit,dx,dy,dr,dp);
1050  Int_t testMnMx = sensor->Intersect(iChanFMn,iChanBMx,xMnMx,yMnMx,zHit,dx,dy,dr,dp);
1051  Int_t testMxMn = sensor->Intersect(iChanFMx,iChanBMn,xMxMn,yMxMn,zHit,dx,dy,dr,dp);
1052  Double_t dx1 = -1;
1053  Double_t dx2 = -1;
1054  Double_t dy1 = -1;
1055  Double_t dy2 = -1;
1056  if( testMnMn != -1 && testMxMx != -1 ) {
1057  dx1 = TMath::Abs( xMnMn - xMxMx );
1058  dy1 = TMath::Abs( yMnMn - yMxMx );
1059  if( testMnMx != -1 && testMxMn != -1 ) {
1060  dx2 = TMath::Abs( xMnMx - xMxMn );
1061  dy2 = TMath::Abs( yMnMx - yMxMn );
1062  sigmaX = TMath::Max( dx1, dx2 ) / TMath::Sqrt(12.);
1063  sigmaY = TMath::Max( dy1, dy2 ) / TMath::Sqrt(12.);
1064  if ( fVerbose > 1 ) {
1065  cout << "fUseCluster is TRUE: dx dy of clusters " << sigmaX << " " << sigmaY << " are calculated for Hits"<< endl;
1066  }
1067  }
1068  }
1069  }
1070  dpos.SetXYZ(sigmaX, sigmaY, sigmaZ);
1071 
1072 
1073  Int_t refIndex = -1;
1074  for ( Int_t irf1 = digiF->GetNIndices()-1 ; irf1 >= 0 ; irf1-- ) {
1075  if ( digiF->GetLink(irf1).GetType() != fMCPointBranchId ) continue;
1076  for ( Int_t irf2 = digiB->GetNIndices()-1 ; irf2 >= 0 ; irf2-- ) {
1077  if ( digiB->GetLink(irf2).GetType() != fMCPointBranchId ) continue;
1078  if ( digiF->GetIndex(irf1) == digiB->GetIndex(irf2) )
1079  refIndex = digiF->GetIndex(irf1);
1080  }
1081  }
1082 
1083  // Hit bit kGemHit will be put when the Hit matches on both sensor.. or problem?
1084  Int_t hitDetId = sensorDetId;
1085  // Int_t hitDetId = sensorDetId | kGemHit << 21;
1086 
1087  // How critical are they??? may I open more??
1088  // SO FAR COMMENTED !! ok??
1089  // if ( TMath::Abs(digiF->GetCharge()-digiB->GetCharge()) > 5. ||
1090  // TMath::Abs(digiF->GetTimeStamp()-digiB->GetTimeStamp()) > 5. )
1091  // continue;
1092 
1093  TString fromStr = "GEMDigi";
1094  Double_t hitTimeStamp = (digiF->GetTimeStamp()+digiB->GetTimeStamp())/2.;
1095  Double_t hitCharge = (digiF->GetCharge() +digiB->GetCharge() )/2.;
1096  if ( fUseClusters ) {
1097  fromStr = "GEMCluster";
1098  hitTimeStamp = (clusterF->GetTimeStamp()+clusterB->GetTimeStamp())/2.;
1099  hitCharge = (clusterF->GetCharge() +clusterB->GetCharge() )/2.;
1100  }
1101  new ((*fHitsTemp)[nHitsTemp++]) PndGemHit(hitDetId, pos, dpos,
1102  hitCharge,
1103  hitTimeStamp,
1104  iDigiF, iDigiB,
1105  sigmaR, sigmaP, refIndex, fromStr);
1106 
1107  fTNofHitsTemp++;
1108  nHitsInSensor++;
1109  }
1110  }
1111 
1112  // return 0;
1113  return nHitsInSensor;
1114 }
1115 // -------------------------------------------------------------------------
1116 
1117 
1118 // ----- Public method Finish ------------------------------------------
1120  if ( fHits ) fHits->Clear();
1121 
1122 // fChargeDiffHist->Draw("colz");
1123 
1124  cout << "-------------------- " << fName.Data() << " : Summary -----------------------" << endl;
1125  cout << " Events: " << setw(10) << fTNofEvents << endl;
1126  cout << " Digis: " << setw(10) << fTNofDigis << " ( " << (Double_t)fTNofDigis/((Double_t)fTNofEvents) << " per event )" << endl;
1127  cout << " HitsTemp: " << setw(10) << fTNofHitsTemp << " ( " << (Double_t)fTNofHitsTemp /((Double_t)fTNofEvents) << " per event )" << endl;
1128  cout << " Hits: " << setw(10) << fTNofHits << " ( " << (Double_t)fTNofHits /((Double_t)fTNofEvents) << " per event )" << endl;
1129  cout << " --> ( " << (Double_t)fTNofHits /((Double_t)fTNofEvents)/((Double_t)fDigiPar->GetNSensors()) << " per sensor )" << endl;
1130  cout << " --> ( " << (Double_t)fTNofHits /((Double_t)fTNofDigis ) << " per digi )" << endl;
1131  cout << "---------------------------------------------------------------------" << endl;
1132 
1133  cout << " >>> HF >>> prep time = " << fPrepTime << "s (get data from input)" << endl;
1134  cout << " >>> HF >>> sort time = " << fSortTime-fPrepTime << "s (sort clusters, " << fSortTime << ")" << endl;
1135  cout << " >>> HF >>> create time = " << fCreateTime-fSortTime << "s (create hits, " << fCreateTime << ")" << endl;
1136  cout << " >>> HF >>> confirm time = " << fConfirmTime-fCreateTime << "s (confirm hits, " << fConfirmTime << ")" << endl;
1137  cout << " >>> HF >>> activ. time = " << fActivateTime-fConfirmTime << "s (activate digis, " << fActivateTime << ")" << endl;
1138  cout << " >>> HF >>> all time = " << fAllTime-fActivateTime << "s (all time spent in Exec, " << fAllTime << ")" << endl;
1139  cout << "---------------------------------------------------------------------" << endl;
1140 
1141 }
1142 // -------------------------------------------------------------------------
1143 
1144 
1146 
1147 
1148 
1149 
1150 
1151 
1152 
1153 
1154 
1155 
1156 
1157 
1158 
TVector3 pos
Double_t fActivateTime
int fVerbose
Definition: poormantracks.C:24
TClonesArray * digi
Int_t GetSensorNr() const
Definition: PndGemCluster.h:90
double dy
Int_t GetClusterBeg() const
Definition: PndGemCluster.h:95
Int_t CreateSensorMonitor(const PndGemSensor &tempSensor)
Int_t run
Definition: autocutx.C:47
Int_t Intersect(Double_t iFStrip, Double_t iBStrip, Double_t &xCross, Double_t &yCross, Double_t &zCross)
TTree * b
Int_t GetStationNr() const
Definition: PndGemDigi.h:84
TString GetDetectorName() const
Definition: PndGemSensor.h:86
PndTransMap * map
Definition: sim_emc_apd.C:99
Int_t GetStationNr() const
Definition: PndGemStation.h:57
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
Double_t fAllTime
Int_t GetSide() const
Definition: PndGemCluster.h:92
PndGemMonitor * fMonitor
TClonesArray * fDigis
void SetPersistency(Bool_t val=kTRUE)
TGeoVolume * sensor
Digitization Parameter Class for GEM part.
Definition: PndGemDigiPar.h:31
virtual void SetParContainers()
virtual void Exec(Option_t *opt)
PndGemSensor * GetSensor(Int_t iSensor)
Definition: PndGemStation.h:63
Double_t ChannelLastActiveAt(Int_t statNr, Int_t sensNr, Int_t sideId, Int_t chanNr)
Int_t GetStationNr() const
Definition: PndGemSensor.h:93
std::map< PndGemSensor *, std::set< Int_t > > fDigiMapF
Int_t GetNIndices()
Definition: PndGemDigi.h:102
Int_t GetSensorNr() const
Definition: PndGemHit.h:83
std::map< PndGemSensor *, std::set< Int_t > > fDigiMapB
static T Abs(const T &x)
Definition: PndCAMath.h:39
Int_t GetSensorNr() const
Definition: PndGemDigi.h:86
Int_t refIndex
Definition: checkmomentum.C:30
Int_t a
Definition: anaLmdDigi.C:126
PndGemSensor * GetSensor(Int_t stationNr, Int_t sensorNr)
Int_t GetNSensors() const
Definition: PndGemStation.h:60
int nHits
Definition: RiemannTest.C:16
TClonesArray * fHitsTemp
Int_t GetChannel(Double_t x, Double_t y, Int_t iSide)
Double_t
void EnableCluster(Int_t eventNr, Int_t clusterNr, PndGemCluster *tempCluster)
Int_t FindHits2(PndGemSensor *sensor, std::set< Int_t > &fSet, std::set< Int_t > &bSet)
virtual void Finish()
static PndGemMonitor * Instance()
Double_t GetDp() const
Definition: PndGemHit.h:76
Double_t GetChannelNr() const
Definition: PndGemDigi.h:80
Double_t fPrepTime
Double_t GetCharge() const
Definition: PndGemHit.h:69
Int_t GetNStations()
Definition: PndGemDigiPar.h:45
TString name
Int_t GetStationNr() const
Definition: PndGemHit.h:81
Int_t GetClusterEnd() const
Definition: PndGemCluster.h:96
Int_t GetNSensors()
Definition: PndGemDigiPar.h:46
double dx
Int_t fMCPointBranchId
Double_t fSortTime
static T Max(const T &x, const T &y)
Definition: PndCAMath.h:36
Double_t GetD() const
Definition: PndGemSensor.h:104
Double_t fConfirmTime
Int_t GetSensorNr() const
Definition: PndGemSensor.h:95
Int_t GetIndex(int i=0) const
Definition: PndGemDigi.h:103
TClonesArray * fHits
ClassImp(PndAnaContFact)
Int_t GetDigiNr(Int_t iside) const
Definition: PndGemHit.h:77
Int_t iVerbose
Double_t GetCharge() const
Definition: PndGemDigi.h:91
Int_t GetDetectorId() const
Definition: PndGemSensor.h:89
Double_t GetZ0() const
Definition: PndGemSensor.h:100
Double_t GetOuterRadius() const
Definition: PndGemSensor.h:103
PndSdsMCPoint * hit
Definition: anasim.C:70
Double_t GetDr() const
Definition: PndGemHit.h:75
Int_t FindHits(PndGemSensor *sensor, std::set< Int_t > &fSet, std::set< Int_t > &bSet)
virtual ~PndGemFindHits()
PndGemDigiPar * fDigiPar
/** GEM monitor **/
virtual InitStatus Init()
TStopwatch fTimer
Int_t GetSide() const
Definition: PndGemDigi.h:88
Double_t GetCharge() const
Definition: PndGemCluster.h:97
Int_t GetStationNr() const
Definition: PndGemCluster.h:88
Double_t fCreateTime
virtual InitStatus ReInit()
PndGemStation * GetStation(Int_t iStation)
Bool_t fTimeOrderedDigi
Double_t fHitWindow