596a855f |
1 | /************************************************************************** |
2 | * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * |
3 | * * |
4 | * Author: The ALICE Off-line Project. * |
5 | * Contributors are mentioned in the code where appropriate. * |
6 | * * |
7 | * Permission to use, copy, modify and distribute this software and its * |
8 | * documentation strictly for non-commercial purposes is hereby granted * |
9 | * without fee, provided that the above copyright notice appears in all * |
10 | * copies and that both the copyright notice and this permission notice * |
11 | * appear in the supporting documentation. The authors make no claims * |
12 | * about the suitability of this software for any purpose. It is * |
13 | * provided "as is" without express or implied warranty. * |
14 | **************************************************************************/ |
15 | |
16 | /* $Id$ */ |
17 | |
18 | /////////////////////////////////////////////////////////////////////////////// |
19 | // // |
20 | // class for running the reconstruction // |
21 | // // |
22 | // Clusters and tracks are created for all detectors and all events by // |
23 | // typing: // |
24 | // // |
25 | // AliReconstruction rec; // |
26 | // rec.Run(); // |
27 | // // |
28 | // The Run method returns kTRUE in case of successful execution. // |
29 | // The name of the galice file can be changed from the default // |
e583c30d |
30 | // "galice.root" by passing it as argument to the AliReconstruction // |
31 | // constructor or by // |
596a855f |
32 | // // |
33 | // rec.SetGAliceFile("..."); // |
34 | // // |
59697224 |
35 | // The local reconstruction can be switched on or off for individual // |
36 | // detectors by // |
596a855f |
37 | // // |
59697224 |
38 | // rec.SetRunLocalReconstruction("..."); // |
596a855f |
39 | // // |
40 | // The argument is a (case sensitive) string with the names of the // |
41 | // detectors separated by a space. The special string "ALL" selects all // |
42 | // available detectors. This is the default. // |
43 | // // |
44 | // The tracking in ITS, TPC and TRD and the creation of ESD tracks can be // |
45 | // switched off by // |
46 | // // |
47 | // rec.SetRunTracking(kFALSE); // |
48 | // // |
49 | // The filling of additional ESD information can be steered by // |
50 | // // |
51 | // rec.SetFillESD("..."); // |
52 | // // |
53 | // Again, the string specifies the list of detectors. The default is "ALL". // |
54 | // // |
55 | // The reconstruction requires digits as input. For the creation of digits // |
56 | // have a look at the class AliSimulation. // |
57 | // // |
24f7a148 |
58 | // For debug purposes the method SetCheckPointLevel can be used. If the // |
59 | // argument is greater than 0, files with ESD events will be written after // |
60 | // selected steps of the reconstruction for each event: // |
61 | // level 1: after tracking and after filling of ESD (final) // |
62 | // level 2: in addition after each tracking step // |
63 | // level 3: in addition after the filling of ESD for each detector // |
64 | // If a final check point file exists for an event, this event will be // |
65 | // skipped in the reconstruction. The tracking and the filling of ESD for // |
66 | // a detector will be skipped as well, if the corresponding check point // |
67 | // file exists. The ESD event will then be loaded from the file instead. // |
68 | // // |
596a855f |
69 | /////////////////////////////////////////////////////////////////////////////// |
70 | |
024a7e64 |
71 | #include <TArrayF.h> |
72 | #include <TFile.h> |
73 | #include <TSystem.h> |
74 | #include <TROOT.h> |
75 | #include <TPluginManager.h> |
fd46e2d2 |
76 | #include <TStopwatch.h> |
596a855f |
77 | |
78 | #include "AliReconstruction.h" |
79 | #include "AliRunLoader.h" |
80 | #include "AliRun.h" |
596a855f |
81 | #include "AliTracker.h" |
82 | #include "AliESD.h" |
2257f27e |
83 | #include "AliESDVertex.h" |
84 | #include "AliVertexer.h" |
596a855f |
85 | #include "AliHeader.h" |
86 | #include "AliGenEventHeader.h" |
87 | #include "AliESDpid.h" |
a866ac60 |
88 | #include "AliMagF.h" |
596a855f |
89 | |
90 | ClassImp(AliReconstruction) |
91 | |
92 | |
c757bafd |
93 | //_____________________________________________________________________________ |
94 | const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "RICH", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "START", "VZERO", "CRT"}; |
95 | |
596a855f |
96 | //_____________________________________________________________________________ |
e583c30d |
97 | AliReconstruction::AliReconstruction(const char* gAliceFilename, |
98 | const char* name, const char* title) : |
99 | TNamed(name, title), |
100 | |
59697224 |
101 | fRunLocalReconstruction("ALL"), |
2257f27e |
102 | fRunVertexFinder(kTRUE), |
e583c30d |
103 | fRunTracking(kTRUE), |
104 | fFillESD("ALL"), |
105 | fGAliceFileName(gAliceFilename), |
106 | fStopOnError(kFALSE), |
24f7a148 |
107 | fCheckPointLevel(0), |
e583c30d |
108 | |
109 | fRunLoader(NULL), |
110 | fITSLoader(NULL), |
2257f27e |
111 | fITSVertexer(NULL), |
e583c30d |
112 | fITSTracker(NULL), |
113 | fTPCLoader(NULL), |
114 | fTPCTracker(NULL), |
115 | fTRDLoader(NULL), |
116 | fTRDTracker(NULL), |
117 | fTOFLoader(NULL), |
118 | fTOFTracker(NULL) |
596a855f |
119 | { |
120 | // create reconstruction object with default parameters |
121 | |
596a855f |
122 | } |
123 | |
124 | //_____________________________________________________________________________ |
125 | AliReconstruction::AliReconstruction(const AliReconstruction& rec) : |
e583c30d |
126 | TNamed(rec), |
127 | |
59697224 |
128 | fRunLocalReconstruction(rec.fRunLocalReconstruction), |
2257f27e |
129 | fRunVertexFinder(rec.fRunVertexFinder), |
e583c30d |
130 | fRunTracking(rec.fRunTracking), |
131 | fFillESD(rec.fFillESD), |
132 | fGAliceFileName(rec.fGAliceFileName), |
133 | fStopOnError(rec.fStopOnError), |
24f7a148 |
134 | fCheckPointLevel(0), |
e583c30d |
135 | |
136 | fRunLoader(NULL), |
137 | fITSLoader(NULL), |
2257f27e |
138 | fITSVertexer(NULL), |
e583c30d |
139 | fITSTracker(NULL), |
140 | fTPCLoader(NULL), |
141 | fTPCTracker(NULL), |
142 | fTRDLoader(NULL), |
143 | fTRDTracker(NULL), |
144 | fTOFLoader(NULL), |
145 | fTOFTracker(NULL) |
596a855f |
146 | { |
147 | // copy constructor |
148 | |
596a855f |
149 | } |
150 | |
151 | //_____________________________________________________________________________ |
152 | AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec) |
153 | { |
154 | // assignment operator |
155 | |
156 | this->~AliReconstruction(); |
157 | new(this) AliReconstruction(rec); |
158 | return *this; |
159 | } |
160 | |
161 | //_____________________________________________________________________________ |
162 | AliReconstruction::~AliReconstruction() |
163 | { |
164 | // clean up |
165 | |
e583c30d |
166 | CleanUp(); |
596a855f |
167 | } |
168 | |
169 | |
170 | //_____________________________________________________________________________ |
171 | void AliReconstruction::SetGAliceFile(const char* fileName) |
172 | { |
173 | // set the name of the galice file |
174 | |
175 | fGAliceFileName = fileName; |
176 | } |
177 | |
178 | |
179 | //_____________________________________________________________________________ |
180 | Bool_t AliReconstruction::Run() |
181 | { |
182 | // run the reconstruction |
183 | |
184 | // open the run loader |
596a855f |
185 | fRunLoader = AliRunLoader::Open(fGAliceFileName.Data()); |
186 | if (!fRunLoader) { |
187 | Error("Run", "no run loader found in file %s", |
188 | fGAliceFileName.Data()); |
e583c30d |
189 | CleanUp(); |
596a855f |
190 | return kFALSE; |
191 | } |
192 | fRunLoader->LoadgAlice(); |
e583c30d |
193 | AliRun* aliRun = fRunLoader->GetAliRun(); |
194 | if (!aliRun) { |
596a855f |
195 | Error("Run", "no gAlice object found in file %s", |
196 | fGAliceFileName.Data()); |
e583c30d |
197 | CleanUp(); |
596a855f |
198 | return kFALSE; |
199 | } |
e583c30d |
200 | gAlice = aliRun; |
596a855f |
201 | |
59697224 |
202 | // load the reconstructor objects |
203 | TPluginManager* pluginManager = gROOT->GetPluginManager(); |
c757bafd |
204 | for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) { |
205 | TString detName = fgkDetectorName[iDet]; |
206 | TString recName = "Ali" + detName + "Reconstructor"; |
207 | if (!gAlice->GetDetector(detName)) continue; |
208 | |
209 | AliReconstructor* reconstructor = NULL; |
210 | // first check if a plugin is defined for the reconstructor |
59697224 |
211 | TPluginHandler* pluginHandler = |
c757bafd |
212 | pluginManager->FindHandler("AliReconstructor", detName); |
213 | // if not, but the reconstructor class is implemented, add a plugin for it |
214 | if (!pluginHandler && gROOT->GetClass(recName.Data())) { |
215 | Info("Run", "defining plugin for %s", recName.Data()); |
216 | pluginManager->AddHandler("AliReconstructor", detName, |
217 | recName, detName, recName + "()"); |
218 | pluginHandler = pluginManager->FindHandler("AliReconstructor", detName); |
219 | } |
220 | if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) { |
221 | reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0); |
222 | } |
223 | // if there is no reconstructor class for the detector use the dummy one |
224 | if (!reconstructor) { |
225 | Info("Run", "using dummy reconstructor for %s", detName.Data()); |
226 | reconstructor = new AliDummyReconstructor(gAlice->GetDetector(detName)); |
227 | } |
59697224 |
228 | if (reconstructor) fReconstructors.Add(reconstructor); |
229 | } |
230 | |
596a855f |
231 | // local reconstruction |
59697224 |
232 | if (!fRunLocalReconstruction.IsNull()) { |
233 | if (!RunLocalReconstruction(fRunLocalReconstruction)) { |
e583c30d |
234 | if (fStopOnError) {CleanUp(); return kFALSE;} |
596a855f |
235 | } |
236 | } |
2257f27e |
237 | if (!fRunVertexFinder && !fRunTracking && fFillESD.IsNull()) return kTRUE; |
238 | |
239 | // get vertexer |
240 | if (fRunVertexFinder && !CreateVertexer()) { |
241 | if (fStopOnError) { |
242 | CleanUp(); |
243 | return kFALSE; |
244 | } |
245 | } |
596a855f |
246 | |
24f7a148 |
247 | // get loaders and trackers |
248 | if (fRunTracking && !CreateTrackers()) { |
249 | if (fStopOnError) { |
250 | CleanUp(); |
251 | return kFALSE; |
252 | } |
596a855f |
253 | } |
24f7a148 |
254 | |
596a855f |
255 | // create the ESD output file |
256 | TFile* file = TFile::Open("AliESDs.root", "RECREATE"); |
257 | if (!file->IsOpen()) { |
258 | Error("Run", "opening AliESDs.root failed"); |
e583c30d |
259 | if (fStopOnError) {CleanUp(file); return kFALSE;} |
596a855f |
260 | } |
261 | |
262 | // loop over events |
263 | for (Int_t iEvent = 0; iEvent < fRunLoader->GetNumberOfEvents(); iEvent++) { |
264 | Info("Run", "processing event %d", iEvent); |
596a855f |
265 | fRunLoader->GetEvent(iEvent); |
24f7a148 |
266 | |
267 | char fileName[256]; |
268 | sprintf(fileName, "ESD_%d.%d_final.root", |
269 | aliRun->GetRunNumber(), aliRun->GetEvNumber()); |
270 | if (!gSystem->AccessPathName(fileName)) continue; |
271 | |
272 | AliESD* esd = new AliESD; |
e583c30d |
273 | esd->SetRunNumber(aliRun->GetRunNumber()); |
274 | esd->SetEventNumber(aliRun->GetEvNumber()); |
275 | esd->SetMagneticField(aliRun->Field()->SolenoidField()); |
596a855f |
276 | |
2257f27e |
277 | // vertex finder |
278 | if (fRunVertexFinder) { |
279 | if (!ReadESD(esd, "vertex")) { |
280 | if (!RunVertexFinder(esd)) { |
281 | if (fStopOnError) {CleanUp(file); return kFALSE;} |
282 | } |
283 | if (fCheckPointLevel > 0) WriteESD(esd, "vertex"); |
284 | } |
285 | } |
286 | |
596a855f |
287 | // barrel tracking |
288 | if (fRunTracking) { |
24f7a148 |
289 | if (!ReadESD(esd, "tracking")) { |
290 | if (!RunTracking(esd)) { |
291 | if (fStopOnError) {CleanUp(file); return kFALSE;} |
292 | } |
293 | if (fCheckPointLevel > 0) WriteESD(esd, "tracking"); |
596a855f |
294 | } |
295 | } |
296 | |
297 | // fill ESD |
298 | if (!fFillESD.IsNull()) { |
299 | if (!FillESD(esd, fFillESD)) { |
e583c30d |
300 | if (fStopOnError) {CleanUp(file); return kFALSE;} |
596a855f |
301 | } |
302 | } |
303 | |
304 | // combined PID |
305 | AliESDpid::MakePID(esd); |
24f7a148 |
306 | if (fCheckPointLevel > 1) WriteESD(esd, "PID"); |
596a855f |
307 | |
308 | // write ESD |
309 | char name[100]; |
310 | sprintf(name, "ESD%d", iEvent); |
311 | file->cd(); |
312 | if (!esd->Write(name)) { |
313 | Error("Run", "writing ESD failed"); |
e583c30d |
314 | if (fStopOnError) {CleanUp(file); return kFALSE;} |
596a855f |
315 | } |
24f7a148 |
316 | file->Flush(); |
317 | |
318 | if (fCheckPointLevel > 0) WriteESD(esd, "final"); |
d9b8978b |
319 | delete esd; |
596a855f |
320 | } |
321 | |
e583c30d |
322 | CleanUp(file); |
596a855f |
323 | |
324 | return kTRUE; |
325 | } |
326 | |
327 | |
328 | //_____________________________________________________________________________ |
59697224 |
329 | Bool_t AliReconstruction::RunLocalReconstruction(const TString& detectors) |
596a855f |
330 | { |
59697224 |
331 | // run the local reconstruction |
596a855f |
332 | |
030b532d |
333 | TStopwatch stopwatch; |
334 | stopwatch.Start(); |
335 | |
596a855f |
336 | TString detStr = detectors; |
c757bafd |
337 | for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) { |
338 | AliReconstructor* reconstructor = |
339 | (AliReconstructor*) fReconstructors[iDet]; |
340 | TString detName = reconstructor->GetDetectorName(); |
341 | if (IsSelected(detName, detStr)) { |
596a855f |
342 | Info("RunReconstruction", "running reconstruction for %s", |
c757bafd |
343 | detName.Data()); |
030b532d |
344 | TStopwatch stopwatchDet; |
345 | stopwatchDet.Start(); |
c757bafd |
346 | reconstructor->Reconstruct(fRunLoader); |
347 | Info("RunReconstruction", "execution time for %s:", detName.Data()); |
030b532d |
348 | stopwatchDet.Print(); |
596a855f |
349 | } |
350 | } |
351 | |
352 | if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) { |
353 | Error("RunReconstruction", "the following detectors were not found: %s", |
354 | detStr.Data()); |
355 | if (fStopOnError) return kFALSE; |
356 | } |
357 | |
030b532d |
358 | Info("RunReconstruction", "execution time:"); |
359 | stopwatch.Print(); |
360 | |
596a855f |
361 | return kTRUE; |
362 | } |
363 | |
364 | //_____________________________________________________________________________ |
2257f27e |
365 | Bool_t AliReconstruction::RunVertexFinder(AliESD*& esd) |
596a855f |
366 | { |
367 | // run the barrel tracking |
368 | |
030b532d |
369 | TStopwatch stopwatch; |
370 | stopwatch.Start(); |
371 | |
2257f27e |
372 | AliESDVertex* vertex = NULL; |
373 | Double_t vtxPos[3] = {0, 0, 0}; |
374 | Double_t vtxErr[3] = {0.07, 0.07, 0.1}; |
375 | TArrayF mcVertex(3); |
376 | fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex); |
377 | for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i]; |
378 | |
379 | if (fITSVertexer) { |
380 | Info("RunVertexFinder", "running the ITS vertex finder"); |
381 | fITSVertexer->SetDebug(1); |
382 | vertex = fITSVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber()); |
383 | if(!vertex){ |
384 | Warning("RunVertexFinder","Vertex not found \n"); |
c710f220 |
385 | vertex = new AliESDVertex(); |
2257f27e |
386 | } |
387 | else { |
388 | vertex->SetTruePos(vtxPos); // store also the vertex from MC |
389 | } |
390 | |
391 | } else { |
392 | Info("RunVertexFinder", "getting the primary vertex from MC"); |
393 | vertex = new AliESDVertex(vtxPos, vtxErr); |
394 | } |
395 | |
396 | if (vertex) { |
397 | vertex->GetXYZ(vtxPos); |
398 | vertex->GetSigmaXYZ(vtxErr); |
399 | } else { |
400 | Warning("RunVertexFinder", "no vertex reconstructed"); |
401 | vertex = new AliESDVertex(vtxPos, vtxErr); |
402 | } |
403 | esd->SetVertex(vertex); |
24f7a148 |
404 | if (fITSTracker) fITSTracker->SetVertex(vtxPos, vtxErr); |
405 | if (fTPCTracker) fTPCTracker->SetVertex(vtxPos, vtxErr); |
406 | if (fTRDTracker) fTRDTracker->SetVertex(vtxPos, vtxErr); |
2257f27e |
407 | delete vertex; |
408 | |
409 | Info("RunVertexFinder", "execution time:"); |
410 | stopwatch.Print(); |
411 | |
412 | return kTRUE; |
413 | } |
414 | |
415 | //_____________________________________________________________________________ |
416 | Bool_t AliReconstruction::RunTracking(AliESD*& esd) |
417 | { |
418 | // run the barrel tracking |
419 | |
420 | TStopwatch stopwatch; |
421 | stopwatch.Start(); |
24f7a148 |
422 | |
423 | if (!fTPCTracker) { |
424 | Error("RunTracking", "no TPC tracker"); |
425 | return kFALSE; |
426 | } |
596a855f |
427 | |
85555ca8 |
428 | // TPC tracking |
429 | Info("RunTracking", "TPC tracking"); |
596a855f |
430 | fTPCLoader->LoadRecPoints("read"); |
431 | TTree* tpcTree = fTPCLoader->TreeR(); |
432 | if (!tpcTree) { |
433 | Error("RunTracking", "Can't get the TPC cluster tree"); |
434 | return kFALSE; |
24f7a148 |
435 | } |
596a855f |
436 | fTPCTracker->LoadClusters(tpcTree); |
437 | if (fTPCTracker->Clusters2Tracks(esd) != 0) { |
85555ca8 |
438 | Error("RunTracking", "TPC Clusters2Tracks failed"); |
596a855f |
439 | return kFALSE; |
440 | } |
24f7a148 |
441 | if (fCheckPointLevel > 1) WriteESD(esd, "TPC.tracking"); |
442 | |
443 | if (!fITSTracker) { |
444 | Warning("RunTracking", "no ITS tracker"); |
445 | } else { |
446 | |
447 | fRunLoader->GetAliRun()->GetDetector("TPC")->FillESD(esd); // preliminary |
448 | AliESDpid::MakePID(esd); // PID for the ITS tracker |
449 | |
450 | // ITS tracking |
451 | Info("RunTracking", "ITS tracking"); |
452 | fITSLoader->LoadRecPoints("read"); |
453 | TTree* itsTree = fITSLoader->TreeR(); |
454 | if (!itsTree) { |
455 | Error("RunTracking", "Can't get the ITS cluster tree"); |
456 | return kFALSE; |
457 | } |
458 | fITSTracker->LoadClusters(itsTree); |
459 | if (fITSTracker->Clusters2Tracks(esd) != 0) { |
460 | Error("RunTracking", "ITS Clusters2Tracks failed"); |
461 | return kFALSE; |
462 | } |
463 | if (fCheckPointLevel > 1) WriteESD(esd, "ITS.tracking"); |
596a855f |
464 | |
24f7a148 |
465 | if (!fTRDTracker) { |
466 | Warning("RunTracking", "no TRD tracker"); |
467 | } else { |
468 | // ITS back propagation |
469 | Info("RunTracking", "ITS back propagation"); |
470 | if (fITSTracker->PropagateBack(esd) != 0) { |
471 | Error("RunTracking", "ITS backward propagation failed"); |
472 | return kFALSE; |
473 | } |
474 | if (fCheckPointLevel > 1) WriteESD(esd, "ITS.back"); |
596a855f |
475 | |
24f7a148 |
476 | // TPC back propagation |
477 | Info("RunTracking", "TPC back propagation"); |
478 | if (fTPCTracker->PropagateBack(esd) != 0) { |
479 | Error("RunTracking", "TPC backward propagation failed"); |
480 | return kFALSE; |
481 | } |
482 | if (fCheckPointLevel > 1) WriteESD(esd, "TPC.back"); |
483 | |
484 | // TRD back propagation |
485 | Info("RunTracking", "TRD back propagation"); |
486 | fTRDLoader->LoadRecPoints("read"); |
487 | TTree* trdTree = fTRDLoader->TreeR(); |
488 | if (!trdTree) { |
489 | Error("RunTracking", "Can't get the TRD cluster tree"); |
490 | return kFALSE; |
491 | } |
492 | fTRDTracker->LoadClusters(trdTree); |
493 | if (fTRDTracker->PropagateBack(esd) != 0) { |
494 | Error("RunTracking", "TRD backward propagation failed"); |
495 | return kFALSE; |
496 | } |
497 | if (fCheckPointLevel > 1) WriteESD(esd, "TRD.back"); |
498 | |
499 | if (!fTOFTracker) { |
500 | Warning("RunTracking", "no TOF tracker"); |
501 | } else { |
502 | // TOF back propagation |
503 | Info("RunTracking", "TOF back propagation"); |
504 | fTOFLoader->LoadDigits("read"); |
505 | TTree* tofTree = fTOFLoader->TreeD(); |
506 | if (!tofTree) { |
507 | Error("RunTracking", "Can't get the TOF digits tree"); |
508 | return kFALSE; |
509 | } |
510 | fTOFTracker->LoadClusters(tofTree); |
511 | if (fTOFTracker->PropagateBack(esd) != 0) { |
512 | Error("RunTracking", "TOF backward propagation failed"); |
513 | return kFALSE; |
514 | } |
515 | if (fCheckPointLevel > 1) WriteESD(esd, "TOF.back"); |
516 | fTOFTracker->UnloadClusters(); |
517 | fTOFLoader->UnloadDigits(); |
518 | } |
596a855f |
519 | |
24f7a148 |
520 | // TRD inward refit |
521 | Info("RunTracking", "TRD inward refit"); |
522 | if (fTRDTracker->RefitInward(esd) != 0) { |
523 | Error("RunTracking", "TRD inward refit failed"); |
524 | return kFALSE; |
525 | } |
526 | if (fCheckPointLevel > 1) WriteESD(esd, "TRD.refit"); |
527 | fTRDTracker->UnloadClusters(); |
528 | fTRDLoader->UnloadRecPoints(); |
529 | |
530 | // TPC inward refit |
531 | Info("RunTracking", "TPC inward refit"); |
532 | if (fTPCTracker->RefitInward(esd) != 0) { |
533 | Error("RunTracking", "TPC inward refit failed"); |
534 | return kFALSE; |
535 | } |
536 | if (fCheckPointLevel > 1) WriteESD(esd, "TPC.refit"); |
537 | |
538 | // ITS inward refit |
539 | Info("RunTracking", "ITS inward refit"); |
540 | if (fITSTracker->RefitInward(esd) != 0) { |
541 | Error("RunTracking", "ITS inward refit failed"); |
542 | return kFALSE; |
543 | } |
544 | if (fCheckPointLevel > 1) WriteESD(esd, "ITS.refit"); |
596a855f |
545 | |
24f7a148 |
546 | } // if TRD tracker |
547 | fITSTracker->UnloadClusters(); |
548 | fITSLoader->UnloadRecPoints(); |
596a855f |
549 | |
24f7a148 |
550 | } // if ITS tracker |
596a855f |
551 | fTPCTracker->UnloadClusters(); |
552 | fTPCLoader->UnloadRecPoints(); |
596a855f |
553 | |
030b532d |
554 | Info("RunTracking", "execution time:"); |
555 | stopwatch.Print(); |
556 | |
596a855f |
557 | return kTRUE; |
558 | } |
559 | |
560 | //_____________________________________________________________________________ |
24f7a148 |
561 | Bool_t AliReconstruction::FillESD(AliESD*& esd, const TString& detectors) |
596a855f |
562 | { |
563 | // fill the event summary data |
564 | |
030b532d |
565 | TStopwatch stopwatch; |
566 | stopwatch.Start(); |
567 | |
596a855f |
568 | TString detStr = detectors; |
c757bafd |
569 | for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) { |
570 | AliReconstructor* reconstructor = |
571 | (AliReconstructor*) fReconstructors[iDet]; |
572 | TString detName = reconstructor->GetDetectorName(); |
573 | if (IsSelected(detName, detStr)) { |
574 | if (!ReadESD(esd, detName.Data())) { |
575 | Info("FillESD", "filling ESD for %s", detName.Data()); |
576 | reconstructor->FillESD(fRunLoader, esd); |
577 | if (fCheckPointLevel > 2) WriteESD(esd, detName.Data()); |
24f7a148 |
578 | } |
596a855f |
579 | } |
580 | } |
581 | |
582 | if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) { |
583 | Error("FillESD", "the following detectors were not found: %s", |
584 | detStr.Data()); |
585 | if (fStopOnError) return kFALSE; |
586 | } |
587 | |
030b532d |
588 | Info("FillESD", "execution time:"); |
589 | stopwatch.Print(); |
590 | |
596a855f |
591 | return kTRUE; |
592 | } |
593 | |
594 | |
595 | //_____________________________________________________________________________ |
596 | Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const |
597 | { |
598 | // check whether detName is contained in detectors |
599 | // if yes, it is removed from detectors |
600 | |
601 | // check if all detectors are selected |
602 | if ((detectors.CompareTo("ALL") == 0) || |
603 | detectors.BeginsWith("ALL ") || |
604 | detectors.EndsWith(" ALL") || |
605 | detectors.Contains(" ALL ")) { |
606 | detectors = "ALL"; |
607 | return kTRUE; |
608 | } |
609 | |
610 | // search for the given detector |
611 | Bool_t result = kFALSE; |
612 | if ((detectors.CompareTo(detName) == 0) || |
613 | detectors.BeginsWith(detName+" ") || |
614 | detectors.EndsWith(" "+detName) || |
615 | detectors.Contains(" "+detName+" ")) { |
616 | detectors.ReplaceAll(detName, ""); |
617 | result = kTRUE; |
618 | } |
619 | |
620 | // clean up the detectors string |
621 | while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " "); |
622 | while (detectors.BeginsWith(" ")) detectors.Remove(0, 1); |
623 | while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1); |
624 | |
625 | return result; |
626 | } |
e583c30d |
627 | |
c757bafd |
628 | //_____________________________________________________________________________ |
629 | AliReconstructor* AliReconstruction::GetReconstructor(const char* detName) const |
630 | { |
631 | // get the reconstructor object for a detector |
632 | |
633 | for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) { |
634 | AliReconstructor* reconstructor = |
635 | (AliReconstructor*) fReconstructors[iDet]; |
636 | if (strcmp(reconstructor->GetDetectorName(), detName) == 0) { |
637 | return reconstructor; |
638 | } |
639 | } |
640 | return NULL; |
641 | } |
642 | |
2257f27e |
643 | //_____________________________________________________________________________ |
644 | Bool_t AliReconstruction::CreateVertexer() |
645 | { |
646 | // create the vertexer |
647 | |
648 | fITSVertexer = NULL; |
c757bafd |
649 | AliReconstructor* itsReconstructor = GetReconstructor("ITS"); |
59697224 |
650 | if (itsReconstructor) { |
651 | fITSVertexer = itsReconstructor->CreateVertexer(fRunLoader); |
2257f27e |
652 | } |
653 | if (!fITSVertexer) { |
654 | Warning("CreateVertexer", "couldn't create a vertexer for ITS"); |
655 | if (fStopOnError) return kFALSE; |
656 | } |
657 | |
658 | return kTRUE; |
659 | } |
660 | |
24f7a148 |
661 | //_____________________________________________________________________________ |
662 | Bool_t AliReconstruction::CreateTrackers() |
663 | { |
664 | // get the loaders and create the trackers |
665 | |
24f7a148 |
666 | fITSTracker = NULL; |
667 | fITSLoader = fRunLoader->GetLoader("ITSLoader"); |
668 | if (!fITSLoader) { |
669 | Warning("CreateTrackers", "no ITS loader found"); |
670 | if (fStopOnError) return kFALSE; |
671 | } else { |
c757bafd |
672 | AliReconstructor* itsReconstructor = GetReconstructor("ITS"); |
59697224 |
673 | if (itsReconstructor) { |
674 | fITSTracker = itsReconstructor->CreateTracker(fRunLoader); |
24f7a148 |
675 | } |
676 | if (!fITSTracker) { |
677 | Warning("CreateTrackers", "couldn't create a tracker for ITS"); |
678 | if (fStopOnError) return kFALSE; |
679 | } |
680 | } |
681 | |
682 | fTPCTracker = NULL; |
683 | fTPCLoader = fRunLoader->GetLoader("TPCLoader"); |
684 | if (!fTPCLoader) { |
685 | Error("CreateTrackers", "no TPC loader found"); |
686 | if (fStopOnError) return kFALSE; |
687 | } else { |
c757bafd |
688 | AliReconstructor* tpcReconstructor = GetReconstructor("TPC"); |
59697224 |
689 | if (tpcReconstructor) { |
690 | fTPCTracker = tpcReconstructor->CreateTracker(fRunLoader); |
24f7a148 |
691 | } |
692 | if (!fTPCTracker) { |
693 | Error("CreateTrackers", "couldn't create a tracker for TPC"); |
694 | if (fStopOnError) return kFALSE; |
695 | } |
696 | } |
697 | |
698 | fTRDTracker = NULL; |
699 | fTRDLoader = fRunLoader->GetLoader("TRDLoader"); |
700 | if (!fTRDLoader) { |
701 | Warning("CreateTrackers", "no TRD loader found"); |
702 | if (fStopOnError) return kFALSE; |
703 | } else { |
c757bafd |
704 | AliReconstructor* trdReconstructor = GetReconstructor("TRD"); |
59697224 |
705 | if (trdReconstructor) { |
706 | fTRDTracker = trdReconstructor->CreateTracker(fRunLoader); |
24f7a148 |
707 | } |
708 | if (!fTRDTracker) { |
709 | Warning("CreateTrackers", "couldn't create a tracker for TRD"); |
710 | if (fStopOnError) return kFALSE; |
711 | } |
712 | } |
713 | |
714 | fTOFTracker = NULL; |
715 | fTOFLoader = fRunLoader->GetLoader("TOFLoader"); |
716 | if (!fTOFLoader) { |
717 | Warning("CreateTrackers", "no TOF loader found"); |
718 | if (fStopOnError) return kFALSE; |
719 | } else { |
c757bafd |
720 | AliReconstructor* tofReconstructor = GetReconstructor("TOF"); |
59697224 |
721 | if (tofReconstructor) { |
722 | fTOFTracker = tofReconstructor->CreateTracker(fRunLoader); |
24f7a148 |
723 | } |
724 | if (!fTOFTracker) { |
725 | Warning("CreateTrackers", "couldn't create a tracker for TOF"); |
726 | if (fStopOnError) return kFALSE; |
727 | } |
728 | } |
729 | |
730 | return kTRUE; |
731 | } |
732 | |
e583c30d |
733 | //_____________________________________________________________________________ |
734 | void AliReconstruction::CleanUp(TFile* file) |
735 | { |
736 | // delete trackers and the run loader and close and delete the file |
737 | |
59697224 |
738 | fReconstructors.Delete(); |
739 | |
2257f27e |
740 | delete fITSVertexer; |
741 | fITSVertexer = NULL; |
e583c30d |
742 | delete fITSTracker; |
743 | fITSTracker = NULL; |
744 | delete fTPCTracker; |
745 | fTPCTracker = NULL; |
746 | delete fTRDTracker; |
747 | fTRDTracker = NULL; |
748 | delete fTOFTracker; |
749 | fTOFTracker = NULL; |
750 | |
751 | delete fRunLoader; |
752 | fRunLoader = NULL; |
753 | |
754 | if (file) { |
755 | file->Close(); |
756 | delete file; |
757 | } |
758 | } |
24f7a148 |
759 | |
760 | |
761 | //_____________________________________________________________________________ |
762 | Bool_t AliReconstruction::ReadESD(AliESD*& esd, const char* recStep) const |
763 | { |
764 | // read the ESD event from a file |
765 | |
766 | if (!esd) return kFALSE; |
767 | char fileName[256]; |
768 | sprintf(fileName, "ESD_%d.%d_%s.root", |
769 | esd->GetRunNumber(), esd->GetEventNumber(), recStep); |
770 | if (gSystem->AccessPathName(fileName)) return kFALSE; |
771 | |
772 | Info("ReadESD", "reading ESD from file %s", fileName); |
773 | TFile* file = TFile::Open(fileName); |
774 | if (!file || !file->IsOpen()) { |
775 | Error("ReadESD", "opening %s failed", fileName); |
776 | delete file; |
777 | return kFALSE; |
778 | } |
779 | |
780 | gROOT->cd(); |
781 | delete esd; |
782 | esd = (AliESD*) file->Get("ESD"); |
783 | file->Close(); |
784 | delete file; |
785 | return kTRUE; |
786 | } |
787 | |
788 | //_____________________________________________________________________________ |
789 | void AliReconstruction::WriteESD(AliESD* esd, const char* recStep) const |
790 | { |
791 | // write the ESD event to a file |
792 | |
793 | if (!esd) return; |
794 | char fileName[256]; |
795 | sprintf(fileName, "ESD_%d.%d_%s.root", |
796 | esd->GetRunNumber(), esd->GetEventNumber(), recStep); |
797 | |
798 | Info("WriteESD", "writing ESD to file %s", fileName); |
799 | TFile* file = TFile::Open(fileName, "recreate"); |
800 | if (!file || !file->IsOpen()) { |
801 | Error("WriteESD", "opening %s failed", fileName); |
802 | } else { |
803 | esd->Write("ESD"); |
804 | file->Close(); |
805 | } |
806 | delete file; |
807 | } |