BmnRoot
Loading...
Searching...
No Matches
BmnSiProfHitMaker.cxx
Go to the documentation of this file.
1#include "BmnSiProfHitMaker.h"
2
3#include "BmnEventHeader.h"
4#include "FairRunAna.h"
5#include "TSystem.h"
6
7#include <TChain.h>
8#include <TRandom.h>
9#include <TStopwatch.h>
10#include <zlib.h>
11
12static Double_t workTime = 0.0;
13
15 : fHitMatching(kTRUE)
16{
17
18 fInputPointsBranchName = "SiProfPoint";
19 fInputDigitsBranchName = "BmnSiProfDigit";
20 fInputDigitMatchesBranchName = "BmnSiProfDigitMatch";
21
22 fOutputHitsBranchName = "BmnSiProfHit";
23 fOutputHitMatchesBranchName = "BmnSiProfHitMatch";
24
25 fVerbose = 1;
26 fField = NULL;
27
28 fCurrentConfig = BmnSiProfConfiguration::None;
29 StationSet = nullptr;
30 TransfSet = nullptr;
31}
32
33BmnSiProfHitMaker::BmnSiProfHitMaker(Int_t run_period, Int_t run_number, Bool_t isExp, TString alignFile)
34 : fHitMatching(kTRUE)
35{
36
37 fInputPointsBranchName = "SiProfPoint";
38 fInputDigitsBranchName = (!isExp) ? "BmnSiProfDigit" : "SiProf";
39 fIsExp = isExp;
40
41 fInputDigitMatchesBranchName = "BmnSiProfDigitMatch";
42
43 fOutputHitsBranchName = "BmnSiProfHit";
44 fOutputHitMatchesBranchName = "BmnSiProfHitMatch";
45
46 fBmnEvQualityBranchName = "BmnEventQuality";
47
48 fVerbose = 1;
49 fField = NULL;
50
51 fCurrentConfig = BmnSiProfConfiguration::None;
52 StationSet = nullptr;
53 TransfSet = nullptr;
54
55 switch (run_period) {
56 case 8: // BM@N RUN-8
57 fCurrentConfig = BmnSiProfConfiguration::Run8;
58 break;
59 case 9: // BM@N RUN-9
60 fCurrentConfig = BmnSiProfConfiguration::Run9;
61 break;
62
63 default:
64 fCurrentConfig = BmnSiProfConfiguration::None;
65 }
66}
67
69{
70 if (StationSet) {
71 delete StationSet;
72 }
73
74 if (TransfSet) {
75 delete TransfSet;
76 }
77}
78
79void BmnSiProfHitMaker::LoadDetectorConfiguration()
80{
81 TString gPathSiProfConfig = gSystem->Getenv("VMCWORKDIR");
82 gPathSiProfConfig += "/parameters/profilometer/XMLConfigs/";
83
84 switch (fCurrentConfig) {
86 XMLConfigFile = gPathSiProfConfig + "ProfRun8.xml";
87 if (fVerbose)
88 cout << " Current SiProf Configuration : Run8"
89 << "\n";
90 break;
91
93 XMLConfigFile = gPathSiProfConfig + "ProfRun9.xml";
94 if (fVerbose)
95 cout << " Current SiProf Configuration : Run9"
96 << "\n";
97 break;
98
99 default:;
100 }
101
102 if (!gSystem->AccessPathName(XMLConfigFile)) { // returns FALSE if the file exists
103 if (fVerbose)
104 cout << " XMLConfigFile : " << XMLConfigFile << "\n";
105 StationSet = new BmnSiProfStationSet(XMLConfigFile);
106 TransfSet = new BmnSiProfTransform();
107 TransfSet->LoadFromXMLFile(XMLConfigFile);
108 } else {
109 Fatal("BmnSiProfHitMaker::LoadDetectorConfiguration()", " !!! Current configuration is not set !!! ");
110 }
111}
112
114{
115 if (fVerbose > 1)
116 cout << "=================== BmnSiProfHitMaker::Init() started ====================" << endl;
117
118 FairRootManager* ioman = FairRootManager::Instance();
119
120 fBmnSiProfDigitsArray = (TClonesArray*)ioman->GetObject(fInputDigitsBranchName);
121 if (!fBmnSiProfDigitsArray) {
122 cout << "BmnSiProfHitMaker::Init(): branch " << fInputDigitsBranchName << " not found! Task will be deactivated"
123 << endl;
124 SetActive(kFALSE);
125 return kERROR;
126 }
127
128 fBmnSiProfDigitMatchesArray = (TClonesArray*)ioman->GetObject(fInputDigitMatchesBranchName);
129
130 if (fVerbose > 1) {
131 if (fBmnSiProfDigitMatchesArray)
132 cout << " Strip matching information exists!\n";
133 else
134 cout << " Strip matching information doesn`t exist!\n";
135 }
136
137 fBmnSiProfHitsArray = new TClonesArray(fOutputHitsBranchName);
138 ioman->Register(fOutputHitsBranchName, "SiProf", fBmnSiProfHitsArray, kTRUE);
139
140 if (fHitMatching && fBmnSiProfDigitMatchesArray) {
141 fBmnSiProfHitMatchesArray = new TClonesArray("BmnMatch");
142 ioman->Register(fOutputHitMatchesBranchName, "SiProf", fBmnSiProfHitMatchesArray, kTRUE);
143 } else {
144 fBmnSiProfHitMatchesArray = 0;
145 }
146
147 // Create SiProf detector ------------------------------------------------------
148 LoadDetectorConfiguration();
149
150 fField = FairRunAna::Instance()->GetField();
151 if (!fField)
152 Fatal("Init", "No Magnetic Field found");
153
154 //--------------------------------------------------------------------------
155
156 fBmnEvQuality = (TClonesArray*)ioman->GetObject(fBmnEvQualityBranchName);
157
158 if (fVerbose > 1)
159 cout << "=================== BmnSiProfHitMaker::Init() finished ===================" << endl;
160
161 return kSUCCESS;
162}
163
164void BmnSiProfHitMaker::Exec(Option_t* opt)
165{
166 TStopwatch sw;
167 sw.Start();
168
169 if (!IsActive())
170 return;
171
172 // Event separation by triggers ...
173 if (fIsExp && fBmnEvQuality) {
174 BmnEventQuality* evQual = (BmnEventQuality*)fBmnEvQuality->UncheckedAt(0);
175 if (!evQual->GetIsGoodEvent())
176 return;
177 }
178 fBmnSiProfHitsArray->Delete();
179
180 if (fHitMatching && fBmnSiProfHitMatchesArray) {
181 fBmnSiProfHitMatchesArray->Delete();
182 }
183
184 if (fVerbose > 1)
185 cout << "=================== BmnSiProfHitMaker::Exec() started ====================" << endl;
186
187 fField = FairRunAna::Instance()->GetField();
188
189 if (fVerbose > 1)
190 cout << " BmnSiProfHitMaker::Exec(), Number of BmnSiProfDigits = " << fBmnSiProfDigitsArray->GetEntriesFast()
191 << "\n";
192
194
195 if (fVerbose > 1)
196 cout << "=================== BmnSiProfHitMaker::Exec() finished ===================" << endl;
197
198 sw.Stop();
199 workTime += sw.RealTime();
200}
201
203{
204 // FairMCPoint* MCPoint;
205 BmnSiProfDigit* digit;
206 BmnMatch* strip_match;
207
208 BmnSiProfStation* station;
209 BmnSiProfModule* module;
210
211 // Loading digits ---------------------------------------------------------------
212 Int_t AddedDigits = 0;
213 Int_t AddedStripDigitMatches = 0;
214
215 for (Int_t idigit = 0; idigit < fBmnSiProfDigitsArray->GetEntriesFast(); idigit++) {
216 digit = (BmnSiProfDigit*)fBmnSiProfDigitsArray->At(idigit);
217 if (!digit->IsGoodDigit())
218 continue;
219 station = StationSet->GetStation(digit->GetStation());
220 module = station->GetModule(digit->GetModule());
221
222 if (module->SetStripSignalInLayer(digit->GetStripLayer(), digit->GetStripNumber(), digit->GetStripSignal()))
223 AddedDigits++;
224
225 if (fBmnSiProfDigitMatchesArray) {
226 strip_match = (BmnMatch*)fBmnSiProfDigitMatchesArray->At(idigit);
227 if (module->SetStripMatchInLayer(digit->GetStripLayer(), digit->GetStripNumber(), *strip_match))
228 AddedStripDigitMatches++;
229 }
230 }
231
232 if (fVerbose > 1)
233 cout << " Processed strip digits : " << AddedDigits << "\n";
234 if (fVerbose > 1 && fBmnSiProfDigitMatchesArray)
235 cout << " Added strip digit matches : " << AddedStripDigitMatches << "\n";
236 //------------------------------------------------------------------------------
237
238 // Processing digits
239 StationSet->ProcessPointsInDetector();
240
241 Int_t NCalculatedPoints = StationSet->CountNProcessedPointsInDetector();
242 if (fVerbose > 1)
243 cout << " Calculated points : " << NCalculatedPoints << "\n";
244 if (fVerbose == 1)
245 cout << "BmnSiProfHitMaker: " << NCalculatedPoints << " hits\n";
246
247 Int_t clear_matched_points_cnt = 0; // points with the only one match-index
248
249 for (Int_t iStation = 0; iStation < StationSet->GetNStations(); ++iStation) {
250 station = StationSet->GetStation(iStation);
251
252 for (Int_t iModule = 0; iModule < station->GetNModules(); ++iModule) {
253 module = station->GetModule(iModule);
254
255 Int_t NIntersectionPointsInModule = module->GetNIntersectionPoints();
256
257 for (Int_t iPoint = 0; iPoint < NIntersectionPointsInModule; ++iPoint) {
258
259 Double_t sigL = module->GetIntersectionPoint_LowerLayerSripTotalSignal(iPoint);
260 Double_t sigU = module->GetIntersectionPoint_UpperLayerSripTotalSignal(iPoint);
261
262 Double_t x = module->GetIntersectionPointX(iPoint);
263 Double_t y = module->GetIntersectionPointY(iPoint);
264 Double_t z = module->GetZPositionRegistered();
265
266 Double_t x_err = module->GetIntersectionPointXError(iPoint);
267 Double_t y_err = module->GetIntersectionPointYError(iPoint);
268 Double_t z_err = 0.0;
269
270 // Transform hit coordinates from local coordinate system of GEM-planes to global
271 if (TransfSet) {
272 Plane3D::Point glob_point = TransfSet->ApplyTransforms(Plane3D::Point(-x, y, z), iStation, iModule);
273 x = -glob_point.X();
274 y = glob_point.Y();
275 z = glob_point.Z();
276 }
277
278 Int_t RefMCIndex = 0;
279
280 // hit matching (define RefMCIndex)) ----------------------------
281 BmnMatch match = module->GetIntersectionPointMatch(iPoint);
282
283 Int_t most_probably_index = -1;
284 Double_t max_weight = 0;
285
286 Int_t n_links = match.GetNofLinks();
287 if (n_links == 1)
288 clear_matched_points_cnt++;
289 for (Int_t ilink = 0; ilink < n_links; ilink++) {
290 Int_t index = match.GetLink(ilink).GetIndex();
291 Double_t weight = match.GetLink(ilink).GetWeight();
292 if (weight > max_weight) {
293 max_weight = weight;
294 most_probably_index = index;
295 }
296 }
297
298 RefMCIndex = most_probably_index;
299 //--------------------------------------------------------------
300
301 // Add hit ------------------------------------------------------
302 x *= -1; // invert to global X //Temporary switched off
303
304 new ((*fBmnSiProfHitsArray)[fBmnSiProfHitsArray->GetEntriesFast()])
305 BmnSiProfHit(0, TVector3(x, y, z), TVector3(x_err, y_err, z_err), RefMCIndex);
306
307 BmnSiProfHit* hit = (BmnSiProfHit*)fBmnSiProfHitsArray->At(fBmnSiProfHitsArray->GetEntriesFast() - 1);
308 hit->SetStation(iStation);
309 hit->SetModule(iModule);
310 hit->SetIndex(fBmnSiProfHitsArray->GetEntriesFast() - 1);
312 module->GetIntersectionPoint_LowerLayerClusterSize(iPoint)); // cluster size (lower layer |||)
313 hit->SetClusterSizeInUpperLayer(module->GetIntersectionPoint_UpperLayerClusterSize(
314 iPoint)); // cluster size (upper layer ///or\\\‍)
316 module->GetIntersectionPoint_LowerLayerSripPosition(iPoint)); // strip position (lower layer |||)
317 hit->SetStripPositionInUpperLayer(module->GetIntersectionPoint_UpperLayerSripPosition(
318 iPoint)); // strip position (upper layer ///or\\\‍)
321
322 if (fVerbose) {
323 cout << " glob(x:y:z) = ( " << x << " : " << y << " : " << z << "\n";
324 cout << " hit(x:y:z) = ( " << hit->GetX() << " : " << hit->GetY() << " : " << hit->GetZ() << "\n";
325 cout << "\n";
326 }
327 //--------------------------------------------------------------
328
329 // hit matching -------------------------------------------------
330 if (fHitMatching && fBmnSiProfHitMatchesArray) {
331 new ((*fBmnSiProfHitMatchesArray)[fBmnSiProfHitMatchesArray->GetEntriesFast()])
332 BmnMatch(module->GetIntersectionPointMatch(iPoint));
333 }
334 //--------------------------------------------------------------
335 }
336 }
337 }
338 if (fVerbose > 1)
339 cout << " N clear matches with MC-points = " << clear_matched_points_cnt << "\n";
340 //------------------------------------------------------------------------------
341 StationSet->Reset();
342}
343
344void BmnSiProfHitMaker::OnlineWrite(const std::unique_ptr<TTree>& resultTree)
345{
346 if (!IsActive())
347 return;
348
349 resultTree->Branch(fOutputHitsBranchName, &fBmnSiProfHitsArray);
350 resultTree->Fill();
351}
352
354{
355 if (StationSet) {
356 delete StationSet;
357 StationSet = nullptr;
358 }
359
360 if (TransfSet) {
361 delete TransfSet;
362 TransfSet = nullptr;
363 }
364
365 printf("Work time of BmnSiProfHitMaker: %4.2f sec.\n", workTime);
366}
Bool_t GetIsGoodEvent()
void SetModule(Int_t mod)
Definition BmnHit.h:73
void SetIndex(Int_t id)
Definition BmnHit.h:57
void SetStation(Short_t st)
Definition BmnHit.h:69
virtual InitStatus Init()
virtual void Exec(Option_t *opt)
virtual void OnlineWrite(const std::unique_ptr< TTree > &resultTree)
Write task resul to tree.
void SetClusterSizeInLowerLayer(Int_t csize)
void SetStripTotalSignalInLowerLayer(Double_t sig)
void SetStripPositionInUpperLayer(Double_t spos)
void SetStripPositionInLowerLayer(Double_t spos)
void SetStripTotalSignalInUpperLayer(Double_t sig)
void SetClusterSizeInUpperLayer(Int_t csize)
BmnSiProfStation * GetStation(Int_t station_num)
Plane3D::Point ApplyTransforms(Plane3D::Point point, Int_t station, Int_t module)
Bool_t LoadFromXMLFile(TString xml_config_file)
Bool_t IsGoodDigit()
Int_t GetStripNumber()
Int_t GetStripLayer()
Int_t GetStation()
Double_t GetStripSignal()