1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // class for running the reconstruction //
22 // Clusters and tracks are created for all detectors and all events by //
25 // AliReconstruction rec; //
28 // The Run method returns kTRUE in case of successful execution. //
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 //
33 // rec.SetInput("..."); //
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 //
41 // By default all events are reconstructed. The reconstruction can be //
42 // limited to a range of events by giving the index of the first and the //
43 // last event as an argument to the Run method or by calling //
45 // rec.SetEventRange(..., ...); //
47 // The index -1 (default) can be used for the last event to indicate no //
48 // upper limit of the event range. //
50 // In case of raw-data reconstruction the user can modify the default //
51 // number of events per digits/clusters/tracks file. In case the option //
52 // is not used the number is set 1. In case the user provides 0, than //
53 // the number of events is equal to the number of events inside the //
54 // raw-data file (i.e. one digits/clusters/tracks file): //
56 // rec.SetNumberOfEventsPerFile(...); //
59 // The name of the galice file can be changed from the default //
60 // "galice.root" by passing it as argument to the AliReconstruction //
61 // constructor or by //
63 // rec.SetGAliceFile("..."); //
65 // The local reconstruction can be switched on or off for individual //
68 // rec.SetRunLocalReconstruction("..."); //
70 // The argument is a (case sensitive) string with the names of the //
71 // detectors separated by a space. The special string "ALL" selects all //
72 // available detectors. This is the default. //
74 // The reconstruction of the primary vertex position can be switched off by //
76 // rec.SetRunVertexFinder(kFALSE); //
78 // The tracking and the creation of ESD tracks can be switched on for //
79 // selected detectors by //
81 // rec.SetRunTracking("..."); //
83 // Uniform/nonuniform field tracking switches (default: uniform field) //
85 // rec.SetUniformFieldTracking(); ( rec.SetUniformFieldTracking(kFALSE); ) //
87 // The filling of additional ESD information can be steered by //
89 // rec.SetFillESD("..."); //
91 // Again, for both methods the string specifies the list of detectors. //
92 // The default is "ALL". //
94 // The call of the shortcut method //
96 // rec.SetRunReconstruction("..."); //
98 // is equivalent to calling SetRunLocalReconstruction, SetRunTracking and //
99 // SetFillESD with the same detector selecting string as argument. //
101 // The reconstruction requires digits or raw data as input. For the creation //
102 // of digits and raw data have a look at the class AliSimulation. //
104 // The input data of a detector can be replaced by the corresponding HLT //
105 // data by calling (usual detector string) //
106 // SetUseHLTData("..."); //
108 // For debug purposes the method SetCheckPointLevel can be used. If the //
109 // argument is greater than 0, files with ESD events will be written after //
110 // selected steps of the reconstruction for each event: //
111 // level 1: after tracking and after filling of ESD (final) //
112 // level 2: in addition after each tracking step //
113 // level 3: in addition after the filling of ESD for each detector //
114 // If a final check point file exists for an event, this event will be //
115 // skipped in the reconstruction. The tracking and the filling of ESD for //
116 // a detector will be skipped as well, if the corresponding check point //
117 // file exists. The ESD event will then be loaded from the file instead. //
119 ///////////////////////////////////////////////////////////////////////////////
126 #include <TPluginManager.h>
127 #include <TGeoManager.h>
128 #include <TLorentzVector.h>
131 #include <TObjArray.h>
133 #include "AliReconstruction.h"
134 #include "AliCodeTimer.h"
135 #include "AliReconstructor.h"
137 #include "AliRunLoader.h"
139 #include "AliRawReaderFile.h"
140 #include "AliRawReaderDate.h"
141 #include "AliRawReaderRoot.h"
142 #include "AliRawEventHeaderBase.h"
143 #include "AliESDEvent.h"
144 #include "AliESDMuonTrack.h"
145 #include "AliESDfriend.h"
146 #include "AliESDVertex.h"
147 #include "AliESDcascade.h"
148 #include "AliESDkink.h"
149 #include "AliESDtrack.h"
150 #include "AliESDCaloCluster.h"
151 #include "AliESDCaloCells.h"
152 #include "AliMultiplicity.h"
153 #include "AliTracker.h"
154 #include "AliVertexer.h"
155 #include "AliVertexerTracks.h"
156 #include "AliV0vertexer.h"
157 #include "AliCascadeVertexer.h"
158 #include "AliHeader.h"
159 #include "AliGenEventHeader.h"
161 #include "AliESDpid.h"
162 #include "AliESDtrack.h"
163 #include "AliESDPmdTrack.h"
165 #include "AliESDTagCreator.h"
166 #include "AliAODTagCreator.h"
168 #include "AliGeomManager.h"
169 #include "AliTrackPointArray.h"
170 #include "AliCDBManager.h"
171 #include "AliCDBStorage.h"
172 #include "AliCDBEntry.h"
173 #include "AliAlignObj.h"
175 #include "AliCentralTrigger.h"
176 #include "AliTriggerConfiguration.h"
177 #include "AliTriggerClass.h"
178 #include "AliCTPRawStream.h"
180 #include "AliQADataMakerRec.h"
181 #include "AliGlobalQADataMaker.h"
183 #include "AliQADataMakerSteer.h"
185 #include "AliPlaneEff.h"
187 #include "AliSysInfo.h" // memory snapshots
188 #include "AliRawHLTManager.h"
191 ClassImp(AliReconstruction)
194 //_____________________________________________________________________________
195 const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
197 //_____________________________________________________________________________
198 AliReconstruction::AliReconstruction(const char* gAliceFilename,
199 const char* name, const char* title) :
202 fUniformField(kTRUE),
203 fRunVertexFinder(kTRUE),
204 fRunVertexFinderTracks(kTRUE),
205 fRunHLTTracking(kFALSE),
206 fRunMuonTracking(kFALSE),
208 fRunCascadeFinder(kTRUE),
209 fStopOnError(kFALSE),
210 fWriteAlignmentData(kFALSE),
211 fWriteESDfriend(kFALSE),
213 fFillTriggerESD(kTRUE),
221 fRunLocalReconstruction("ALL"),
224 fUseTrackingErrorsForAlignment(""),
225 fGAliceFileName(gAliceFilename),
230 fNumberOfEventsPerFile(1),
233 fLoadAlignFromCDB(kTRUE),
234 fLoadAlignData("ALL"),
240 fParentRawReader(NULL),
243 fDiamondProfile(NULL),
244 fDiamondProfileTPC(NULL),
245 fMeanVertexConstraint(kTRUE),
249 fAlignObjArray(NULL),
252 fInitCDBCalled(kFALSE),
253 fSetRunNumberFromDataCalled(kFALSE),
257 fSameQACycle(kFALSE),
259 fRunPlaneEff(kFALSE),
271 fIsNewRunLoader(kFALSE)
273 // create reconstruction object with default parameters
275 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
276 fReconstructor[iDet] = NULL;
277 fLoader[iDet] = NULL;
278 fTracker[iDet] = NULL;
279 fQADataMaker[iDet] = NULL;
280 fQACycles[iDet] = 999999;
282 fQADataMaker[fgkNDetectors]=NULL; //Global QA
286 //_____________________________________________________________________________
287 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
290 fUniformField(rec.fUniformField),
291 fRunVertexFinder(rec.fRunVertexFinder),
292 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
293 fRunHLTTracking(rec.fRunHLTTracking),
294 fRunMuonTracking(rec.fRunMuonTracking),
295 fRunV0Finder(rec.fRunV0Finder),
296 fRunCascadeFinder(rec.fRunCascadeFinder),
297 fStopOnError(rec.fStopOnError),
298 fWriteAlignmentData(rec.fWriteAlignmentData),
299 fWriteESDfriend(rec.fWriteESDfriend),
300 fWriteAOD(rec.fWriteAOD),
301 fFillTriggerESD(rec.fFillTriggerESD),
303 fCleanESD(rec.fCleanESD),
304 fV0DCAmax(rec.fV0DCAmax),
305 fV0CsPmin(rec.fV0CsPmin),
309 fRunLocalReconstruction(rec.fRunLocalReconstruction),
310 fRunTracking(rec.fRunTracking),
311 fFillESD(rec.fFillESD),
312 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
313 fGAliceFileName(rec.fGAliceFileName),
315 fEquipIdMap(rec.fEquipIdMap),
316 fFirstEvent(rec.fFirstEvent),
317 fLastEvent(rec.fLastEvent),
318 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
321 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
322 fLoadAlignData(rec.fLoadAlignData),
323 fESDPar(rec.fESDPar),
324 fUseHLTData(rec.fUseHLTData),
328 fParentRawReader(NULL),
331 fDiamondProfile(NULL),
332 fDiamondProfileTPC(NULL),
333 fMeanVertexConstraint(rec.fMeanVertexConstraint),
337 fAlignObjArray(rec.fAlignObjArray),
338 fCDBUri(rec.fCDBUri),
340 fInitCDBCalled(rec.fInitCDBCalled),
341 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
343 fRunGlobalQA(rec.fRunGlobalQA),
344 fInLoopQA(rec.fInLoopQA),
345 fSameQACycle(rec.fSameQACycle),
346 fRunPlaneEff(rec.fRunPlaneEff),
358 fIsNewRunLoader(rec.fIsNewRunLoader)
362 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
363 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
365 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
366 fReconstructor[iDet] = NULL;
367 fLoader[iDet] = NULL;
368 fTracker[iDet] = NULL;
369 fQADataMaker[iDet] = NULL;
370 fQACycles[iDet] = rec.fQACycles[iDet];
372 fQADataMaker[fgkNDetectors]=NULL; //Global QA
373 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
374 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
378 //_____________________________________________________________________________
379 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
381 // assignment operator
383 this->~AliReconstruction();
384 new(this) AliReconstruction(rec);
388 //_____________________________________________________________________________
389 AliReconstruction::~AliReconstruction()
395 fSpecCDBUri.Delete();
397 AliCodeTimer::Instance()->Print();
400 //_____________________________________________________________________________
401 void AliReconstruction::InitCDB()
403 // activate a default CDB storage
404 // First check if we have any CDB storage set, because it is used
405 // to retrieve the calibration and alignment constants
407 if (fInitCDBCalled) return;
408 fInitCDBCalled = kTRUE;
410 AliCDBManager* man = AliCDBManager::Instance();
411 if (man->IsDefaultStorageSet())
413 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
414 AliWarning("Default CDB storage has been already set !");
415 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
416 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
417 fCDBUri = man->GetDefaultStorage()->GetURI();
420 if (fCDBUri.Length() > 0)
422 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
423 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
424 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
426 fCDBUri="local://$ALICE_ROOT";
427 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
428 AliWarning("Default CDB storage not yet set !!!!");
429 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
430 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
433 man->SetDefaultStorage(fCDBUri);
436 // Now activate the detector specific CDB storage locations
437 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
438 TObject* obj = fSpecCDBUri[i];
440 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
441 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
442 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
443 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
448 //_____________________________________________________________________________
449 void AliReconstruction::SetDefaultStorage(const char* uri) {
450 // Store the desired default CDB storage location
451 // Activate it later within the Run() method
457 //_____________________________________________________________________________
458 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
459 // Store a detector-specific CDB storage location
460 // Activate it later within the Run() method
462 AliCDBPath aPath(calibType);
463 if(!aPath.IsValid()){
464 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
465 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
466 if(!strcmp(calibType, fgkDetectorName[iDet])) {
467 aPath.SetPath(Form("%s/*", calibType));
468 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
472 if(!aPath.IsValid()){
473 AliError(Form("Not a valid path or detector: %s", calibType));
478 // // check that calibType refers to a "valid" detector name
479 // Bool_t isDetector = kFALSE;
480 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
481 // TString detName = fgkDetectorName[iDet];
482 // if(aPath.GetLevel0() == detName) {
483 // isDetector = kTRUE;
489 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
493 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
494 if (obj) fSpecCDBUri.Remove(obj);
495 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
499 //_____________________________________________________________________________
500 Bool_t AliReconstruction::SetRunNumberFromData()
502 // The method is called in Run() in order
503 // to set a correct run number.
504 // In case of raw data reconstruction the
505 // run number is taken from the raw data header
507 if (fSetRunNumberFromDataCalled) return kTRUE;
508 fSetRunNumberFromDataCalled = kTRUE;
510 AliCDBManager* man = AliCDBManager::Instance();
512 if(man->GetRun() > 0) {
513 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
517 AliError("No run loader is found !");
520 // read run number from gAlice
521 if(fRunLoader->GetAliRun())
522 AliCDBManager::Instance()->SetRun(fRunLoader->GetHeader()->GetRun());
525 if(fRawReader->NextEvent()) {
526 AliCDBManager::Instance()->SetRun(fRawReader->GetRunNumber());
527 fRawReader->RewindEvents();
530 if(man->GetRun() > 0) {
531 AliWarning("No raw events is found ! Using settings in AliCDBManager !");
536 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
542 AliError("Neither gAlice nor RawReader objects are found !");
552 //_____________________________________________________________________________
553 void AliReconstruction::SetCDBLock() {
554 // Set CDB lock: from now on it is forbidden to reset the run number
555 // or the default storage or to activate any further storage!
557 AliCDBManager::Instance()->SetLock(1);
560 //_____________________________________________________________________________
561 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
563 // Read the alignment objects from CDB.
564 // Each detector is supposed to have the
565 // alignment objects in DET/Align/Data CDB path.
566 // All the detector objects are then collected,
567 // sorted by geometry level (starting from ALIC) and
568 // then applied to the TGeo geometry.
569 // Finally an overlaps check is performed.
571 // Load alignment data from CDB and fill fAlignObjArray
572 if(fLoadAlignFromCDB){
574 TString detStr = detectors;
575 TString loadAlObjsListOfDets = "";
577 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
578 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
579 loadAlObjsListOfDets += fgkDetectorName[iDet];
580 loadAlObjsListOfDets += " ";
581 } // end loop over detectors
582 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
583 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
585 // Check if the array with alignment objects was
586 // provided by the user. If yes, apply the objects
587 // to the present TGeo geometry
588 if (fAlignObjArray) {
589 if (gGeoManager && gGeoManager->IsClosed()) {
590 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
591 AliError("The misalignment of one or more volumes failed!"
592 "Compare the list of simulated detectors and the list of detector alignment data!");
597 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
603 delete fAlignObjArray; fAlignObjArray=0;
608 //_____________________________________________________________________________
609 void AliReconstruction::SetGAliceFile(const char* fileName)
611 // set the name of the galice file
613 fGAliceFileName = fileName;
616 //_____________________________________________________________________________
617 void AliReconstruction::SetInput(const char* input)
619 // In case the input string starts with 'mem://', we run in an online mode
620 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
621 // file is assumed. One can give as an input:
622 // mem://: - events taken from DAQ monitoring libs online
624 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
628 //_____________________________________________________________________________
629 void AliReconstruction::SetOption(const char* detector, const char* option)
631 // set options for the reconstruction of a detector
633 TObject* obj = fOptions.FindObject(detector);
634 if (obj) fOptions.Remove(obj);
635 fOptions.Add(new TNamed(detector, option));
638 //_____________________________________________________________________________
639 Bool_t AliReconstruction::Run(const char* input)
642 AliCodeTimerAuto("");
644 if (!InitRun(input)) return kFALSE;
646 //******* The loop over events
648 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
649 (fRawReader && fRawReader->NextEvent())) {
650 if (!RunEvent(iEvent)) return kFALSE;
654 if (!FinishRun()) return kFALSE;
659 //_____________________________________________________________________________
660 Bool_t AliReconstruction::InitRun(const char* input)
662 // Initialize all the stuff before
663 // going into the event loop
664 // If the second argument is given, the first one is ignored and
665 // the reconstruction works in an online mode
666 AliCodeTimerAuto("");
668 // Overwrite the previous setting
669 if (input) fInput = input;
671 // set the input in case of raw data
672 fRawReader = AliRawReader::Create(fInput.Data());
674 AliInfo("Reconstruction will run over digits");
676 if (!fEquipIdMap.IsNull() && fRawReader)
677 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
679 if (!fUseHLTData.IsNull()) {
680 // create the RawReaderHLT which performs redirection of HLT input data for
681 // the specified detectors
682 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
684 fParentRawReader=fRawReader;
685 fRawReader=pRawReader;
687 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
691 AliSysInfo::AddStamp("Start");
692 // get the run loader
693 if (!InitRunLoader()) return kFALSE;
694 AliSysInfo::AddStamp("LoadLoader");
696 // Initialize the CDB storage
699 AliSysInfo::AddStamp("LoadCDB");
701 // Set run number in CDBManager (if it is not already set by the user)
702 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
704 // Set CDB lock: from now on it is forbidden to reset the run number
705 // or the default storage or to activate any further storage!
708 // Import ideal TGeo geometry and apply misalignment
710 TString geom(gSystem->DirName(fGAliceFileName));
711 geom += "/geometry.root";
712 AliGeomManager::LoadGeometry(geom.Data());
713 if (!gGeoManager) if (fStopOnError) return kFALSE;
716 if (!MisalignGeometry(fLoadAlignData)) if (fStopOnError) return kFALSE;
717 AliSysInfo::AddStamp("LoadGeom");
721 // local reconstruction
722 if (!fRunLocalReconstruction.IsNull()) {
723 if (!RunLocalReconstruction(fRunLocalReconstruction)) {
724 if (fStopOnError) {CleanUp(); return kFALSE;}
730 if (fRunVertexFinder && !CreateVertexer()) {
736 AliSysInfo::AddStamp("Vertexer");
739 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
745 AliSysInfo::AddStamp("LoadTrackers");
747 // get the possibly already existing ESD file and tree
748 fesd = new AliESDEvent(); fhltesd = new AliESDEvent();
749 if (!gSystem->AccessPathName("AliESDs.root")){
750 gSystem->CopyFile("AliESDs.root", "AliESDs.old.root", kTRUE);
751 ffileOld = TFile::Open("AliESDs.old.root");
752 if (ffileOld && ffileOld->IsOpen()) {
753 ftreeOld = (TTree*) ffileOld->Get("esdTree");
754 if (ftreeOld)fesd->ReadFromTree(ftreeOld);
755 fhlttreeOld = (TTree*) ffileOld->Get("HLTesdTree");
756 if (fhlttreeOld) fhltesd->ReadFromTree(fhlttreeOld);
760 // create the ESD output file and tree
761 ffile = TFile::Open("AliESDs.root", "RECREATE");
762 ffile->SetCompressionLevel(2);
763 if (!ffile->IsOpen()) {
764 AliError("opening AliESDs.root failed");
765 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
768 ftree = new TTree("esdTree", "Tree with ESD objects");
769 fesd = new AliESDEvent();
770 fesd->CreateStdContent();
771 fesd->WriteToTree(ftree);
773 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
774 fhltesd = new AliESDEvent();
775 fhltesd->CreateStdContent();
776 fhltesd->WriteToTree(fhlttree);
779 delete esd; delete hltesd;
780 esd = NULL; hltesd = NULL;
782 // create the branch with ESD additions
786 if (fWriteESDfriend) {
787 fesdf = new AliESDfriend();
788 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
789 br->SetFile("AliESDfriends.root");
790 fesd->AddObject(fesdf);
794 // Get the GRP CDB entry
795 AliCDBEntry* entryGRP = AliCDBManager::Instance()->Get("GRP/GRP/Data");
798 fGRPList = dynamic_cast<TList*> (entryGRP->GetObject());
800 AliError("No GRP entry found in OCDB!");
803 // Get the diamond profile from OCDB
804 AliCDBEntry* entry = AliCDBManager::Instance()
805 ->Get("GRP/Calib/MeanVertex");
808 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
810 AliError("No diamond profile found in OCDB!");
814 entry = AliCDBManager::Instance()
815 ->Get("GRP/Calib/MeanVertexTPC");
818 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
820 AliError("No diamond profile found in OCDB!");
823 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
824 if(fDiamondProfile && fMeanVertexConstraint) ftVertexer->SetVtxStart(fDiamondProfile);
826 if (fRawReader) fRawReader->RewindEvents();
829 gSystem->GetProcInfo(&ProcInfo);
830 AliInfo(Form("Current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual));
833 //Initialize the QA and start of cycle for out-of-cycle QA
835 AliQADataMakerSteer qas ;
836 if (fRunQA && fRawReader) {
837 qas.SetEventRange(fFirstEvent, fLastEvent) ;
838 qas.Run(fRunLocalReconstruction, fRawReader) ;
839 fSameQACycle = kTRUE ;
841 // checking the QA of previous steps
844 TString detStr(fFillESD);
845 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
846 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
847 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
849 AliInfo(Form("Initializing the QA data maker for %s",
850 fgkDetectorName[iDet]));
852 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
853 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
855 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
856 qadm->StartOfCycle(AliQA::kESDS,"same");
860 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
861 AliInfo(Form("Initializing the global QA data maker"));
863 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
864 AliTracker::SetResidualsArray(arr);
865 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
867 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
868 qadm->StartOfCycle(AliQA::kESDS, "same");
872 fSameQACycle = kTRUE;
875 //Initialize the Plane Efficiency framework
876 if (fRunPlaneEff && !InitPlaneEff()) {
877 if(fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
883 //_____________________________________________________________________________
884 Bool_t AliReconstruction::RunEvent(Int_t iEvent)
886 // run the reconstruction over a single event
887 // The event loop is steered in Run method
889 AliCodeTimerAuto("");
891 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
892 fRunLoader->SetEventNumber(iEvent);
893 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
895 //?? fRunLoader->MakeTree("H");
896 fRunLoader->TreeE()->Fill();
899 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
900 // copy old ESD to the new one
902 fesd->ReadFromTree(ftreeOld);
903 ftreeOld->GetEntry(iEvent);
907 fesd->ReadFromTree(fhlttreeOld);
908 fhlttreeOld->GetEntry(iEvent);
914 AliInfo(Form("processing event %d", iEvent));
916 //Start of cycle for the in-loop QA
919 TString detStr(fFillESD);
920 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
921 if (!IsSelected(fgkDetectorName[iDet], detStr))
923 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
926 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
927 qadm->StartOfCycle(AliQA::kESDS, "same") ;
930 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
931 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
932 qadm->StartOfCycle(AliQA::kESDS, "same");
937 fRunLoader->GetEvent(iEvent);
940 sprintf(aFileName, "ESD_%d.%d_final.root",
941 fRunLoader->GetHeader()->GetRun(),
942 fRunLoader->GetHeader()->GetEventNrInRun());
943 if (!gSystem->AccessPathName(aFileName)) return kTRUE;
945 // local signle event reconstruction
946 if (!fRunLocalReconstruction.IsNull()) {
947 TString detectors="HLT";
948 // run HLT event reconstruction first
949 if (IsSelected(detectors, fRunLocalReconstruction) &&
950 !RunLocalEventReconstruction(detectors)) {
951 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
953 detectors=fRunLocalReconstruction;
954 detectors.ReplaceAll("HLT", "");
955 if (!RunLocalEventReconstruction(detectors)) {
956 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
960 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
961 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
962 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
963 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
965 // Set magnetic field from the tracker
966 fesd->SetMagneticField(AliTracker::GetBz());
967 fhltesd->SetMagneticField(AliTracker::GetBz());
971 // Fill raw-data error log into the ESD
972 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
975 if (fRunVertexFinder) {
976 if (!ReadESD(fesd, "vertex")) {
977 if (!RunVertexFinder(fesd)) {
978 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
980 if (fCheckPointLevel > 0) WriteESD(fesd, "vertex");
985 if (!fRunTracking.IsNull()) {
986 if (fRunHLTTracking) {
987 fhltesd->SetPrimaryVertexSPD(fesd->GetVertex());
988 if (!RunHLTTracking(fhltesd)) {
989 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
995 if (!fRunTracking.IsNull()) {
996 if (fRunMuonTracking) {
997 if (!RunMuonTracking(fesd)) {
998 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1004 if (!fRunTracking.IsNull()) {
1005 if (!ReadESD(fesd, "tracking")) {
1006 if (!RunTracking(fesd)) {
1007 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1009 if (fCheckPointLevel > 0) WriteESD(fesd, "tracking");
1014 if (!fFillESD.IsNull()) {
1015 if (!FillESD(fesd, fFillESD)) {
1016 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1020 // fill Event header information from the RawEventHeader
1021 if (fRawReader){FillRawEventHeaderESD(fesd);}
1024 AliESDpid::MakePID(fesd);
1025 if (fCheckPointLevel > 1) WriteESD(fesd, "PID");
1027 if (fFillTriggerESD) {
1028 if (!ReadESD(fesd, "trigger")) {
1029 if (!FillTriggerESD(fesd)) {
1030 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1032 if (fCheckPointLevel > 1) WriteESD(fesd, "trigger");
1039 // Propagate track to the beam pipe (if not laready done by ITS)
1041 const Int_t ntracks = fesd->GetNumberOfTracks();
1042 const Double_t kBz = fesd->GetMagneticField();
1043 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1046 UShort_t *selectedIdx=new UShort_t[ntracks];
1048 for (Int_t itrack=0; itrack<ntracks; itrack++){
1049 const Double_t kMaxStep = 5; //max step over the material
1052 AliESDtrack *track = fesd->GetTrack(itrack);
1053 if (!track) continue;
1055 AliExternalTrackParam *tpcTrack =
1056 (AliExternalTrackParam *)track->GetTPCInnerParam();
1060 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kTRUE);
1065 Int_t n=trkArray.GetEntriesFast();
1066 selectedIdx[n]=track->GetID();
1067 trkArray.AddLast(tpcTrack);
1071 //Tracks refitted by ITS should already be at the SPD vertex
1072 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1076 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kTRUE);
1077 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kVeryBig);
1082 // Improve the reconstructed primary vertex position using the tracks
1084 TObject *obj = fOptions.FindObject("ITS");
1086 TString optITS = obj->GetTitle();
1087 if (optITS.Contains("cosmics") || optITS.Contains("COSMICS"))
1088 fRunVertexFinderTracks=kFALSE;
1090 if (fRunVertexFinderTracks) {
1091 // TPC + ITS primary vertex
1092 ftVertexer->SetITSrefitRequired();
1093 if(fDiamondProfile && fMeanVertexConstraint) {
1094 ftVertexer->SetVtxStart(fDiamondProfile);
1096 ftVertexer->SetConstraintOff();
1098 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1100 if (pvtx->GetStatus()) {
1101 fesd->SetPrimaryVertex(pvtx);
1102 for (Int_t i=0; i<ntracks; i++) {
1103 AliESDtrack *t = fesd->GetTrack(i);
1104 t->RelateToVertex(pvtx, kBz, kVeryBig);
1109 // TPC-only primary vertex
1110 ftVertexer->SetITSrefitNotRequired();
1111 if(fDiamondProfileTPC && fMeanVertexConstraint) {
1112 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1114 ftVertexer->SetConstraintOff();
1116 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1118 if (pvtx->GetStatus()) {
1119 fesd->SetPrimaryVertexTPC(pvtx);
1120 Int_t nsel=trkArray.GetEntriesFast();
1121 for (Int_t i=0; i<nsel; i++) {
1122 AliExternalTrackParam *t =
1123 (AliExternalTrackParam *)trkArray.UncheckedAt(i);
1124 t->PropagateToDCA(pvtx, kBz, kRadius);
1130 delete[] selectedIdx;
1132 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1137 AliV0vertexer vtxer;
1138 vtxer.Tracks2V0vertices(fesd);
1140 if (fRunCascadeFinder) {
1142 AliCascadeVertexer cvtxer;
1143 cvtxer.V0sTracks2CascadeVertices(fesd);
1148 if (fCleanESD) CleanESD(fesd);
1152 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1153 if (qadm) qadm->Exec(AliQA::kESDS, fesd);
1157 if (fWriteESDfriend) {
1158 fesdf->~AliESDfriend();
1159 new (fesdf) AliESDfriend(); // Reset...
1160 fesd->GetESDfriend(fesdf);
1167 if (fCheckPointLevel > 0) WriteESD(fesd, "final");
1170 if (fWriteESDfriend) {
1171 fesdf->~AliESDfriend();
1172 new (fesdf) AliESDfriend(); // Reset...
1175 ProcInfo_t ProcInfo;
1176 gSystem->GetProcInfo(&ProcInfo);
1177 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1180 // End of cycle for the in-loop QA
1183 RunQA(fFillESD.Data(), fesd);
1184 TString detStr(fFillESD);
1185 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1186 if (!IsSelected(fgkDetectorName[iDet], detStr))
1188 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1191 qadm->EndOfCycle(AliQA::kRECPOINTS);
1192 qadm->EndOfCycle(AliQA::kESDS);
1197 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1199 qadm->EndOfCycle(AliQA::kRECPOINTS);
1200 qadm->EndOfCycle(AliQA::kESDS);
1209 //_____________________________________________________________________________
1210 Bool_t AliReconstruction::FinishRun()
1213 // Called after the exit
1214 // from the event loop
1215 AliCodeTimerAuto("");
1217 if (fIsNewRunLoader) { // galice.root didn't exist
1218 fRunLoader->WriteHeader("OVERWRITE");
1219 fRunLoader->CdGAFile();
1220 fRunLoader->Write(0, TObject::kOverwrite);
1223 ftree->GetUserInfo()->Add(fesd);
1224 fhlttree->GetUserInfo()->Add(fhltesd);
1226 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1227 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1229 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1230 cdbMapCopy->SetOwner(1);
1231 cdbMapCopy->SetName("cdbMap");
1232 TIter iter(cdbMap->GetTable());
1235 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1236 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1237 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1238 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1241 TList *cdbListCopy = new TList();
1242 cdbListCopy->SetOwner(1);
1243 cdbListCopy->SetName("cdbList");
1245 TIter iter2(cdbList);
1248 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1249 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1252 ftree->GetUserInfo()->Add(cdbMapCopy);
1253 ftree->GetUserInfo()->Add(cdbListCopy);
1256 if(fESDPar.Contains("ESD.par")){
1257 AliInfo("Attaching ESD.par to Tree");
1258 TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1259 ftree->GetUserInfo()->Add(fn);
1265 if (fWriteESDfriend)
1266 ftree->SetBranchStatus("ESDfriend*",0);
1267 // we want to have only one tree version number
1268 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1271 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1272 if (fRunPlaneEff && !FinishPlaneEff()) {
1273 AliWarning("Finish PlaneEff evaluation failed");
1277 CleanUp(ffile, ffileOld);
1280 AliWarning("AOD creation not supported anymore during reconstruction. See ANALYSIS/AliAnalysisTaskESDfilter.cxx instead.");
1283 // Create tags for the events in the ESD tree (the ESD tree is always present)
1284 // In case of empty events the tags will contain dummy values
1285 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1286 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPList);
1288 AliWarning("AOD tag creation not supported anymore during reconstruction.");
1291 //Finish QA and end of cycle for out-of-loop QA
1294 AliQADataMakerSteer qas;
1295 qas.SetEventRange(fFirstEvent, fLastEvent) ;
1296 qas.Run(fRunLocalReconstruction.Data(), AliQA::kRECPOINTS, fSameQACycle);
1298 qas.Run(fRunTracking.Data(), AliQA::kESDS, fSameQACycle);
1300 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1302 qadm->EndOfCycle(AliQA::kRECPOINTS);
1303 qadm->EndOfCycle(AliQA::kESDS);
1310 // Cleanup of CDB manager: cache and active storages!
1311 AliCDBManager::Instance()->ClearCache();
1317 //_____________________________________________________________________________
1318 Bool_t AliReconstruction::RunLocalReconstruction(const TString& /*detectors*/)
1320 // run the local reconstruction
1321 static Int_t eventNr=0;
1322 AliCodeTimerAuto("")
1324 // AliCDBManager* man = AliCDBManager::Instance();
1325 // Bool_t origCache = man->GetCacheFlag();
1327 // TString detStr = detectors;
1328 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1329 // if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1330 // AliReconstructor* reconstructor = GetReconstructor(iDet);
1331 // if (!reconstructor) continue;
1332 // if (reconstructor->HasLocalReconstruction()) continue;
1334 // AliCodeTimerStart(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1335 // AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1337 // AliCodeTimerStart(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1338 // AliInfo(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1340 // man->SetCacheFlag(kTRUE);
1341 // TString calibPath = Form("%s/Calib/*", fgkDetectorName[iDet]);
1342 // man->GetAll(calibPath); // entries are cached!
1344 // AliCodeTimerStop(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1346 // if (fRawReader) {
1347 // fRawReader->RewindEvents();
1348 // reconstructor->Reconstruct(fRunLoader, fRawReader);
1350 // reconstructor->Reconstruct(fRunLoader);
1353 // AliCodeTimerStop(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1354 // AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr));
1356 // // unload calibration data
1357 // man->UnloadFromCache(calibPath);
1358 // //man->ClearCache();
1361 // man->SetCacheFlag(origCache);
1363 // if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1364 // AliError(Form("the following detectors were not found: %s",
1366 // if (fStopOnError) return kFALSE;
1373 //_____________________________________________________________________________
1374 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1376 // run the local reconstruction
1378 static Int_t eventNr=0;
1379 AliCodeTimerAuto("")
1381 TString detStr = detectors;
1382 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1383 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1384 AliReconstructor* reconstructor = GetReconstructor(iDet);
1385 if (!reconstructor) continue;
1386 AliLoader* loader = fLoader[iDet];
1388 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1391 // conversion of digits
1392 if (fRawReader && reconstructor->HasDigitConversion()) {
1393 AliInfo(Form("converting raw data digits into root objects for %s",
1394 fgkDetectorName[iDet]));
1395 AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1396 fgkDetectorName[iDet]));
1397 loader->LoadDigits("update");
1398 loader->CleanDigits();
1399 loader->MakeDigitsContainer();
1400 TTree* digitsTree = loader->TreeD();
1401 reconstructor->ConvertDigits(fRawReader, digitsTree);
1402 loader->WriteDigits("OVERWRITE");
1403 loader->UnloadDigits();
1405 // local reconstruction
1406 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1407 AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1408 loader->LoadRecPoints("update");
1409 loader->CleanRecPoints();
1410 loader->MakeRecPointsContainer();
1411 TTree* clustersTree = loader->TreeR();
1412 if (fRawReader && !reconstructor->HasDigitConversion()) {
1413 reconstructor->Reconstruct(fRawReader, clustersTree);
1415 loader->LoadDigits("read");
1416 TTree* digitsTree = loader->TreeD();
1418 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1419 if (fStopOnError) return kFALSE;
1421 reconstructor->Reconstruct(digitsTree, clustersTree);
1423 loader->UnloadDigits();
1426 // In-loop QA for local reconstrucion
1427 if (fRunQA && fInLoopQA) {
1428 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1431 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1433 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1435 qadm->Exec(AliQA::kRECPOINTS, clustersTree) ;
1438 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1442 loader->WriteRecPoints("OVERWRITE");
1443 loader->UnloadRecPoints();
1444 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1447 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1448 AliError(Form("the following detectors were not found: %s",
1450 if (fStopOnError) return kFALSE;
1456 //_____________________________________________________________________________
1457 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1459 // run the barrel tracking
1461 AliCodeTimerAuto("")
1463 AliESDVertex* vertex = NULL;
1464 Double_t vtxPos[3] = {0, 0, 0};
1465 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1466 TArrayF mcVertex(3);
1467 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1468 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1469 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1473 if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1474 AliInfo("running the ITS vertex finder");
1475 if (fLoader[0]) fLoader[0]->LoadRecPoints();
1476 vertex = fVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber());
1477 if (fLoader[0]) fLoader[0]->UnloadRecPoints();
1479 AliWarning("Vertex not found");
1480 vertex = new AliESDVertex();
1481 vertex->SetName("default");
1484 vertex->SetName("reconstructed");
1488 AliInfo("getting the primary vertex from MC");
1489 vertex = new AliESDVertex(vtxPos, vtxErr);
1493 vertex->GetXYZ(vtxPos);
1494 vertex->GetSigmaXYZ(vtxErr);
1496 AliWarning("no vertex reconstructed");
1497 vertex = new AliESDVertex(vtxPos, vtxErr);
1499 esd->SetPrimaryVertexSPD(vertex);
1500 // if SPD multiplicity has been determined, it is stored in the ESD
1501 AliMultiplicity *mult = fVertexer->GetMultiplicity();
1502 if(mult)esd->SetMultiplicity(mult);
1504 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1505 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1512 //_____________________________________________________________________________
1513 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
1515 // run the HLT barrel tracking
1517 AliCodeTimerAuto("")
1520 AliError("Missing runLoader!");
1524 AliInfo("running HLT tracking");
1526 // Get a pointer to the HLT reconstructor
1527 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
1528 if (!reconstructor) return kFALSE;
1531 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1532 TString detName = fgkDetectorName[iDet];
1533 AliDebug(1, Form("%s HLT tracking", detName.Data()));
1534 reconstructor->SetOption(detName.Data());
1535 AliTracker *tracker = reconstructor->CreateTracker();
1537 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
1538 if (fStopOnError) return kFALSE;
1542 Double_t vtxErr[3]={0.005,0.005,0.010};
1543 const AliESDVertex *vertex = esd->GetVertex();
1544 vertex->GetXYZ(vtxPos);
1545 tracker->SetVertex(vtxPos,vtxErr);
1547 fLoader[iDet]->LoadRecPoints("read");
1548 TTree* tree = fLoader[iDet]->TreeR();
1550 AliError(Form("Can't get the %s cluster tree", detName.Data()));
1553 tracker->LoadClusters(tree);
1555 if (tracker->Clusters2Tracks(esd) != 0) {
1556 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
1560 tracker->UnloadClusters();
1568 //_____________________________________________________________________________
1569 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
1571 // run the muon spectrometer tracking
1573 AliCodeTimerAuto("")
1576 AliError("Missing runLoader!");
1579 Int_t iDet = 7; // for MUON
1581 AliInfo("is running...");
1583 // Get a pointer to the MUON reconstructor
1584 AliReconstructor *reconstructor = GetReconstructor(iDet);
1585 if (!reconstructor) return kFALSE;
1588 TString detName = fgkDetectorName[iDet];
1589 AliDebug(1, Form("%s tracking", detName.Data()));
1590 AliTracker *tracker = reconstructor->CreateTracker();
1592 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1597 fLoader[iDet]->LoadRecPoints("read");
1599 tracker->LoadClusters(fLoader[iDet]->TreeR());
1601 Int_t rv = tracker->Clusters2Tracks(esd);
1605 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1609 fLoader[iDet]->UnloadRecPoints();
1611 tracker->UnloadClusters();
1619 //_____________________________________________________________________________
1620 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
1622 // run the barrel tracking
1623 static Int_t eventNr=0;
1624 AliCodeTimerAuto("")
1626 AliInfo("running tracking");
1628 //Fill the ESD with the T0 info (will be used by the TOF)
1629 if (fReconstructor[11] && fLoader[11]) {
1630 fLoader[11]->LoadRecPoints("READ");
1631 TTree *treeR = fLoader[11]->TreeR();
1632 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
1635 // pass 1: TPC + ITS inwards
1636 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1637 if (!fTracker[iDet]) continue;
1638 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
1641 fLoader[iDet]->LoadRecPoints("read");
1642 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
1643 TTree* tree = fLoader[iDet]->TreeR();
1645 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1648 fTracker[iDet]->LoadClusters(tree);
1649 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1651 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
1652 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1655 if (fCheckPointLevel > 1) {
1656 WriteESD(esd, Form("%s.tracking", fgkDetectorName[iDet]));
1658 // preliminary PID in TPC needed by the ITS tracker
1660 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1661 AliESDpid::MakePID(esd);
1663 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
1666 // pass 2: ALL backwards
1668 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1669 if (!fTracker[iDet]) continue;
1670 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
1673 if (iDet > 1) { // all except ITS, TPC
1675 fLoader[iDet]->LoadRecPoints("read");
1676 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
1677 tree = fLoader[iDet]->TreeR();
1679 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1682 fTracker[iDet]->LoadClusters(tree);
1683 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1687 if (iDet>1) // start filling residuals for the "outer" detectors
1688 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1690 if (fTracker[iDet]->PropagateBack(esd) != 0) {
1691 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
1694 if (fCheckPointLevel > 1) {
1695 WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
1699 if (iDet > 2) { // all except ITS, TPC, TRD
1700 fTracker[iDet]->UnloadClusters();
1701 fLoader[iDet]->UnloadRecPoints();
1703 // updated PID in TPC needed by the ITS tracker -MI
1705 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1706 AliESDpid::MakePID(esd);
1708 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1710 //stop filling residuals for the "outer" detectors
1711 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1713 // write space-points to the ESD in case alignment data output
1715 if (fWriteAlignmentData)
1716 WriteAlignmentData(esd);
1718 // pass 3: TRD + TPC + ITS refit inwards
1720 for (Int_t iDet = 2; iDet >= 0; iDet--) {
1721 if (!fTracker[iDet]) continue;
1722 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
1725 if (iDet<2) // start filling residuals for TPC and ITS
1726 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1728 if (fTracker[iDet]->RefitInward(esd) != 0) {
1729 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
1732 // run postprocessing
1733 if (fTracker[iDet]->PostProcess(esd) != 0) {
1734 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
1737 if (fCheckPointLevel > 1) {
1738 WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
1740 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1742 fTracker[iDet]->UnloadClusters();
1743 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
1744 fLoader[iDet]->UnloadRecPoints();
1745 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
1747 // stop filling residuals for TPC and ITS
1748 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1754 //_____________________________________________________________________________
1755 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
1757 // Remove the data which are not needed for the physics analysis.
1760 Int_t nTracks=esd->GetNumberOfTracks();
1761 Int_t nV0s=esd->GetNumberOfV0s();
1763 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
1765 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
1766 Bool_t rc=esd->Clean(cleanPars);
1768 nTracks=esd->GetNumberOfTracks();
1769 nV0s=esd->GetNumberOfV0s();
1771 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
1776 //_____________________________________________________________________________
1777 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
1779 // fill the event summary data
1781 AliCodeTimerAuto("")
1782 static Int_t eventNr=0;
1783 TString detStr = detectors;
1785 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1786 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1787 AliReconstructor* reconstructor = GetReconstructor(iDet);
1788 if (!reconstructor) continue;
1789 if (!ReadESD(esd, fgkDetectorName[iDet])) {
1790 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
1791 TTree* clustersTree = NULL;
1792 if (fLoader[iDet]) {
1793 fLoader[iDet]->LoadRecPoints("read");
1794 clustersTree = fLoader[iDet]->TreeR();
1795 if (!clustersTree) {
1796 AliError(Form("Can't get the %s clusters tree",
1797 fgkDetectorName[iDet]));
1798 if (fStopOnError) return kFALSE;
1801 if (fRawReader && !reconstructor->HasDigitConversion()) {
1802 reconstructor->FillESD(fRawReader, clustersTree, esd);
1804 TTree* digitsTree = NULL;
1805 if (fLoader[iDet]) {
1806 fLoader[iDet]->LoadDigits("read");
1807 digitsTree = fLoader[iDet]->TreeD();
1809 AliError(Form("Can't get the %s digits tree",
1810 fgkDetectorName[iDet]));
1811 if (fStopOnError) return kFALSE;
1814 reconstructor->FillESD(digitsTree, clustersTree, esd);
1815 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
1817 if (fLoader[iDet]) {
1818 fLoader[iDet]->UnloadRecPoints();
1821 if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
1825 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1826 AliError(Form("the following detectors were not found: %s",
1828 if (fStopOnError) return kFALSE;
1830 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
1835 //_____________________________________________________________________________
1836 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
1838 // Reads the trigger decision which is
1839 // stored in Trigger.root file and fills
1840 // the corresponding esd entries
1842 AliCodeTimerAuto("")
1844 AliInfo("Filling trigger information into the ESD");
1846 AliCentralTrigger *aCTP = NULL;
1849 AliCTPRawStream input(fRawReader);
1850 if (!input.Next()) {
1851 AliError("No valid CTP (trigger) DDL raw data is found ! The trigger information is not stored in the ESD !");
1854 esd->SetTriggerMask(input.GetClassMask());
1855 esd->SetTriggerCluster(input.GetClusterMask());
1857 aCTP = new AliCentralTrigger();
1858 TString configstr("");
1859 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
1860 AliError("No trigger configuration found in OCDB! The trigger classes information will no be stored in ESD!");
1866 AliRunLoader *runloader = AliRunLoader::GetRunLoader();
1868 if (!runloader->LoadTrigger()) {
1869 aCTP = runloader->GetTrigger();
1870 esd->SetTriggerMask(aCTP->GetClassMask());
1871 esd->SetTriggerCluster(aCTP->GetClusterMask());
1874 AliWarning("No trigger can be loaded! The trigger information is not stored in the ESD !");
1879 AliError("No run loader is available! The trigger information is not stored in the ESD !");
1884 // Now fill the trigger class names into AliESDRun object
1885 AliTriggerConfiguration *config = aCTP->GetConfiguration();
1887 AliError("No trigger configuration has been found! The trigger classes information will no be stored in ESD!");
1888 if (fRawReader) delete aCTP;
1892 const TObjArray& classesArray = config->GetClasses();
1893 Int_t nclasses = classesArray.GetEntriesFast();
1894 for( Int_t j=0; j<nclasses; j++ ) {
1895 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At( j );
1896 Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
1897 esd->SetTriggerClass(trclass->GetName(),trindex);
1900 if (fRawReader) delete aCTP;
1908 //_____________________________________________________________________________
1909 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
1912 // Filling information from RawReader Header
1915 AliInfo("Filling information from RawReader Header");
1916 esd->SetBunchCrossNumber(0);
1917 esd->SetOrbitNumber(0);
1918 esd->SetPeriodNumber(0);
1919 esd->SetTimeStamp(0);
1920 esd->SetEventType(0);
1921 const AliRawEventHeaderBase * eventHeader = fRawReader->GetEventHeader();
1924 const UInt_t *id = eventHeader->GetP("Id");
1925 esd->SetBunchCrossNumber((id)[1]&0x00000fff);
1926 esd->SetOrbitNumber((((id)[0]<<20)&0xf00000)|(((id)[1]>>12)&0xfffff));
1927 esd->SetPeriodNumber(((id)[0]>>4)&0x0fffffff);
1929 esd->SetTimeStamp((eventHeader->Get("Timestamp")));
1930 esd->SetEventType((eventHeader->Get("Type")));
1937 //_____________________________________________________________________________
1938 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
1940 // check whether detName is contained in detectors
1941 // if yes, it is removed from detectors
1943 // check if all detectors are selected
1944 if ((detectors.CompareTo("ALL") == 0) ||
1945 detectors.BeginsWith("ALL ") ||
1946 detectors.EndsWith(" ALL") ||
1947 detectors.Contains(" ALL ")) {
1952 // search for the given detector
1953 Bool_t result = kFALSE;
1954 if ((detectors.CompareTo(detName) == 0) ||
1955 detectors.BeginsWith(detName+" ") ||
1956 detectors.EndsWith(" "+detName) ||
1957 detectors.Contains(" "+detName+" ")) {
1958 detectors.ReplaceAll(detName, "");
1962 // clean up the detectors string
1963 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1964 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1965 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1970 //_____________________________________________________________________________
1971 Bool_t AliReconstruction::InitRunLoader()
1973 // get or create the run loader
1975 if (gAlice) delete gAlice;
1978 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
1979 // load all base libraries to get the loader classes
1980 TString libs = gSystem->GetLibraries();
1981 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1982 TString detName = fgkDetectorName[iDet];
1983 if (detName == "HLT") continue;
1984 if (libs.Contains("lib" + detName + "base.so")) continue;
1985 gSystem->Load("lib" + detName + "base.so");
1987 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
1989 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1993 fRunLoader->CdGAFile();
1994 if (gFile->GetKey(AliRunLoader::GetGAliceName())) {
1995 if (fRunLoader->LoadgAlice() == 0) {
1996 gAlice = fRunLoader->GetAliRun();
1997 AliTracker::SetFieldMap(gAlice->Field(),fUniformField);
2000 if (!gAlice && !fRawReader) {
2001 AliError(Form("no gAlice object found in file %s",
2002 fGAliceFileName.Data()));
2007 //PH This is a temporary fix to give access to the kinematics
2008 //PH that is needed for the labels of ITS clusters
2009 fRunLoader->LoadHeader();
2010 fRunLoader->LoadKinematics();
2012 } else { // galice.root does not exist
2014 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2018 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2019 AliConfig::GetDefaultEventFolderName(),
2022 AliError(Form("could not create run loader in file %s",
2023 fGAliceFileName.Data()));
2027 fIsNewRunLoader = kTRUE;
2028 fRunLoader->MakeTree("E");
2030 if (fNumberOfEventsPerFile > 0)
2031 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2033 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2039 //_____________________________________________________________________________
2040 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2042 // get the reconstructor object and the loader for a detector
2044 if (fReconstructor[iDet]) return fReconstructor[iDet];
2046 // load the reconstructor object
2047 TPluginManager* pluginManager = gROOT->GetPluginManager();
2048 TString detName = fgkDetectorName[iDet];
2049 TString recName = "Ali" + detName + "Reconstructor";
2050 if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT")) return NULL;
2052 AliReconstructor* reconstructor = NULL;
2053 // first check if a plugin is defined for the reconstructor
2054 TPluginHandler* pluginHandler =
2055 pluginManager->FindHandler("AliReconstructor", detName);
2056 // if not, add a plugin for it
2057 if (!pluginHandler) {
2058 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2059 TString libs = gSystem->GetLibraries();
2060 if (libs.Contains("lib" + detName + "base.so") ||
2061 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2062 pluginManager->AddHandler("AliReconstructor", detName,
2063 recName, detName + "rec", recName + "()");
2065 pluginManager->AddHandler("AliReconstructor", detName,
2066 recName, detName, recName + "()");
2068 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2070 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2071 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2073 if (reconstructor) {
2074 TObject* obj = fOptions.FindObject(detName.Data());
2075 if (obj) reconstructor->SetOption(obj->GetTitle());
2076 reconstructor->Init();
2077 fReconstructor[iDet] = reconstructor;
2080 // get or create the loader
2081 if (detName != "HLT") {
2082 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2083 if (!fLoader[iDet]) {
2084 AliConfig::Instance()
2085 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2087 // first check if a plugin is defined for the loader
2089 pluginManager->FindHandler("AliLoader", detName);
2090 // if not, add a plugin for it
2091 if (!pluginHandler) {
2092 TString loaderName = "Ali" + detName + "Loader";
2093 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2094 pluginManager->AddHandler("AliLoader", detName,
2095 loaderName, detName + "base",
2096 loaderName + "(const char*, TFolder*)");
2097 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2099 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2101 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2102 fRunLoader->GetEventFolder());
2104 if (!fLoader[iDet]) { // use default loader
2105 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2107 if (!fLoader[iDet]) {
2108 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2109 if (fStopOnError) return NULL;
2111 fRunLoader->AddLoader(fLoader[iDet]);
2112 fRunLoader->CdGAFile();
2113 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2114 fRunLoader->Write(0, TObject::kOverwrite);
2119 return reconstructor;
2122 //_____________________________________________________________________________
2123 Bool_t AliReconstruction::CreateVertexer()
2125 // create the vertexer
2128 AliReconstructor* itsReconstructor = GetReconstructor(0);
2129 if (itsReconstructor) {
2130 fVertexer = itsReconstructor->CreateVertexer();
2133 AliWarning("couldn't create a vertexer for ITS");
2134 if (fStopOnError) return kFALSE;
2140 //_____________________________________________________________________________
2141 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2143 // create the trackers
2145 TString detStr = detectors;
2146 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2147 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2148 AliReconstructor* reconstructor = GetReconstructor(iDet);
2149 if (!reconstructor) continue;
2150 TString detName = fgkDetectorName[iDet];
2151 if (detName == "HLT") {
2152 fRunHLTTracking = kTRUE;
2155 if (detName == "MUON") {
2156 fRunMuonTracking = kTRUE;
2161 fTracker[iDet] = reconstructor->CreateTracker();
2162 if (!fTracker[iDet] && (iDet < 7)) {
2163 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2164 if (fStopOnError) return kFALSE;
2166 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2172 //_____________________________________________________________________________
2173 void AliReconstruction::CleanUp(TFile* file, TFile* fileOld)
2175 // delete trackers and the run loader and close and delete the file
2177 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2178 delete fReconstructor[iDet];
2179 fReconstructor[iDet] = NULL;
2180 fLoader[iDet] = NULL;
2181 delete fTracker[iDet];
2182 fTracker[iDet] = NULL;
2183 // delete fQADataMaker[iDet];
2184 // fQADataMaker[iDet] = NULL;
2189 if (ftVertexer) delete ftVertexer;
2192 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2193 delete fDiamondProfile;
2194 fDiamondProfile = NULL;
2195 delete fDiamondProfileTPC;
2196 fDiamondProfileTPC = NULL;
2206 if (fParentRawReader) delete fParentRawReader;
2207 fParentRawReader=NULL;
2217 gSystem->Unlink("AliESDs.old.root");
2222 //_____________________________________________________________________________
2224 Bool_t AliReconstruction::ReadESD(AliESDEvent*& esd, const char* recStep) const
2226 // read the ESD event from a file
2228 if (!esd) return kFALSE;
2230 sprintf(fileName, "ESD_%d.%d_%s.root",
2231 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2232 if (gSystem->AccessPathName(fileName)) return kFALSE;
2234 AliInfo(Form("reading ESD from file %s", fileName));
2235 AliDebug(1, Form("reading ESD from file %s", fileName));
2236 TFile* file = TFile::Open(fileName);
2237 if (!file || !file->IsOpen()) {
2238 AliError(Form("opening %s failed", fileName));
2245 esd = (AliESDEvent*) file->Get("ESD");
2254 //_____________________________________________________________________________
2255 void AliReconstruction::WriteESD(AliESDEvent* esd, const char* recStep) const
2257 // write the ESD event to a file
2261 sprintf(fileName, "ESD_%d.%d_%s.root",
2262 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2264 AliDebug(1, Form("writing ESD to file %s", fileName));
2265 TFile* file = TFile::Open(fileName, "recreate");
2266 if (!file || !file->IsOpen()) {
2267 AliError(Form("opening %s failed", fileName));
2276 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2278 // Write space-points which are then used in the alignment procedures
2279 // For the moment only ITS, TRD and TPC
2281 // Load TOF clusters
2283 fLoader[3]->LoadRecPoints("read");
2284 TTree* tree = fLoader[3]->TreeR();
2286 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[3]));
2289 fTracker[3]->LoadClusters(tree);
2291 Int_t ntracks = esd->GetNumberOfTracks();
2292 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2294 AliESDtrack *track = esd->GetTrack(itrack);
2297 for (Int_t iDet = 3; iDet >= 0; iDet--)
2298 nsp += track->GetNcls(iDet);
2300 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2301 track->SetTrackPointArray(sp);
2303 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2304 AliTracker *tracker = fTracker[iDet];
2305 if (!tracker) continue;
2306 Int_t nspdet = track->GetNcls(iDet);
2307 if (nspdet <= 0) continue;
2308 track->GetClusters(iDet,idx);
2312 while (isp2 < nspdet) {
2314 TString dets = fgkDetectorName[iDet];
2315 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2316 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2317 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2318 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2319 isvalid = tracker->GetTrackPointTrackingError(idx[isp2],p,track);
2321 isvalid = tracker->GetTrackPoint(idx[isp2],p);
2324 const Int_t kNTPCmax = 159;
2325 if (iDet==1 && isp2>kNTPCmax) break; // to be fixed
2326 if (!isvalid) continue;
2327 sp->AddPoint(isptrack,&p); isptrack++; isp++;
2333 fTracker[3]->UnloadClusters();
2334 fLoader[3]->UnloadRecPoints();
2338 //_____________________________________________________________________________
2339 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2341 // The method reads the raw-data error log
2342 // accumulated within the rawReader.
2343 // It extracts the raw-data errors related to
2344 // the current event and stores them into
2345 // a TClonesArray inside the esd object.
2347 if (!fRawReader) return;
2349 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2351 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2353 if (iEvent != log->GetEventNumber()) continue;
2355 esd->AddRawDataErrorLog(log);
2360 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString fName){
2361 // Dump a file content into a char in TNamed
2363 in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2364 Int_t kBytes = (Int_t)in.tellg();
2365 printf("Size: %d \n",kBytes);
2368 char* memblock = new char [kBytes];
2369 in.seekg (0, ios::beg);
2370 in.read (memblock, kBytes);
2372 TString fData(memblock,kBytes);
2373 fn = new TNamed(fName,fData);
2374 printf("fData Size: %d \n",fData.Sizeof());
2375 printf("fName Size: %d \n",fName.Sizeof());
2376 printf("fn Size: %d \n",fn->Sizeof());
2380 AliInfo(Form("Could not Open %s\n",fPath.Data()));
2386 void AliReconstruction::TNamedToFile(TTree* fTree, TString fName){
2387 // This is not really needed in AliReconstruction at the moment
2388 // but can serve as a template
2390 TList *fList = fTree->GetUserInfo();
2391 TNamed *fn = (TNamed*)fList->FindObject(fName.Data());
2392 printf("fn Size: %d \n",fn->Sizeof());
2394 TString fTmp(fn->GetName()); // to be 100% sure in principle fName also works
2395 const char* cdata = fn->GetTitle();
2396 printf("fTmp Size %d\n",fTmp.Sizeof());
2398 int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2399 printf("calculated size %d\n",size);
2400 ofstream out(fName.Data(),ios::out | ios::binary);
2401 out.write(cdata,size);
2406 //_____________________________________________________________________________
2407 AliQADataMakerRec * AliReconstruction::GetQADataMaker(Int_t iDet)
2409 // get the quality assurance data maker object and the loader for a detector
2411 if (fQADataMaker[iDet])
2412 return fQADataMaker[iDet];
2414 AliQADataMakerRec * qadm = NULL;
2415 if (iDet == fgkNDetectors) { //Global QA
2416 qadm = new AliGlobalQADataMaker();
2417 fQADataMaker[iDet] = qadm;
2421 // load the QA data maker object
2422 TPluginManager* pluginManager = gROOT->GetPluginManager();
2423 TString detName = fgkDetectorName[iDet];
2424 TString qadmName = "Ali" + detName + "QADataMakerRec";
2425 if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT"))
2428 // first check if a plugin is defined for the quality assurance data maker
2429 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2430 // if not, add a plugin for it
2431 if (!pluginHandler) {
2432 AliDebug(1, Form("defining plugin for %s", qadmName.Data()));
2433 TString libs = gSystem->GetLibraries();
2434 if (libs.Contains("lib" + detName + "base.so") ||
2435 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2436 pluginManager->AddHandler("AliQADataMakerRec", detName,
2437 qadmName, detName + "qadm", qadmName + "()");
2439 pluginManager->AddHandler("AliQADataMakerRec", detName,
2440 qadmName, detName, qadmName + "()");
2442 pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2444 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2445 qadm = (AliQADataMakerRec *) pluginHandler->ExecPlugin(0);
2448 fQADataMaker[iDet] = qadm;
2453 //_____________________________________________________________________________
2454 Bool_t AliReconstruction::RunQA(const char* detectors, AliESDEvent *& esd)
2456 // run the Quality Assurance data producer
2458 AliCodeTimerAuto("")
2459 TString detStr = detectors;
2460 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2461 if (!IsSelected(fgkDetectorName[iDet], detStr))
2463 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
2466 AliCodeTimerStart(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2467 AliInfo(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2469 qadm->Exec(AliQA::kESDS, esd) ;
2472 AliCodeTimerStop(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2474 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2475 AliError(Form("the following detectors were not found: %s",
2485 //_____________________________________________________________________________
2486 void AliReconstruction::CheckQA()
2488 // check the QA of SIM for this run and remove the detectors
2489 // with status Fatal
2491 TString newRunLocalReconstruction ;
2492 TString newRunTracking ;
2493 TString newFillESD ;
2495 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
2496 TString detName(AliQA::GetDetName(iDet)) ;
2497 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX_t(iDet)) ;
2498 if ( qa->IsSet(AliQA::DETECTORINDEX_t(iDet), AliQA::kSIM, AliQA::kFATAL)) {
2499 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
2501 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
2502 fRunLocalReconstruction.Contains("ALL") ) {
2503 newRunLocalReconstruction += detName ;
2504 newRunLocalReconstruction += " " ;
2506 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
2507 fRunTracking.Contains("ALL") ) {
2508 newRunTracking += detName ;
2509 newRunTracking += " " ;
2511 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
2512 fFillESD.Contains("ALL") ) {
2513 newFillESD += detName ;
2518 fRunLocalReconstruction = newRunLocalReconstruction ;
2519 fRunTracking = newRunTracking ;
2520 fFillESD = newFillESD ;
2523 //_____________________________________________________________________________
2524 Int_t AliReconstruction::GetDetIndex(const char* detector)
2526 // return the detector index corresponding to detector
2528 for (index = 0; index < fgkNDetectors ; index++) {
2529 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2534 //_____________________________________________________________________________
2535 Bool_t AliReconstruction::FinishPlaneEff() {
2537 // Here execute all the necessary operationis, at the end of the tracking phase,
2538 // in case that evaluation of PlaneEfficiencies was required for some detector.
2539 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
2541 // This Preliminary version works only FOR ITS !!!!!
2542 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2545 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2548 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2549 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
2550 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2551 if(fTracker[iDet]) {
2552 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
2553 ret=planeeff->WriteIntoCDB();
2554 if(planeeff->GetCreateHistos()) {
2555 TString name="PlaneEffHisto";
2556 name+=fgkDetectorName[iDet];
2558 ret*=planeeff->WriteHistosToFile(name,"RECREATE");
2564 //_____________________________________________________________________________
2565 Bool_t AliReconstruction::InitPlaneEff() {
2567 // Here execute all the necessary operations, before of the tracking phase,
2568 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
2569 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
2570 // which should be updated/recalculated.
2572 // This Preliminary version will work only FOR ITS !!!!!
2573 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2576 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2578 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));