FairRoot/PandaRoot
PndDrcHitFinder.cxx
Go to the documentation of this file.
1 // -------------------------------------------------------------------------
2 // ----- PndDrcHitFinder source file -----
3 // -------------------------------------------------------------------------
4 
5 
6 #include "TClonesArray.h"
7 #include "TArrayD.h"
8 #include "TGeoManager.h"
9 
10 #include "FairRootManager.h"
11 #include "PndDrcHitFinder.h"
12 #include "FairRun.h"
13 #include "FairRuntimeDb.h"
14 #include "FairGeoNode.h"
15 #include "FairGeoVector.h"
16 #include "FairRunAna.h"
17 #include "FairEventHeader.h"
18 
19 //#include "PndDrcDigi.h"
20 
21 #include "PndDetectorList.h"
22 
23 // ----- Default constructor -------------------------------------------
25 PndPersistencyTask("DrcHitFinder", 1)
26 {
27  SetPersistency(kTRUE);
28  fPixelHits = 0;
29  fEventNr = 0;
30  fPixelFactor = 1;
31 
32  if(fVerbose>0) Info("PndDrcHitFinder","DRC Hit Finder created, Parameters will be taken from RTDB");
33  fGeoH=NULL;
34  fGeo = new PndGeoDrc();
35  fDigiArray = NULL;
36  fPdHitArray = NULL;
37  fMCEventHeader = NULL;
38  fStopFunctor = NULL;
39 
40  fPixelSize = fGeo->PixelSize(); //pixel size;
41  fNpix = fGeo->Npixels(); //pixel columns
44  fPixelStep = fMcpActiveArea/fNpix;//fPixelSize + 0.5*fPixelGap;
45 
46 }
47 // -------------------------------------------------------------------------
48 
50  PndPersistencyTask("DrcHitFinder", iVerbose)
51 {
52  SetPersistency(kTRUE);
53  fPixelHits = 0;
54  fEventNr = 0;
55  fPixelFactor = 1;
56  fGeoH = NULL;
57  fGeo = new PndGeoDrc();
58  fDigiPixelMCInfo = kFALSE;
59  if(fVerbose>0) Info("PndDrcHitFinder","DrcHitFinder created, Parameters will be taken from RTDB");
60 
61  fDigiArray = NULL;
62  fPdHitArray = NULL;
63  fMCEventHeader = NULL;
64  fStopFunctor = NULL;
65 
66  fPixelSize = fGeo->PixelSize(); //pixel size;
67  fNpix = fGeo->Npixels(); //pixel rows in one FE
70  fPixelStep = fPixelSize + 0.5*fPixelGap;
71 }
72 
74 PndPersistencyTask(name, iVerbose)
75 {
76  fPixelHits = 0;
77  fEventNr = 0;
78  fPixelFactor = 1;
79  fGeoH = NULL;
80  fGeo = new PndGeoDrc();
81  fDigiPixelMCInfo = kFALSE;
82  if(fVerbose>0) Info("PndDrcHitFinder","%s created, Parameters will be taken from RTDB",name);
83 
84  fDigiArray = NULL;
85  fPdHitArray = NULL;
86  fMCEventHeader = NULL;
87 
88  fPixelSize = fGeo->PixelSize(); //pixel size;
89  fNpix = fGeo->Npixels(); //pixel rows in one FE
92  fPixelStep = fPixelSize + 0.5*fPixelGap;
93 }
94 
95 
96 // ----- Destructor ----------------------------------------------------
98  if (fGeo) delete fGeo;
99  if (fGeoH) delete fGeoH;
100 }
101 // -------------------------------------------------------------------------
102 
103 // ----- Initialization of Parameter Containers -------------------------
105  if ( fGeoH == NULL ) fGeoH = PndGeoHandling::Instance();
108  if(fVerbose>1) Info("SetParContainers","done.");
109  return;
110 }
111 
114  return kSUCCESS;
115 }
116 
117 // ----- Public method Init --------------------------------------------
119  //FairRun* ana = FairRun::Instance(); //[R.K. 01/2017] unused variable?
120  FairRootManager* ioman = FairRootManager::Instance();
121  if ( ! ioman )
122  {
123  std::cout << "-E- PndDrcHitFinder::Init: "
124  << "RootManager not instantiated!" << std::endl;
125  return kFATAL;
126  }
127 
128  fInBranchName = "DrcDigi";
129  // Get input array
130  fDigiArray = (TClonesArray*) ioman->GetObject(fInBranchName);
131  if ( ! fDigiArray )
132  {
133  std::cout << "-W- PndDrcHitFinder::Init: "
134  << "No DrcDigi array!" << std::endl;
135  return kERROR;
136  }
137 
138  // Create and register output array
139  fPdHitArray = new TClonesArray("PndDrcPDHit");
140  ioman->Register("DrcPDHit", "Drc", fPdHitArray, GetPersistency());
141 
142  fGapFunctor = new TimeGap();
143  fStopFunctor = new StopTime();
144  return kSUCCESS;
145 }
146 
147 // ----- Public method Exec --------------------------------------------
148 void PndDrcHitFinder::Exec(Option_t*){
149 
150  if(fVerbose>3) Info("Exec","Start");
151  if (!fPdHitArray) Fatal("Exec", "No PdHitArray");
152  fPdHitArray->Clear();
153 
154  Double_t etime = FairRootManager::Instance()->GetEventTime();
155  if (FairRunAna::Instance()->IsTimeStamp()){
156  fDigiArray = FairRootManager::Instance()->GetData(fInBranchName, fStopFunctor, etime, fStopFunctor, etime+100);
157  }
158 
159  Int_t nDigis = fDigiArray->GetEntriesFast();
160  if(fVerbose>1) std::cout<<"-I- PndDrcHitFinder: Event # "<< fEventNr<<" has "<<nDigis<<" digis."<< std::endl;
161  else if(fVerbose==1 && fEventNr%1000==0) std::cout<<"-I- PndDrcHitFinder: Event # "<< fEventNr<<" has "<<nDigis<<" digis."<< std::endl;
162 
163  Int_t detID = 0, mcpID = 0, pixelID = 0;
164  TVector3 HitPosGlobal, HitPosLocal, dPosHit;
165  Double_t hitTime = 0.;
166 
167  for (Int_t iDigi = 0; iDigi < nDigis; iDigi++){
168  fDigi = (PndDrcDigi*) fDigiArray->At(iDigi);
169  detID = fDigi->GetDetectorId();
170  pixelID = detID - 100*(Int_t)TMath::Floor((Double_t)detID/100.);
171  mcpID = detID/100;
172  hitTime = fDigi->GetTime();
173 
174  // the pixel number shows local coordinates of the hit:
175  HitPosLocal.SetXYZ(fPixelStep*((Double_t)(pixelID % fNpix) - (Double_t)(fNpix/2) + 0.5),
176  fPixelStep*(TMath::Floor(((Double_t)pixelID)/((Double_t)fNpix)) - (Double_t)(fNpix/2) + 0.5), 0.);
177  Int_t sensorId = fDigi->GetSensorId()/fPixelFactor;
178  if(fPixelFactor==2) { //double pixels
179  sensorId++;
180  if(fDigi->GetSensorId()%2) HitPosLocal.SetX(HitPosLocal.X()+fPixelSize/2.);
181  else HitPosLocal.SetX(HitPosLocal.X()-fPixelSize/2.);
182  }
183 
184  // local coordinates of the hit on the MCP are translated into global ones as the following:
185  HitPosGlobal = fGeoH->LocalToMasterShortId(HitPosLocal, mcpID);
186  dPosHit.SetXYZ(fPixelSize/2., fPixelSize/2., 0.);
187  if(fPixelFactor==2) dPosHit.SetXYZ(fPixelSize, fPixelSize/2., 0.);
188 
189 
190  if(fDigi->GetTimeStamp()!=etime+hitTime) hitTime = fDigi->GetTimeStamp() - etime;
191  PndDrcPDHit pdhit = PndDrcPDHit(detID, sensorId , HitPosGlobal, dPosHit, hitTime, 0., iDigi);
192  pdhit.SetPdgCode(fDigi->GetPdgCode());
193  pdhit.SetTimeStamp(fDigi->GetTimeStamp());
194  pdhit.SetTimeAtBar(fDigi->GetTimeAtBar());
195  pdhit.SetLink(fDigi->GetLink(0)); // MCTrack
196  pdhit.AddLink(fDigi->GetLink(1)); // DrcPDPoint
197  pdhit.AddLink(FairLink(-1,fEventNr, "DrcDigi", iDigi));
198  //((FairMultiLinkedData)pdhit).Print(); std::cout<<std::endl;
199  new((*fPdHitArray)[fPdHitArray->GetEntriesFast()]) PndDrcPDHit(pdhit);
200  }
201 
202  fEventNr++;
203  if(fVerbose>3) Info("Exec","Loop MC points");
204 
205  fPdHitArray->Sort();
206  fDigiArray->Delete();
207 }
208 
209 
210 // -------------------------------------------------------------------------
212 {
213  FinishEvents();
214 }
215 
216 // -------------------------------------------------------------------------
218 {
219 }
220 
void SetTimeAtBar(Double_t TimeAtBar)
Definition: PndDrcPDHit.cxx:75
int fVerbose
Definition: poormantracks.C:24
PndGeoDrc * fGeo
virtual void Exec(Option_t *opt)
virtual void SetParContainers()
virtual InitStatus Init()
BinaryFunctor * fGapFunctor
TClonesArray * fPdHitArray
void SetPersistency(Bool_t val=kTRUE)
virtual void FinishEvent()
FairMCEventHeader * fMCEventHeader
TClonesArray * fDigiArray
virtual void FinishTask()
virtual InitStatus ReInit()
Double_t
Int_t GetPdgCode() const
Definition: PndDrcDigi.h:81
Double_t McpActiveArea()
Definition: PndGeoDrc.h:166
BinaryFunctor * fStopFunctor
PndGeoHandling * fGeoH
PndDrcDigi * fDigi
static PndGeoHandling * Instance()
TString name
void SetVerbose(Int_t v)
Int_t GetSensorId() const
Definition: PndDrcDigi.h:91
ClassImp(PndAnaContFact)
TVector3 LocalToMasterShortId(const TVector3 &local, const Int_t &shortId)
Int_t iVerbose
Double_t fMcpActiveArea
virtual ~PndDrcHitFinder()
Int_t Npixels()
Definition: PndGeoDrc.h:172
Int_t GetDetectorId() const
Definition: PndDrcDigi.h:90
void SetPdgCode(Int_t Pdg)
Definition: PndDrcPDHit.cxx:66
virtual void SetParContainers()
Double_t PixelSize()
Definition: PndGeoDrc.h:175
Double_t GetTimeAtBar() const
Definition: PndDrcDigi.h:84
Double_t GetTime() const
Definition: PndDrcDigi.h:96