CbmRoot
Loading...
Searching...
No Matches
CbmRecoStsPixel.cxx
Go to the documentation of this file.
1/* Copyright (C) 2021 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Sergey Gorbunov [committer] */
4
9
10#include "CbmRecoStsPixel.h"
11
12#include "CbmAddress.h"
13#include "CbmDigiManager.h"
14#include "CbmEvent.h"
15#include "CbmMCDataArray.h"
16#include "CbmMCDataManager.h"
17#include "CbmStsDigi.h"
18#include "CbmStsModule.h"
19#include "CbmStsParSetModule.h"
20#include "CbmStsParSetSensor.h"
22#include "CbmStsParSim.h"
23#include "CbmStsPoint.h"
24#include "CbmStsRecoModule.h"
25#include "CbmStsSetup.h"
26#include "CbmStsStation.h"
27
28#include <FairField.h>
29#include <FairRootManager.h>
30#include <FairRun.h>
31#include <FairRuntimeDb.h>
32
33#include <TClonesArray.h>
34#include <TGeoBBox.h>
35#include <TGeoPhysicalNode.h>
36#include <TRandom.h>
37
38#include <iomanip>
39
40using std::fixed;
41using std::left;
42using std::right;
43using std::setprecision;
44using std::setw;
45using std::stringstream;
46using std::vector;
47
48
50
51
52 // ----- Constructor ---------------------------------------------------
54 : FairTask("RecoStsPixel", 1)
55 , fMode(mode)
56{
57}
58// -------------------------------------------------------------------------
59
60
61// ----- Destructor ----------------------------------------------------
63// -------------------------------------------------------------------------
64
65
66// ----- End-of-run action ---------------------------------------------
68// -------------------------------------------------------------------------
69
70
71// ----- Initialisation ------------------------------------------------
73{
74
75 // --- Something for the screen
76
77 LOG(info) << "==========================================================";
78 LOG(info) << GetName() << ": Initialising ";
79
80 // --- Check IO-Manager
81 FairRootManager* ioman = FairRootManager::Instance();
82 assert(ioman);
83
84 // --- In event mode: get input array (CbmEvent)
86 LOG(info) << GetName() << ": Using event-by-event mode";
87 fEvents = dynamic_cast<TClonesArray*>(ioman->GetObject("CbmEvent"));
88 if (nullptr == fEvents) {
89 LOG(warn) << GetName() << ": Event mode selected but no event array found!";
90 return kFATAL;
91 } //? Event branch not present
92 } //? Event mode
93 else {
94 LOG(info) << GetName() << ": Using time-based mode";
95 }
96
97
98 // --- Digi Manager
100 fDigiManager->Init();
101
102 // --- Check input array (StsDigis)
103 if (!fDigiManager->IsPresent(ECbmModuleId::kSts)) {
104 LOG(fatal) << GetName() << ": No StsDigi branch in input!";
105 return kERROR;
106 }
107
108 if (!fDigiManager->IsMatchPresent(ECbmModuleId::kSts)) {
109 LOG(error) << GetName() << " sts digi matches are not present";
110 return kERROR;
111 }
112
113 CbmMCDataManager* mcManager = (CbmMCDataManager*) ioman->GetObject("MCDataManager");
114 if (!mcManager) {
115 LOG(error) << GetName() << ": No CbmMCDataManager!";
116 return kERROR;
117 }
118
119 fStsPoints = mcManager->InitBranch("StsPoint");
120
121 if (!fStsPoints) {
122 LOG(fatal) << GetName() << ": No StsPoint branch in input!";
123 return kERROR;
124 }
125
126
127 // --- Register output array
128 fClusters = new TClonesArray("CbmStsCluster", 1);
129 ioman->Register("StsCluster", "Clusters in STS", fClusters, IsOutputBranchPersistent("StsCluster"));
130
131 // --- Register output array
132 fHits = new TClonesArray("CbmStsHit", 1);
133 ioman->Register("StsHit", "Hits in STS", fHits, IsOutputBranchPersistent("StsHit"));
134
135 // --- Simulation settings
136 assert(fParSim);
137 LOG(info) << GetName() << ": Sim settings " << fParSim->ToString();
138
139 // --- Module parameters
140 assert(fParSetModule);
141 LOG(info) << GetName() << ": Module parameters " << fParSetModule->ToString();
142
143 // --- Sensor parameters
144 assert(fParSetSensor);
145 LOG(info) << GetName() << ": Sensor parameters " << fParSetModule->ToString();
146
147 // --- Sensor conditions
148 assert(fParSetCond);
149 //assert(fParSetCond->IsSet());
150 LOG(info) << GetName() << ": Sensor conditions " << fParSetCond->ToString();
151
152 // --- Initialise STS setup
154 if (!fSetup->IsInit()) {
155 fSetup->Init(nullptr);
156 }
157 if (!fSetup->IsModuleParsInit()) {
158 fSetup->SetModuleParameters(fParSetModule);
159 }
160 if (!fSetup->IsSensorParsInit()) {
161 fSetup->SetSensorParameters(fParSetSensor);
162 }
163 if (!fSetup->IsSensorCondInit()) {
164 fSetup->SetSensorConditions(fParSetCond);
165 }
166
167 // make sure that the STS digis were produced by the experimental Pixel digitiser
168 if (strcmp(fParSetSensor->getDescription(), "Experimental STS Pixels")) {
169 LOG(error) << GetName() << " STS digis must be produced by the CbmStsDigitizePixel digitizer";
170 return kERROR;
171 }
172
173 LOG(info) << GetName() << ": Initialisation successful.";
174 LOG(info) << "==========================================================";
175
176 return kSUCCESS;
177}
178// -------------------------------------------------------------------------
179
180
181// ----- Task execution ------------------------------------------------
183{
184
185 // --- Clear hit output array
186 fHits->Delete();
187
188 // --- Reset cluster output array
189 fClusters->Delete();
190
192 // --- Time-slice mode: process entire array
193 ProcessData(nullptr);
194 }
195 else {
196 // --- Event mode: loop over events
197 assert(fEvents);
198 for (Int_t iEvent = 0; iEvent < fEvents->GetEntriesFast(); iEvent++) {
199 CbmEvent* event = dynamic_cast<CbmEvent*>(fEvents->At(iEvent));
200 assert(event);
201 ProcessData(event);
202 } //# events
203 } //? event mode
204}
205// -------------------------------------------------------------------------
206
207
208// ----- Process one time slice or event -------------------------------
210{
211 if (!fDigiManager->IsMatchPresent(ECbmModuleId::kSts)) {
212 LOG(error) << GetName() << ": No StsDigi branch in input!";
213 return;
214 }
215
216 Int_t nDigis = (event ? event->GetNofData(ECbmDataType::kStsDigi) : fDigiManager->GetNofDigis(ECbmModuleId::kSts));
217 int nHits = 0;
218
219 for (Int_t iDigi = 0; iDigi < nDigis; iDigi++) {
220
221 UInt_t digiIndex = (event ? event->GetIndex(ECbmDataType::kStsDigi, iDigi) : iDigi);
222 const CbmStsDigi* digi = fDigiManager->Get<const CbmStsDigi>(digiIndex);
223 assert(digi);
224
225 // Check system ID. There are pulser digis in which will be ignored here.
226 Int_t systemId = CbmAddress::GetSystemId(digi->GetAddress());
227 if (systemId != ToIntegralType(ECbmModuleId::kSts)) {
228 continue;
229 }
230
231 const CbmMatch* match = fDigiManager->GetMatch(ECbmModuleId::kSts, digiIndex);
232 assert(match);
233
234 // make sure that the digi was produced by CbmStsDigitizePixel task
235 if (digi->GetChannel() != 0 || match->GetNofLinks() != 1) {
236 LOG(error) << GetName() << " sts digis were not produced by CbmStsDigitizePixel task";
237 break;
238 }
239
240 // the digi was produced by one MC point
241
242 const CbmLink& link = match->GetLink(0);
243 CbmStsPoint* pt = dynamic_cast<CbmStsPoint*>(fStsPoints->Get(link.GetFile(), link.GetEntry(), link.GetIndex()));
244
245 // create one cluster and one hit per digi
246
247 UInt_t iCluster = fClusters->GetEntriesFast();
248 CbmStsCluster* cluster = new ((*fClusters)[iCluster]) CbmStsCluster;
249 cluster->AddDigi(iDigi);
250 if (event) event->AddData(ECbmDataType::kStsCluster, iCluster);
251
252 UInt_t iHit = fHits->GetEntriesFast();
253 CbmStsHit* hit = new ((*fHits)[iHit]) CbmStsHit;
254
255 hit->SetAddress(digi->GetAddress());
256 hit->SetFrontClusterId(iCluster);
257 hit->SetBackClusterId(iCluster);
258
260 if (ista < 0 || ista >= CbmStsSetup::Instance()->GetNofStations()) {
261 LOG(error) << "wrong Sts station number " << ista;
262 break;
263 }
264
266 double staZ = station->GetZ();
267
268 if (fabs(pt->GetZ() - staZ) > 1.) {
269 LOG(error) << "Sts point Z " << pt->GetZ() << " is far from the station Z " << staZ;
270 break;
271 }
272
273 hit->SetX(0.5 * (pt->GetXIn() + pt->GetXOut()));
274 hit->SetY(0.5 * (pt->GetYIn() + pt->GetYOut()));
275 hit->SetZ(0.5 * (pt->GetZIn() + pt->GetZOut()));
276
277 Double_t resX = station->GetSensorPitch(0) / sqrt(12.);
278 Double_t resY = station->GetSensorPitch(1) / sqrt(12.);
279
280 assert(resX > 1.e-5);
281 assert(resY > 1.e-5);
282
283 auto gaus = []() {
284 double x = 5;
285 while (fabs(x) > 3.5) {
286 x = gRandom->Gaus(0., 1.);
287 }
288 return x;
289 };
290
291 hit->SetX(hit->GetX() + resX * gaus());
292 hit->SetY(hit->GetY() + resY * gaus());
293 hit->SetDx(resX);
294 hit->SetDy(resY);
295 hit->SetDxy(0.);
296 hit->SetDu(hit->GetDx());
297 hit->SetDv(hit->GetDy());
298
299 hit->SetTime(digi->GetTime());
300 const CbmStsParModule& module = fParSetModule->GetParModule(digi->GetAddress());
301 const CbmStsParAsic& asic = module.GetParAsic(0);
302 hit->SetTimeError(asic.GetTimeResol());
303
304 if (event) event->AddData(ECbmDataType::kStsHit, iHit);
305 nHits++;
306 } // digis
307
308 LOG(info) << GetName() << " Processed " << nDigis << " digis, created " << nHits << " hits";
309}
310// -------------------------------------------------------------------------
311
312
313// ----- Connect parameter container -----------------------------------
315{
316 FairRuntimeDb* db = FairRun::Instance()->GetRuntimeDb();
317 fParSim = dynamic_cast<CbmStsParSim*>(db->getContainer("CbmStsParSim"));
318 fParSetModule = dynamic_cast<CbmStsParSetModule*>(db->getContainer("CbmStsParSetModule"));
319 fParSetSensor = dynamic_cast<CbmStsParSetSensor*>(db->getContainer("CbmStsParSetSensor"));
320 fParSetCond = dynamic_cast<CbmStsParSetSensorCond*>(db->getContainer("CbmStsParSetSensorCond"));
321}
322// -------------------------------------------------------------------------
ECbmRecoMode
Reconstruct the full time slice or event-by-event.
Definition CbmDefs.h:202
XPU_D constexpr auto ToIntegralType(T enumerator) -> typename std::underlying_type< T >::type
Converts an element of enum class to its underlying integral type.
Definition CbmDefs.h:33
@ kSts
Silicon Tracking System.
Definition CbmDefs.h:48
ClassImp(CbmRecoStsPixel) CbmRecoStsPixel
friend fvec sqrt(const fvec &a)
int Int_t
static int32_t GetSystemId(uint32_t address)
Definition CbmAddress.h:47
void AddDigi(int32_t index)
Add digi to cluster.
Definition CbmCluster.h:51
static CbmDigiManager * Instance()
Static instance.
Class characterising one event by a collection of links (indices) to data objects,...
Definition CbmEvent.h:34
void SetTimeError(double error)
Definition CbmHit.h:94
void SetAddress(int32_t address)
Definition CbmHit.h:86
int32_t GetAddress() const
Definition CbmHit.h:77
void SetZ(double z)
Definition CbmHit.h:83
void SetTime(double time)
Definition CbmHit.h:88
Task class creating and managing CbmMCDataArray objects.
CbmMCDataObject * GetObject(const char *name)
CbmMCDataArray * InitBranch(const char *name)
const CbmLink & GetLink(int32_t i) const
Definition CbmMatch.h:39
int32_t GetNofLinks() const
Definition CbmMatch.h:42
void SetDxy(double dxy)
Definition CbmPixelHit.h:96
void SetDx(double dx)
Definition CbmPixelHit.h:94
double GetDy() const
Definition CbmPixelHit.h:76
double GetDx() const
Definition CbmPixelHit.h:75
void SetX(double x)
Definition CbmPixelHit.h:92
double GetY() const
Definition CbmPixelHit.h:74
void SetY(double y)
Definition CbmPixelHit.h:93
double GetX() const
Definition CbmPixelHit.h:73
void SetDy(double dy)
Definition CbmPixelHit.h:95
Task class for local reconstruction in the STS Pixel detector.
void Exec(Option_t *opt)
Task execution.
CbmStsParSetSensorCond * fParSetCond
Sensor conditions.
TClonesArray * fEvents
CbmStsParSetModule * fParSetModule
Module parameters.
TClonesArray * fClusters
void SetParContainers()
Define the needed parameter containers.
CbmMCDataArray * fStsPoints
Interface to digi branch.
InitStatus Init()
Initialisation.
TClonesArray * fHits
Output cluster array.
CbmDigiManager * fDigiManager
Input array of events.
CbmStsParSetSensor * fParSetSensor
Sensor parameters.
CbmRecoStsPixel(ECbmRecoMode mode=ECbmRecoMode::Timeslice)
Constructor.
~CbmRecoStsPixel()
Destructor.
ECbmRecoMode fMode
Time-slice or event.
CbmStsParSim * fParSim
Instance of STS setup.
void ProcessData(CbmEvent *event=nullptr)
Process one time slice or event.
void Finish()
End-of-run action.
CbmStsSetup * fSetup
Output hit array.
Data class for STS clusters.
Data class for a single-channel message in the STS.
Definition CbmStsDigi.h:40
XPU_D uint16_t GetChannel() const
Channel number in module @value Channel number.
Definition CbmStsDigi.h:93
XPU_D int32_t GetAddress() const
Definition CbmStsDigi.h:74
data class for a reconstructed 3-d hit in the STS
Definition CbmStsHit.h:35
void SetDv(Double_t dv)
Error of coordinate across front-side strips @value Coordinate error [cm].
Definition CbmStsHit.h:99
void SetFrontClusterId(int32_t index)
Set the index of the frontside cluster To keep track of the input during matching.
Definition CbmStsHit.h:111
void SetBackClusterId(int32_t index)
Set the index of the backside cluster To keep track of the input during matching.
Definition CbmStsHit.h:117
void SetDu(Double_t du)
Error of coordinate across front-side strips @value Coordinate error [cm].
Definition CbmStsHit.h:82
Parameters of the STS readout ASIC.
double GetTimeResol() const
Time resolution.
Parameters for one STS module.
Parameters container for CbmStsParModule.
Parameters container for CbmStsParSensorCond.
Parameters container for CbmStsParSensor.
Settings for STS simulation (digitizer)
double GetZOut() const
Definition CbmStsPoint.h:78
double GetXOut() const
Definition CbmStsPoint.h:76
double GetYIn() const
Definition CbmStsPoint.h:74
double GetXIn() const
Definition CbmStsPoint.h:73
double GetZIn() const
Definition CbmStsPoint.h:75
double GetYOut() const
Definition CbmStsPoint.h:77
static CbmStsSetup * Instance()
CbmStsStation * GetStation(Int_t stationId) const
Int_t GetStationNumber(Int_t address)
Class representing a station of the StsSystem.
Double_t GetSensorPitch(Int_t iSide) const
Double_t GetZ() const