]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - STEER/AliReconstruction.cxx
Bug correction
[u/mrichter/AliRoot.git] / STEER / AliReconstruction.cxx
... / ...
CommitLineData
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// //
30// If the input to the reconstruction are not simulated digits but raw data, //
31// this can be specified by an argument of the Run method or by the method //
32// //
33// rec.SetInput("..."); //
34// //
35// The input formats and the corresponding argument are: //
36// - DDL raw data files: directory name, ends with "/" //
37// - raw data root file: root file name, extension ".root" //
38// - raw data DATE file: DATE file name, any other non-empty string //
39// - MC root files : empty string, default //
40// //
41// The name of the galice file can be changed from the default //
42// "galice.root" by passing it as argument to the AliReconstruction //
43// constructor or by //
44// //
45// rec.SetGAliceFile("..."); //
46// //
47// The local reconstruction can be switched on or off for individual //
48// detectors by //
49// //
50// rec.SetRunLocalReconstruction("..."); //
51// //
52// The argument is a (case sensitive) string with the names of the //
53// detectors separated by a space. The special string "ALL" selects all //
54// available detectors. This is the default. //
55// //
56// The reconstruction of the primary vertex position can be switched off by //
57// //
58// rec.SetRunVertexFinder(kFALSE); //
59// //
60// The tracking in ITS, TPC and TRD and the creation of ESD tracks can be //
61// switched off by //
62// //
63// rec.SetRunTracking(kFALSE); //
64// //
65// The filling of additional ESD information can be steered by //
66// //
67// rec.SetFillESD("..."); //
68// //
69// Again, the string specifies the list of detectors. The default is "ALL". //
70// //
71// The reconstruction requires digits or raw data as input. For the creation //
72// of digits and raw data have a look at the class AliSimulation. //
73// //
74// For debug purposes the method SetCheckPointLevel can be used. If the //
75// argument is greater than 0, files with ESD events will be written after //
76// selected steps of the reconstruction for each event: //
77// level 1: after tracking and after filling of ESD (final) //
78// level 2: in addition after each tracking step //
79// level 3: in addition after the filling of ESD for each detector //
80// If a final check point file exists for an event, this event will be //
81// skipped in the reconstruction. The tracking and the filling of ESD for //
82// a detector will be skipped as well, if the corresponding check point //
83// file exists. The ESD event will then be loaded from the file instead. //
84// //
85///////////////////////////////////////////////////////////////////////////////
86
87#include <TArrayF.h>
88#include <TFile.h>
89#include <TSystem.h>
90#include <TROOT.h>
91#include <TPluginManager.h>
92#include <TStopwatch.h>
93
94#include "AliReconstruction.h"
95#include "AliLog.h"
96#include "AliRunLoader.h"
97#include "AliRun.h"
98#include "AliRawReaderFile.h"
99#include "AliRawReaderDate.h"
100#include "AliRawReaderRoot.h"
101#include "AliTracker.h"
102#include "AliESD.h"
103#include "AliESDVertex.h"
104#include "AliVertexer.h"
105#include "AliHeader.h"
106#include "AliGenEventHeader.h"
107#include "AliESDpid.h"
108#include "AliMagF.h"
109
110ClassImp(AliReconstruction)
111
112
113//_____________________________________________________________________________
114const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "RICH", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "START", "VZERO", "CRT", "HLT"};
115
116//_____________________________________________________________________________
117AliReconstruction::AliReconstruction(const char* gAliceFilename,
118 const char* name, const char* title) :
119 TNamed(name, title),
120
121 fRunLocalReconstruction("ALL"),
122 fRunVertexFinder(kTRUE),
123 fRunTracking(kTRUE),
124 fFillESD("ALL"),
125 fGAliceFileName(gAliceFilename),
126 fInput(""),
127 fStopOnError(kFALSE),
128 fCheckPointLevel(0),
129
130 fRunLoader(NULL),
131 fRawReader(NULL),
132 fITSLoader(NULL),
133 fITSVertexer(NULL),
134 fITSTracker(NULL),
135 fTPCLoader(NULL),
136 fTPCTracker(NULL),
137 fTRDLoader(NULL),
138 fTRDTracker(NULL),
139 fTOFLoader(NULL),
140 fTOFTracker(NULL),
141 fPHOSLoader(NULL),
142 fPHOSTracker(NULL),
143 fEMCALLoader(NULL),
144 fEMCALTracker(NULL),
145 fRICHLoader(NULL),
146 fRICHTracker(NULL),
147
148 fReconstructors(),
149 fOptions()
150{
151// create reconstruction object with default parameters
152
153}
154
155//_____________________________________________________________________________
156AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
157 TNamed(rec),
158
159 fRunLocalReconstruction(rec.fRunLocalReconstruction),
160 fRunVertexFinder(rec.fRunVertexFinder),
161 fRunTracking(rec.fRunTracking),
162 fFillESD(rec.fFillESD),
163 fGAliceFileName(rec.fGAliceFileName),
164 fInput(rec.fInput),
165 fStopOnError(rec.fStopOnError),
166 fCheckPointLevel(0),
167
168 fRunLoader(NULL),
169 fRawReader(NULL),
170 fITSLoader(NULL),
171 fITSVertexer(NULL),
172 fITSTracker(NULL),
173 fTPCLoader(NULL),
174 fTPCTracker(NULL),
175 fTRDLoader(NULL),
176 fTRDTracker(NULL),
177 fTOFLoader(NULL),
178 fTOFTracker(NULL),
179 fPHOSLoader(NULL),
180 fPHOSTracker(NULL),
181 fEMCALLoader(NULL),
182 fEMCALTracker(NULL),
183 fRICHLoader(NULL),
184 fRICHTracker(NULL),
185
186 fReconstructors(),
187 fOptions()
188{
189// copy constructor
190
191 for (Int_t i = 0; i < fOptions.GetEntriesFast(); i++) {
192 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
193 }
194}
195
196//_____________________________________________________________________________
197AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
198{
199// assignment operator
200
201 this->~AliReconstruction();
202 new(this) AliReconstruction(rec);
203 return *this;
204}
205
206//_____________________________________________________________________________
207AliReconstruction::~AliReconstruction()
208{
209// clean up
210
211 CleanUp();
212 fOptions.Delete();
213}
214
215
216//_____________________________________________________________________________
217void AliReconstruction::SetGAliceFile(const char* fileName)
218{
219// set the name of the galice file
220
221 fGAliceFileName = fileName;
222}
223
224//_____________________________________________________________________________
225void AliReconstruction::SetOption(const char* detector, const char* option)
226{
227// set options for the reconstruction of a detector
228
229 TObject* obj = fOptions.FindObject(detector);
230 if (obj) fOptions.Remove(obj);
231 fOptions.Add(new TNamed(detector, option));
232}
233
234
235//_____________________________________________________________________________
236Bool_t AliReconstruction::Run(const char* input)
237{
238// run the reconstruction
239
240 // set the input
241 if (!input) input = fInput.Data();
242 TString fileName(input);
243 if (fileName.EndsWith("/")) {
244 fRawReader = new AliRawReaderFile(fileName);
245 } else if (fileName.EndsWith(".root")) {
246 fRawReader = new AliRawReaderRoot(fileName);
247 } else if (!fileName.IsNull()) {
248 fRawReader = new AliRawReaderDate(fileName);
249 fRawReader->SelectEvents(7);
250 }
251
252 // open the run loader
253 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
254 if (!fRunLoader) {
255 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
256 CleanUp();
257 return kFALSE;
258 }
259 fRunLoader->LoadgAlice();
260 AliRun* aliRun = fRunLoader->GetAliRun();
261 if (!aliRun) {
262 AliError(Form("no gAlice object found in file %s",
263 fGAliceFileName.Data()));
264 CleanUp();
265 return kFALSE;
266 }
267 gAlice = aliRun;
268 AliTracker::SetFieldMap(gAlice->Field());
269
270 // load the reconstructor objects
271 TPluginManager* pluginManager = gROOT->GetPluginManager();
272 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
273 TString detName = fgkDetectorName[iDet];
274 TString recName = "Ali" + detName + "Reconstructor";
275 if (!gAlice->GetDetector(detName) && detName != "HLT") continue;
276
277 if(detName == "HLT") {
278 if (!gROOT->GetClass("AliLevel3")) {
279 gSystem->Load("libAliL3Src.so");
280 gSystem->Load("libAliL3Misc.so");
281 gSystem->Load("libAliL3Hough.so");
282 gSystem->Load("libAliL3Comp.so");
283 }
284 }
285
286 AliReconstructor* reconstructor = NULL;
287 // first check if a plugin is defined for the reconstructor
288 TPluginHandler* pluginHandler =
289 pluginManager->FindHandler("AliReconstructor", detName);
290 // if not, but the reconstructor class is implemented, add a plugin for it
291 if (!pluginHandler && gROOT->GetClass(recName.Data())) {
292 AliDebug(1, Form("defining plugin for %s", recName.Data()));
293 pluginManager->AddHandler("AliReconstructor", detName,
294 recName, detName, recName + "()");
295 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
296 }
297 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
298 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
299 }
300 // if there is no reconstructor class for the detector use the dummy one
301 if (!reconstructor && gAlice->GetDetector(detName)) {
302 AliDebug(1, Form("using dummy reconstructor for %s", detName.Data()));
303 reconstructor = new AliDummyReconstructor(gAlice->GetDetector(detName));
304 }
305 if (reconstructor) {
306 TObject* obj = fOptions.FindObject(detName.Data());
307 if (obj) reconstructor->SetOption(obj->GetTitle());
308 fReconstructors.Add(reconstructor);
309 }
310 }
311
312 // local reconstruction
313 if (!fRunLocalReconstruction.IsNull()) {
314 if (!RunLocalReconstruction(fRunLocalReconstruction)) {
315 if (fStopOnError) {CleanUp(); return kFALSE;}
316 }
317 }
318 if (!fRunVertexFinder && !fRunTracking && fFillESD.IsNull()) return kTRUE;
319
320 // get vertexer
321 if (fRunVertexFinder && !CreateVertexer()) {
322 if (fStopOnError) {
323 CleanUp();
324 return kFALSE;
325 }
326 }
327
328 // get loaders and trackers
329 if (fRunTracking && !CreateTrackers()) {
330 if (fStopOnError) {
331 CleanUp();
332 return kFALSE;
333 }
334 }
335
336 // create the ESD output file and tree
337 TFile* file = TFile::Open("AliESDs.root", "RECREATE");
338 if (!file->IsOpen()) {
339 AliError("opening AliESDs.root failed");
340 if (fStopOnError) {CleanUp(file); return kFALSE;}
341 }
342 AliESD* esd = new AliESD;
343 TTree* tree = new TTree("esdTree", "Tree with ESD objects");
344 tree->Branch("ESD", "AliESD", &esd);
345 delete esd;
346 gROOT->cd();
347
348 // loop over events
349 if (fRawReader) fRawReader->RewindEvents();
350 for (Int_t iEvent = 0; iEvent < fRunLoader->GetNumberOfEvents(); iEvent++) {
351 AliInfo(Form("processing event %d", iEvent));
352 fRunLoader->GetEvent(iEvent);
353 if (fRawReader) fRawReader->NextEvent();
354
355 char fileName[256];
356 sprintf(fileName, "ESD_%d.%d_final.root",
357 aliRun->GetRunNumber(), aliRun->GetEvNumber());
358 if (!gSystem->AccessPathName(fileName)) continue;
359
360 esd = new AliESD;
361 esd->SetRunNumber(aliRun->GetRunNumber());
362 esd->SetEventNumber(aliRun->GetEvNumber());
363 esd->SetMagneticField(aliRun->Field()->SolenoidField());
364
365 // vertex finder
366 if (fRunVertexFinder) {
367 if (!ReadESD(esd, "vertex")) {
368 if (!RunVertexFinder(esd)) {
369 if (fStopOnError) {CleanUp(file); return kFALSE;}
370 }
371 if (fCheckPointLevel > 0) WriteESD(esd, "vertex");
372 }
373 }
374
375 // barrel tracking
376 if (fRunTracking) {
377 if (!ReadESD(esd, "tracking")) {
378 if (!RunTracking(esd)) {
379 if (fStopOnError) {CleanUp(file); return kFALSE;}
380 }
381 if (fCheckPointLevel > 0) WriteESD(esd, "tracking");
382 }
383 }
384
385 // fill ESD
386 if (!fFillESD.IsNull()) {
387 if (!FillESD(esd, fFillESD)) {
388 if (fStopOnError) {CleanUp(file); return kFALSE;}
389 }
390 }
391
392 // combined PID
393 AliESDpid::MakePID(esd);
394 if (fCheckPointLevel > 1) WriteESD(esd, "PID");
395
396 // write ESD
397 tree->Fill();
398
399 if (fCheckPointLevel > 0) WriteESD(esd, "final");
400 delete esd;
401 }
402
403 file->cd();
404 tree->Write();
405 CleanUp(file);
406
407 return kTRUE;
408}
409
410
411//_____________________________________________________________________________
412Bool_t AliReconstruction::RunLocalReconstruction(const TString& detectors)
413{
414// run the local reconstruction
415
416 TStopwatch stopwatch;
417 stopwatch.Start();
418
419 TString detStr = detectors;
420 for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
421 AliReconstructor* reconstructor =
422 (AliReconstructor*) fReconstructors[iDet];
423 TString detName = reconstructor->GetDetectorName();
424 if (IsSelected(detName, detStr)) {
425 AliInfo(Form("running reconstruction for %s", detName.Data()));
426 TStopwatch stopwatchDet;
427 stopwatchDet.Start();
428 if (fRawReader) {
429 fRawReader->RewindEvents();
430 reconstructor->Reconstruct(fRunLoader, fRawReader);
431 } else {
432 reconstructor->Reconstruct(fRunLoader);
433 }
434 AliInfo(Form("execution time for %s:", detName.Data()));
435 ToAliInfo(stopwatchDet.Print());
436 }
437 }
438
439 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
440 AliError(Form("the following detectors were not found: %s",
441 detStr.Data()));
442 if (fStopOnError) return kFALSE;
443 }
444
445 AliInfo("execution time:");
446 ToAliInfo(stopwatch.Print());
447
448 return kTRUE;
449}
450
451//_____________________________________________________________________________
452Bool_t AliReconstruction::RunVertexFinder(AliESD*& esd)
453{
454// run the barrel tracking
455
456 TStopwatch stopwatch;
457 stopwatch.Start();
458
459 AliESDVertex* vertex = NULL;
460 Double_t vtxPos[3] = {0, 0, 0};
461 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
462 TArrayF mcVertex(3);
463 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
464 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
465 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
466 }
467
468 if (fITSVertexer) {
469 AliInfo("running the ITS vertex finder");
470 vertex = fITSVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber());
471 if(!vertex){
472 AliWarning("Vertex not found");
473 vertex = new AliESDVertex();
474 }
475 else {
476 vertex->SetTruePos(vtxPos); // store also the vertex from MC
477 }
478
479 } else {
480 AliInfo("getting the primary vertex from MC");
481 vertex = new AliESDVertex(vtxPos, vtxErr);
482 }
483
484 if (vertex) {
485 vertex->GetXYZ(vtxPos);
486 vertex->GetSigmaXYZ(vtxErr);
487 } else {
488 AliWarning("no vertex reconstructed");
489 vertex = new AliESDVertex(vtxPos, vtxErr);
490 }
491 esd->SetVertex(vertex);
492 if (fITSTracker) fITSTracker->SetVertex(vtxPos, vtxErr);
493 if (fTPCTracker) fTPCTracker->SetVertex(vtxPos, vtxErr);
494 if (fTRDTracker) fTRDTracker->SetVertex(vtxPos, vtxErr);
495 delete vertex;
496
497 AliInfo("execution time:");
498 ToAliInfo(stopwatch.Print());
499
500 return kTRUE;
501}
502
503//_____________________________________________________________________________
504Bool_t AliReconstruction::RunTracking(AliESD*& esd)
505{
506// run the barrel tracking
507
508 TStopwatch stopwatch;
509 stopwatch.Start();
510
511 if (!fTPCTracker) {
512 AliError("no TPC tracker");
513 return kFALSE;
514 }
515 AliInfo("running tracking");
516
517 // TPC tracking
518 AliDebug(1, "TPC tracking");
519 fTPCLoader->LoadRecPoints("read");
520 TTree* tpcTree = fTPCLoader->TreeR();
521 if (!tpcTree) {
522 AliError("Can't get the TPC cluster tree");
523 return kFALSE;
524 }
525 fTPCTracker->LoadClusters(tpcTree);
526 if (fTPCTracker->Clusters2Tracks(esd) != 0) {
527 AliError("TPC Clusters2Tracks failed");
528 return kFALSE;
529 }
530 if (fCheckPointLevel > 1) WriteESD(esd, "TPC.tracking");
531
532 if (!fITSTracker) {
533 AliWarning("no ITS tracker");
534 } else {
535
536 GetReconstructor("TPC")->FillESD(fRunLoader, esd); // preliminary
537 AliESDpid::MakePID(esd); // PID for the ITS tracker
538
539 // ITS tracking
540 AliDebug(1, "ITS tracking");
541 fITSLoader->LoadRecPoints("read");
542 TTree* itsTree = fITSLoader->TreeR();
543 if (!itsTree) {
544 Error("RunTracking", "Can't get the ITS cluster tree");
545 return kFALSE;
546 }
547 fITSTracker->LoadClusters(itsTree);
548 if (fITSTracker->Clusters2Tracks(esd) != 0) {
549 AliError("ITS Clusters2Tracks failed");
550 return kFALSE;
551 }
552 if (fCheckPointLevel > 1) WriteESD(esd, "ITS.tracking");
553
554 if (fTRDTracker||fTOFTracker||fPHOSTracker||fEMCALTracker||fRICHTracker) {
555 // ITS back propagation
556 AliDebug(1, "ITS back propagation");
557 if (fITSTracker->PropagateBack(esd) != 0) {
558 AliError("ITS backward propagation failed");
559 return kFALSE;
560 }
561 if (fCheckPointLevel > 1) WriteESD(esd, "ITS.back");
562
563 // TPC back propagation
564 AliDebug(1, "TPC back propagation");
565 if (fTPCTracker->PropagateBack(esd) != 0) {
566 AliError("TPC backward propagation failed");
567 return kFALSE;
568 }
569 if (fCheckPointLevel > 1) WriteESD(esd, "TPC.back");
570
571 if (!fTRDTracker) {
572 AliWarning("no TRD tracker");
573 } else {
574 // TRD back propagation
575 AliDebug(1, "TRD back propagation");
576 fTRDLoader->LoadRecPoints("read");
577 TTree* trdTree = fTRDLoader->TreeR();
578 if (!trdTree) {
579 AliError("Can't get the TRD cluster tree");
580 return kFALSE;
581 }
582 fTRDTracker->LoadClusters(trdTree);
583 if (fTRDTracker->PropagateBack(esd) != 0) {
584 AliError("TRD backward propagation failed");
585 return kFALSE;
586 }
587 if (fCheckPointLevel > 1) WriteESD(esd, "TRD.back");
588 }
589
590 if (!fTOFTracker) {
591 AliWarning("no TOF tracker");
592 } else {
593 // TOF back propagation
594 AliDebug(1, "TOF back propagation");
595 fTOFLoader->LoadDigits("read");
596 TTree* tofTree = fTOFLoader->TreeD();
597 if (!tofTree) {
598 AliError("Can't get the TOF digits tree");
599 return kFALSE;
600 }
601 fTOFTracker->LoadClusters(tofTree);
602 if (fTOFTracker->PropagateBack(esd) != 0) {
603 AliError("TOF backward propagation failed");
604 return kFALSE;
605 }
606 if (fCheckPointLevel > 1) WriteESD(esd, "TOF.back");
607 fTOFTracker->UnloadClusters();
608 fTOFLoader->UnloadDigits();
609
610 if (!fPHOSTracker) {
611 AliWarning("no PHOS tracker");
612 } else {
613 // PHOS back propagation
614 AliDebug(1, "PHOS back propagation");
615 fPHOSLoader->LoadRecPoints("read");
616 TTree* phosTree = fPHOSLoader->TreeR();
617 if (!phosTree) {
618 AliError("Can't get the PHOS cluster tree");
619 return kFALSE;
620 }
621 fPHOSTracker->LoadClusters(phosTree);
622 if (fPHOSTracker->PropagateBack(esd) != 0) {
623 AliError("PHOS backward propagation failed");
624 return kFALSE;
625 }
626 if (fCheckPointLevel > 1) WriteESD(esd, "PHOS.back");
627 fPHOSTracker->UnloadClusters();
628 fPHOSLoader->UnloadRecPoints();
629 }
630
631 if (!fEMCALTracker) {
632 AliWarning("no EMCAL tracker");
633 } else {
634 // EMCAL back propagation
635 AliDebug(1, "EMCAL back propagation");
636 fEMCALLoader->LoadRecPoints("read");
637 TTree* emcalTree = fEMCALLoader->TreeR();
638 if (!emcalTree) {
639 AliError("Can't get the EMCAL cluster tree");
640 return kFALSE;
641 }
642 fEMCALTracker->LoadClusters(emcalTree);
643 if (fEMCALTracker->PropagateBack(esd) != 0) {
644 AliError("EMCAL backward propagation failed");
645 return kFALSE;
646 }
647 if (fCheckPointLevel > 1) WriteESD(esd, "EMCAL.back");
648 fEMCALTracker->UnloadClusters();
649 fEMCALLoader->UnloadRecPoints();
650 }
651
652 if (!fRICHTracker) {
653 AliWarning("no RICH tracker");
654 } else {
655 // RICH back propagation
656 AliDebug(1, "RICH back propagation");
657 fRICHLoader->LoadRecPoints("read");
658 TTree* richTree = fRICHLoader->TreeR();
659 if (!richTree) {
660 AliError("Can't get the RICH cluster tree");
661 return kFALSE;
662 }
663 fRICHTracker->LoadClusters(richTree);
664 if (fRICHTracker->PropagateBack(esd) != 0) {
665 AliError("RICH backward propagation failed");
666 return kFALSE;
667 }
668 if (fCheckPointLevel > 1) WriteESD(esd, "RICH.back");
669 fRICHTracker->UnloadClusters();
670 fRICHLoader->UnloadRecPoints();
671 }
672 }
673
674 if (fTRDTracker) {
675 // TRD inward refit
676 AliDebug(1, "TRD inward refit");
677 if (fTRDTracker->RefitInward(esd) != 0) {
678 AliError("TRD inward refit failed");
679 return kFALSE;
680 }
681 if (fCheckPointLevel > 1) WriteESD(esd, "TRD.refit");
682 fTRDTracker->UnloadClusters();
683 fTRDLoader->UnloadRecPoints();
684 }
685
686 // TPC inward refit
687 AliInfo("TPC inward refit");
688 if (fTPCTracker->RefitInward(esd) != 0) {
689 AliError("TPC inward refit failed");
690 return kFALSE;
691 }
692 if (fCheckPointLevel > 1) WriteESD(esd, "TPC.refit");
693
694 // ITS inward refit
695 AliInfo("ITS inward refit");
696 if (fITSTracker->RefitInward(esd) != 0) {
697 AliError("ITS inward refit failed");
698 return kFALSE;
699 }
700 if (fCheckPointLevel > 1) WriteESD(esd, "ITS.refit");
701
702 } // if TRD tracker or TOF tracker or PHOS tracker ...
703 fITSTracker->UnloadClusters();
704 fITSLoader->UnloadRecPoints();
705
706 } // if ITS tracker
707 fTPCTracker->UnloadClusters();
708 fTPCLoader->UnloadRecPoints();
709
710 AliInfo("execution time:");
711 ToAliInfo(stopwatch.Print());
712
713 return kTRUE;
714}
715
716//_____________________________________________________________________________
717Bool_t AliReconstruction::FillESD(AliESD*& esd, const TString& detectors)
718{
719// fill the event summary data
720
721 TStopwatch stopwatch;
722 stopwatch.Start();
723 AliInfo("filling ESD");
724
725 TString detStr = detectors;
726 for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
727 AliReconstructor* reconstructor =
728 (AliReconstructor*) fReconstructors[iDet];
729 TString detName = reconstructor->GetDetectorName();
730 if (IsSelected(detName, detStr)) {
731 if (!ReadESD(esd, detName.Data())) {
732 AliDebug(1, Form("filling ESD for %s", detName.Data()));
733 if (fRawReader) {
734 reconstructor->FillESD(fRunLoader, fRawReader, esd);
735 } else {
736 reconstructor->FillESD(fRunLoader, esd);
737 }
738 if (fCheckPointLevel > 2) WriteESD(esd, detName.Data());
739 }
740 }
741 }
742
743 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
744 AliError(Form("the following detectors were not found: %s",
745 detStr.Data()));
746 if (fStopOnError) return kFALSE;
747 }
748
749 AliInfo("execution time:");
750 ToAliInfo(stopwatch.Print());
751
752 return kTRUE;
753}
754
755
756//_____________________________________________________________________________
757Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
758{
759// check whether detName is contained in detectors
760// if yes, it is removed from detectors
761
762 // check if all detectors are selected
763 if ((detectors.CompareTo("ALL") == 0) ||
764 detectors.BeginsWith("ALL ") ||
765 detectors.EndsWith(" ALL") ||
766 detectors.Contains(" ALL ")) {
767 detectors = "ALL";
768 return kTRUE;
769 }
770
771 // search for the given detector
772 Bool_t result = kFALSE;
773 if ((detectors.CompareTo(detName) == 0) ||
774 detectors.BeginsWith(detName+" ") ||
775 detectors.EndsWith(" "+detName) ||
776 detectors.Contains(" "+detName+" ")) {
777 detectors.ReplaceAll(detName, "");
778 result = kTRUE;
779 }
780
781 // clean up the detectors string
782 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
783 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
784 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
785
786 return result;
787}
788
789//_____________________________________________________________________________
790AliReconstructor* AliReconstruction::GetReconstructor(const char* detName) const
791{
792// get the reconstructor object for a detector
793
794 for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
795 AliReconstructor* reconstructor =
796 (AliReconstructor*) fReconstructors[iDet];
797 if (strcmp(reconstructor->GetDetectorName(), detName) == 0) {
798 return reconstructor;
799 }
800 }
801 return NULL;
802}
803
804//_____________________________________________________________________________
805Bool_t AliReconstruction::CreateVertexer()
806{
807// create the vertexer
808
809 fITSVertexer = NULL;
810 AliReconstructor* itsReconstructor = GetReconstructor("ITS");
811 if (itsReconstructor) {
812 fITSVertexer = itsReconstructor->CreateVertexer(fRunLoader);
813 }
814 if (!fITSVertexer) {
815 AliWarning("couldn't create a vertexer for ITS");
816 if (fStopOnError) return kFALSE;
817 }
818
819 return kTRUE;
820}
821
822//_____________________________________________________________________________
823Bool_t AliReconstruction::CreateTrackers()
824{
825// get the loaders and create the trackers
826
827 fITSTracker = NULL;
828 fITSLoader = fRunLoader->GetLoader("ITSLoader");
829 if (!fITSLoader) {
830 AliWarning("no ITS loader found");
831 if (fStopOnError) return kFALSE;
832 } else {
833 AliReconstructor* itsReconstructor = GetReconstructor("ITS");
834 if (itsReconstructor) {
835 fITSTracker = itsReconstructor->CreateTracker(fRunLoader);
836 }
837 if (!fITSTracker) {
838 AliWarning("couldn't create a tracker for ITS");
839 if (fStopOnError) return kFALSE;
840 }
841 }
842
843 fTPCTracker = NULL;
844 fTPCLoader = fRunLoader->GetLoader("TPCLoader");
845 if (!fTPCLoader) {
846 AliError("no TPC loader found");
847 if (fStopOnError) return kFALSE;
848 } else {
849 AliReconstructor* tpcReconstructor = GetReconstructor("TPC");
850 if (tpcReconstructor) {
851 fTPCTracker = tpcReconstructor->CreateTracker(fRunLoader);
852 }
853 if (!fTPCTracker) {
854 AliError("couldn't create a tracker for TPC");
855 if (fStopOnError) return kFALSE;
856 }
857 }
858
859 fTRDTracker = NULL;
860 fTRDLoader = fRunLoader->GetLoader("TRDLoader");
861 if (!fTRDLoader) {
862 AliWarning("no TRD loader found");
863 if (fStopOnError) return kFALSE;
864 } else {
865 AliReconstructor* trdReconstructor = GetReconstructor("TRD");
866 if (trdReconstructor) {
867 fTRDTracker = trdReconstructor->CreateTracker(fRunLoader);
868 }
869 if (!fTRDTracker) {
870 AliWarning("couldn't create a tracker for TRD");
871 if (fStopOnError) return kFALSE;
872 }
873 }
874
875 fTOFTracker = NULL;
876 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
877 if (!fTOFLoader) {
878 AliWarning("no TOF loader found");
879 if (fStopOnError) return kFALSE;
880 } else {
881 AliReconstructor* tofReconstructor = GetReconstructor("TOF");
882 if (tofReconstructor) {
883 fTOFTracker = tofReconstructor->CreateTracker(fRunLoader);
884 }
885 if (!fTOFTracker) {
886 AliWarning("couldn't create a tracker for TOF");
887 if (fStopOnError) return kFALSE;
888 }
889 }
890
891 fPHOSTracker = NULL;
892 fPHOSLoader = fRunLoader->GetLoader("PHOSLoader");
893 if (!fPHOSLoader) {
894 AliWarning("no PHOS loader found");
895 if (fStopOnError) return kFALSE;
896 } else {
897 AliReconstructor* phosReconstructor = GetReconstructor("PHOS");
898 if (phosReconstructor) {
899 fPHOSTracker = phosReconstructor->CreateTracker(fRunLoader);
900 }
901 if (!fPHOSTracker) {
902 AliWarning("couldn't create a tracker for PHOS");
903 if (fStopOnError) return kFALSE;
904 }
905 }
906
907 fEMCALTracker = NULL;
908 fEMCALLoader = fRunLoader->GetLoader("EMCALLoader");
909 if (!fEMCALLoader) {
910 AliWarning("no EMCAL loader found");
911 if (fStopOnError) return kFALSE;
912 } else {
913 AliReconstructor* emcalReconstructor = GetReconstructor("EMCAL");
914 if (emcalReconstructor) {
915 fEMCALTracker = emcalReconstructor->CreateTracker(fRunLoader);
916 }
917 if (!fEMCALTracker) {
918 AliWarning("couldn't create a tracker for EMCAL");
919 if (fStopOnError) return kFALSE;
920 }
921 }
922
923 fRICHTracker = NULL;
924 fRICHLoader = fRunLoader->GetLoader("RICHLoader");
925 if (!fRICHLoader) {
926 AliWarning("no RICH loader found");
927 if (fStopOnError) return kFALSE;
928 } else {
929 AliReconstructor* tofReconstructor = GetReconstructor("RICH");
930 if (tofReconstructor) {
931 fRICHTracker = tofReconstructor->CreateTracker(fRunLoader);
932 }
933 if (!fRICHTracker) {
934 AliWarning("couldn't create a tracker for RICH");
935 if (fStopOnError) return kFALSE;
936 }
937 }
938
939 return kTRUE;
940}
941
942//_____________________________________________________________________________
943void AliReconstruction::CleanUp(TFile* file)
944{
945// delete trackers and the run loader and close and delete the file
946
947 fReconstructors.Delete();
948
949 delete fITSVertexer;
950 fITSVertexer = NULL;
951 delete fITSTracker;
952 fITSTracker = NULL;
953 delete fTPCTracker;
954 fTPCTracker = NULL;
955 delete fTRDTracker;
956 fTRDTracker = NULL;
957 delete fTOFTracker;
958 fTOFTracker = NULL;
959 delete fPHOSTracker;
960 fPHOSTracker = NULL;
961 delete fEMCALTracker;
962 fEMCALTracker = NULL;
963 delete fRICHTracker;
964 fRICHTracker = NULL;
965
966 delete fRunLoader;
967 fRunLoader = NULL;
968 delete fRawReader;
969 fRawReader = NULL;
970
971 if (file) {
972 file->Close();
973 delete file;
974 }
975}
976
977
978//_____________________________________________________________________________
979Bool_t AliReconstruction::ReadESD(AliESD*& esd, const char* recStep) const
980{
981// read the ESD event from a file
982
983 if (!esd) return kFALSE;
984 char fileName[256];
985 sprintf(fileName, "ESD_%d.%d_%s.root",
986 esd->GetRunNumber(), esd->GetEventNumber(), recStep);
987 if (gSystem->AccessPathName(fileName)) return kFALSE;
988
989 AliDebug(1, Form("reading ESD from file %s", fileName));
990 TFile* file = TFile::Open(fileName);
991 if (!file || !file->IsOpen()) {
992 AliError(Form("opening %s failed", fileName));
993 delete file;
994 return kFALSE;
995 }
996
997 gROOT->cd();
998 delete esd;
999 esd = (AliESD*) file->Get("ESD");
1000 file->Close();
1001 delete file;
1002 return kTRUE;
1003}
1004
1005//_____________________________________________________________________________
1006void AliReconstruction::WriteESD(AliESD* esd, const char* recStep) const
1007{
1008// write the ESD event to a file
1009
1010 if (!esd) return;
1011 char fileName[256];
1012 sprintf(fileName, "ESD_%d.%d_%s.root",
1013 esd->GetRunNumber(), esd->GetEventNumber(), recStep);
1014
1015 AliDebug(1, Form("writing ESD to file %s", fileName));
1016 TFile* file = TFile::Open(fileName, "recreate");
1017 if (!file || !file->IsOpen()) {
1018 AliError(Form("opening %s failed", fileName));
1019 } else {
1020 esd->Write("ESD");
1021 file->Close();
1022 }
1023 delete file;
1024}