FairRoot
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
PixelFindHitsTask.cxx
Go to the documentation of this file.
1 /********************************************************************************
2  * Copyright (C) 2014 GSI Helmholtzzentrum fuer Schwerionenforschung GmbH *
3  * *
4  * This software is distributed under the terms of the *
5  * GNU Lesser General Public Licence (LGPL) version 3, *
6  * copied verbatim in the file "LICENSE" *
7  ********************************************************************************/
8 /*
9  * PixelHit.h
10  *
11  * Created on: 18.02.2016
12  * Author: R. Karabowicz
13  */
14 
15 #include "PixelFindHitsTask.h"
16 
17 #include "FairGeoParSet.h"
18 #include "PixelDigi.h"
19 #include "PixelDigiPar.h"
20 #include "PixelHit.h"
21 
22 #include <FairMQLogger.h>
23 #include <TClonesArray.h>
24 #include <TGeoBBox.h>
25 #include <TGeoManager.h>
26 #include <TGeoNode.h>
27 #include <TGeoVolume.h>
28 #include <TList.h>
29 #include <TMath.h>
30 #include <TString.h>
31 #include <TVector3.h>
32 
34  : fGeoParSet(nullptr)
35  , fNDigis(0)
36  , fNHits(0)
37  , fTNofEvents(0)
38  , fTNofDigis(0)
39  , fTNofHits(0)
40  , fFeCols(0)
41  , fFeRows(0)
42  , fMaxFEperCol(0)
43  , fPitchX(0.)
44  , fPitchY(0.)
45 {
46  // Reset();
47 }
48 
50 {
51  // Reset();
52 }
53 
54 void PixelFindHitsTask::Reset(TClonesArray* hits)
55 {
56  fNDigis = fNHits = 0;
57  if (hits)
58  hits->Delete();
59 }
60 
61 void PixelFindHitsTask::Finish()
62 {
63  // if (fDigis) fDigis->Delete();
64 
65  LOG(info) << "-------------------- PixelFindHitsTask : Summary ------------------------";
66  LOG(info) << " Events: " << fTNofEvents;
67  LOG(info) << " Digis: " << fTNofDigis << " ("
68  << static_cast<Double_t>(fTNofDigis) / (static_cast<Double_t>(fTNofEvents)) << " per event)";
69  LOG(info) << " Hits: " << fTNofHits << " ("
70  << static_cast<Double_t>(fTNofHits) / (static_cast<Double_t>(fTNofEvents)) << " per event)";
71  LOG(info) << "---------------------------------------------------------------------";
72 }
73 
75 {
76  fFeCols = digipar->GetFECols();
77  fFeRows = digipar->GetFERows();
78  fMaxFEperCol = digipar->GetMaxFEperCol();
79  fPitchX = digipar->GetXPitch();
80  fPitchY = digipar->GetYPitch();
81  fGeoParSet = geopar;
82  LOG(info) << "PixelFindHitsTask::Init"
83  << " fFeCols=" << fFeCols << " fFeRows=" << fFeRows << " fMaxFEperCol=" << fMaxFEperCol
84  << " fPitchX=" << fPitchX << " fPitchY=" << fPitchX;
85  LOG(info) << "geopar->printParams()";
86  geopar->printParams();
87  LOG(info) << "fGeoParSet->printParams()";
88  fGeoParSet->printParams();
89 }
90 
91 void PixelFindHitsTask::Exec(TClonesArray* digis, TClonesArray* hits)
92 {
93  Reset(hits);
94 
95  LOG(info) << "PixelFindHits::Exec(TCA) EVENT " << fTNofEvents;
96 
97  fTNofEvents++;
98  LOG(debug) << "PixelFindHits::Exec() ok 0 ";
99  fNDigis = digis->GetEntriesFast();
100  LOG(debug) << "PixelFindHits::Exec() fNDigis = " << fNDigis;
101  fTNofDigis += fNDigis;
102  //*
103  for (Int_t iDigi = 0; iDigi < fNDigis; iDigi++) {
104  PixelDigi* currentDigi = (PixelDigi*)digis->At(iDigi);
105 
106  Int_t detId = currentDigi->GetDetectorID();
107  TString nodeName = Form("/cave/Pixel%d_%d", detId / 256, detId % 256);
108  LOG(debug) << "PixelFindHits::Exec() ok 1 node name = " << nodeName.Data();
109  fGeoParSet->GetGeometry()->cd(nodeName.Data());
110  LOG(debug) << "PixelFindHits::Exec() ok 2 ";
111  TGeoNode* curNode = fGeoParSet->GetGeometry()->GetCurrentNode();
112  LOG(debug) << "PixelFindHits::Exec() ok 3 ";
113  // TGeoMatrix* matrix = curNode->GetMatrix();
114 
115  TGeoVolume* actVolume = fGeoParSet->GetGeometry()->GetCurrentVolume();
116  TGeoBBox* actBox = (TGeoBBox*)(actVolume->GetShape());
117 
118  Int_t feId = currentDigi->GetFeID();
119  Int_t col = currentDigi->GetCol();
120  Int_t row = currentDigi->GetRow();
121 
122  Double_t locPosCalc[3];
123  locPosCalc[0] = (((feId - 1) / fMaxFEperCol) * fFeCols + col + 0.5) * fPitchX;
124  locPosCalc[1] = (((feId - 1) % fMaxFEperCol) * fFeRows + row + 0.5) * fPitchY;
125  locPosCalc[2] = 0.;
126 
127  locPosCalc[0] -= actBox->GetDX();
128  locPosCalc[1] -= actBox->GetDY();
129 
130  Double_t globPos[3];
131 
132  curNode->LocalToMaster(locPosCalc, globPos);
133 
134  LOG(debug) << "HIT ON " << detId << " POSITION: " << locPosCalc[0] << " / " << locPosCalc[1];
135  LOG(debug) << "GLOB HIT " << detId << " POSITION: " << globPos[0] << " / " << globPos[1] << " / "
136  << globPos[2];
137 
138  TVector3 pos(globPos[0], globPos[1], globPos[2]);
139  TVector3 posErr(fPitchX / TMath::Sqrt(12.), fPitchY / TMath::Sqrt(12.), actBox->GetDZ());
140 
141  new ((*hits)[fNHits]) PixelHit(detId, currentDigi->GetIndex(), pos, posErr);
142 
143  fNHits++;
144  }
145 
146  fTNofHits += fNHits;
147  // return fHits;
148 }
149 
150 void PixelFindHitsTask::Exec(TList* list, TClonesArray* hits)
151 {
152  Reset(hits);
153 
154  LOG(info) << "PixelFindHits::Exec(TList) EVENT " << fTNofEvents;
155 
156  TClonesArray* digis = (TClonesArray*)list->FindObject("PixelDigis");
157  if (!digis) {
158  LOG(error) << "OOOPS!!! no digis array";
159  return;
160  }
161 
162  fTNofEvents++;
163  LOG(debug) << "PixelFindHits::Exec() ok 0 ";
164  fNDigis = digis->GetEntriesFast();
165  LOG(debug) << "PixelFindHits::Exec() fNDigis = " << fNDigis;
166  fTNofDigis += fNDigis;
167  //*
168  for (Int_t iDigi = 0; iDigi < fNDigis; iDigi++) {
169  PixelDigi* currentDigi = static_cast<PixelDigi*>(digis->At(iDigi));
170 
171  Int_t detId = currentDigi->GetDetectorID();
172  TString nodeName = Form("/cave/Pixel%d_%d", detId / 256, detId % 256);
173  LOG(debug) << "PixelFindHits::Exec() ok 1 node name = " << nodeName.Data();
174  fGeoParSet->GetGeometry()->cd(nodeName.Data());
175  LOG(debug) << "PixelFindHits::Exec() ok 2 ";
176  TGeoNode* curNode = fGeoParSet->GetGeometry()->GetCurrentNode();
177  LOG(debug) << "PixelFindHits::Exec() ok 3 ";
178  // TGeoMatrix* matrix = curNode->GetMatrix();
179 
180  TGeoVolume* actVolume = fGeoParSet->GetGeometry()->GetCurrentVolume();
181  TGeoBBox* actBox = static_cast<TGeoBBox*>(actVolume->GetShape());
182 
183  Int_t feId = currentDigi->GetFeID();
184  Int_t col = currentDigi->GetCol();
185  Int_t row = currentDigi->GetRow();
186 
187  Double_t locPosCalc[3];
188  locPosCalc[0] = (((feId - 1) / fMaxFEperCol) * fFeCols + col + 0.5) * fPitchX;
189  locPosCalc[1] = (((feId - 1) % fMaxFEperCol) * fFeRows + row + 0.5) * fPitchY;
190  locPosCalc[2] = 0.;
191 
192  locPosCalc[0] -= actBox->GetDX();
193  locPosCalc[1] -= actBox->GetDY();
194 
195  Double_t globPos[3];
196 
197  curNode->LocalToMaster(locPosCalc, globPos);
198 
199  LOG(debug) << "HIT ON " << detId << " POSITION: " << locPosCalc[0] << " / " << locPosCalc[1];
200  LOG(debug) << "GLOB HIT " << detId << " POSITION: " << globPos[0] << " / " << globPos[1] << " / "
201  << globPos[2];
202 
203  TVector3 pos(globPos[0], globPos[1], globPos[2]);
204  TVector3 posErr(fPitchX / TMath::Sqrt(12.), fPitchY / TMath::Sqrt(12.), actBox->GetDZ());
205 
206  new ((*hits)[fNHits]) PixelHit(detId, currentDigi->GetIndex(), pos, posErr);
207 
208  fNHits++;
209  }
210 
211  fTNofHits += fNHits;
212  // return fHits;
213 }
Int_t GetIndex()
Definition: PixelDigi.h:44
TGeoManager * GetGeometry()
Definition: FairGeoParSet.h:71
Double_t GetYPitch() const
Definition: PixelDigiPar.h:38
Int_t GetRow()
Definition: PixelDigi.h:49
Int_t GetFERows() const
Definition: PixelDigiPar.h:43
void Init(PixelDigiPar *digipar, FairGeoParSet *geopar)
Digitization Parameter Class for Pixel detector.
Definition: PixelDigiPar.h:24
void Exec(TClonesArray *digis, TClonesArray *hits)
Int_t GetFECols() const
Definition: PixelDigiPar.h:42
virtual void printParams()
Int_t GetDetectorID()
Definition: PixelDigi.h:45
Int_t GetCol()
Definition: PixelDigi.h:48
Int_t GetFeID()
Definition: PixelDigi.h:46
Double_t GetXPitch() const
Definition: PixelDigiPar.h:37
Int_t GetMaxFEperCol() const
Definition: PixelDigiPar.h:44