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),
260 // create reconstruction object with default parameters
262 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
263 fReconstructor[iDet] = NULL;
264 fLoader[iDet] = NULL;
265 fTracker[iDet] = NULL;
266 fQADataMaker[iDet] = NULL;
267 fQACycles[iDet] = 999999;
269 fQADataMaker[fgkNDetectors]=NULL; //Global QA
273 //_____________________________________________________________________________
274 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
277 fUniformField(rec.fUniformField),
278 fRunVertexFinder(rec.fRunVertexFinder),
279 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
280 fRunHLTTracking(rec.fRunHLTTracking),
281 fRunMuonTracking(rec.fRunMuonTracking),
282 fRunV0Finder(rec.fRunV0Finder),
283 fRunCascadeFinder(rec.fRunCascadeFinder),
284 fStopOnError(rec.fStopOnError),
285 fWriteAlignmentData(rec.fWriteAlignmentData),
286 fWriteESDfriend(rec.fWriteESDfriend),
287 fWriteAOD(rec.fWriteAOD),
288 fFillTriggerESD(rec.fFillTriggerESD),
290 fCleanESD(rec.fCleanESD),
291 fV0DCAmax(rec.fV0DCAmax),
292 fV0CsPmin(rec.fV0CsPmin),
296 fRunLocalReconstruction(rec.fRunLocalReconstruction),
297 fRunTracking(rec.fRunTracking),
298 fFillESD(rec.fFillESD),
299 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
300 fGAliceFileName(rec.fGAliceFileName),
302 fEquipIdMap(rec.fEquipIdMap),
303 fFirstEvent(rec.fFirstEvent),
304 fLastEvent(rec.fLastEvent),
305 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
308 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
309 fLoadAlignData(rec.fLoadAlignData),
310 fESDPar(rec.fESDPar),
311 fUseHLTData(rec.fUseHLTData),
315 fParentRawReader(NULL),
318 fDiamondProfile(NULL),
319 fDiamondProfileTPC(NULL),
320 fMeanVertexConstraint(rec.fMeanVertexConstraint),
324 fAlignObjArray(rec.fAlignObjArray),
325 fCDBUri(rec.fCDBUri),
327 fInitCDBCalled(rec.fInitCDBCalled),
328 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
330 fRunGlobalQA(rec.fRunGlobalQA),
331 fInLoopQA(rec.fInLoopQA),
332 fRunPlaneEff(rec.fRunPlaneEff)
336 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
337 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
339 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
340 fReconstructor[iDet] = NULL;
341 fLoader[iDet] = NULL;
342 fTracker[iDet] = NULL;
343 fQADataMaker[iDet] = NULL;
344 fQACycles[iDet] = rec.fQACycles[iDet];
346 fQADataMaker[fgkNDetectors]=NULL; //Global QA
347 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
348 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
352 //_____________________________________________________________________________
353 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
355 // assignment operator
357 this->~AliReconstruction();
358 new(this) AliReconstruction(rec);
362 //_____________________________________________________________________________
363 AliReconstruction::~AliReconstruction()
369 fSpecCDBUri.Delete();
371 AliCodeTimer::Instance()->Print();
374 //_____________________________________________________________________________
375 void AliReconstruction::InitCDB()
377 // activate a default CDB storage
378 // First check if we have any CDB storage set, because it is used
379 // to retrieve the calibration and alignment constants
381 if (fInitCDBCalled) return;
382 fInitCDBCalled = kTRUE;
384 AliCDBManager* man = AliCDBManager::Instance();
385 if (man->IsDefaultStorageSet())
387 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
388 AliWarning("Default CDB storage has been already set !");
389 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
390 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
391 fCDBUri = man->GetDefaultStorage()->GetURI();
394 if (fCDBUri.Length() > 0)
396 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
397 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
398 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
400 fCDBUri="local://$ALICE_ROOT";
401 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
402 AliWarning("Default CDB storage not yet set !!!!");
403 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
404 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
407 man->SetDefaultStorage(fCDBUri);
410 // Now activate the detector specific CDB storage locations
411 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
412 TObject* obj = fSpecCDBUri[i];
414 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
415 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
416 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
417 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
422 //_____________________________________________________________________________
423 void AliReconstruction::SetDefaultStorage(const char* uri) {
424 // Store the desired default CDB storage location
425 // Activate it later within the Run() method
431 //_____________________________________________________________________________
432 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
433 // Store a detector-specific CDB storage location
434 // Activate it later within the Run() method
436 AliCDBPath aPath(calibType);
437 if(!aPath.IsValid()){
438 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
439 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
440 if(!strcmp(calibType, fgkDetectorName[iDet])) {
441 aPath.SetPath(Form("%s/*", calibType));
442 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
446 if(!aPath.IsValid()){
447 AliError(Form("Not a valid path or detector: %s", calibType));
452 // // check that calibType refers to a "valid" detector name
453 // Bool_t isDetector = kFALSE;
454 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
455 // TString detName = fgkDetectorName[iDet];
456 // if(aPath.GetLevel0() == detName) {
457 // isDetector = kTRUE;
463 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
467 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
468 if (obj) fSpecCDBUri.Remove(obj);
469 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
473 //_____________________________________________________________________________
474 Bool_t AliReconstruction::SetRunNumberFromData()
476 // The method is called in Run() in order
477 // to set a correct run number.
478 // In case of raw data reconstruction the
479 // run number is taken from the raw data header
481 if (fSetRunNumberFromDataCalled) return kTRUE;
482 fSetRunNumberFromDataCalled = kTRUE;
484 AliCDBManager* man = AliCDBManager::Instance();
486 if(man->GetRun() > 0) {
487 AliWarning("Run number is taken from event header! Ignoring settings in AliCDBManager!");
491 AliError("No run loader is found !");
494 // read run number from gAlice
495 if(fRunLoader->GetAliRun())
496 AliCDBManager::Instance()->SetRun(fRunLoader->GetHeader()->GetRun());
499 if(fRawReader->NextEvent()) {
500 AliCDBManager::Instance()->SetRun(fRawReader->GetRunNumber());
501 fRawReader->RewindEvents();
504 AliError("No raw-data events found !");
509 AliError("Neither gAlice nor RawReader objects are found !");
519 //_____________________________________________________________________________
520 void AliReconstruction::SetCDBLock() {
521 // Set CDB lock: from now on it is forbidden to reset the run number
522 // or the default storage or to activate any further storage!
524 AliCDBManager::Instance()->SetLock(1);
527 //_____________________________________________________________________________
528 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
530 // Read the alignment objects from CDB.
531 // Each detector is supposed to have the
532 // alignment objects in DET/Align/Data CDB path.
533 // All the detector objects are then collected,
534 // sorted by geometry level (starting from ALIC) and
535 // then applied to the TGeo geometry.
536 // Finally an overlaps check is performed.
538 // Load alignment data from CDB and fill fAlignObjArray
539 if(fLoadAlignFromCDB){
541 TString detStr = detectors;
542 TString loadAlObjsListOfDets = "";
544 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
545 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
546 loadAlObjsListOfDets += fgkDetectorName[iDet];
547 loadAlObjsListOfDets += " ";
548 } // end loop over detectors
549 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
550 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
552 // Check if the array with alignment objects was
553 // provided by the user. If yes, apply the objects
554 // to the present TGeo geometry
555 if (fAlignObjArray) {
556 if (gGeoManager && gGeoManager->IsClosed()) {
557 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
558 AliError("The misalignment of one or more volumes failed!"
559 "Compare the list of simulated detectors and the list of detector alignment data!");
564 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
570 delete fAlignObjArray; fAlignObjArray=0;
575 //_____________________________________________________________________________
576 void AliReconstruction::SetGAliceFile(const char* fileName)
578 // set the name of the galice file
580 fGAliceFileName = fileName;
583 //_____________________________________________________________________________
584 void AliReconstruction::SetOption(const char* detector, const char* option)
586 // set options for the reconstruction of a detector
588 TObject* obj = fOptions.FindObject(detector);
589 if (obj) fOptions.Remove(obj);
590 fOptions.Add(new TNamed(detector, option));
594 //_____________________________________________________________________________
595 Bool_t AliReconstruction::Run(const char* input, Bool_t IsOnline)
597 // run the reconstruction
603 if (!input) input = fInput.Data();
604 TString fileName(input);
605 if (fileName.EndsWith("/")) {
606 fRawReader = new AliRawReaderFile(fileName);
607 } else if (fileName.EndsWith(".root")) {
608 fRawReader = new AliRawReaderRoot(fileName);
609 } else if (!fileName.IsNull()) {
610 fRawReader = new AliRawReaderDate(fileName);
615 AliError("Null pointer to the event structure!");
618 fRawReader = new AliRawReaderDate((void *)input);
621 if (!fEquipIdMap.IsNull() && fRawReader)
622 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
624 if (!fUseHLTData.IsNull()) {
625 // create the RawReaderHLT which performs redirection of HLT input data for
626 // the specified detectors
627 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
629 fParentRawReader=fRawReader;
630 fRawReader=pRawReader;
632 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
636 AliSysInfo::AddStamp("Start");
637 // get the run loader
638 if (!InitRunLoader()) return kFALSE;
639 AliSysInfo::AddStamp("LoadLoader");
641 // Initialize the CDB storage
644 AliSysInfo::AddStamp("LoadCDB");
646 // Set run number in CDBManager (if it is not already set by the user)
647 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
649 // Set CDB lock: from now on it is forbidden to reset the run number
650 // or the default storage or to activate any further storage!
653 // Import ideal TGeo geometry and apply misalignment
655 TString geom(gSystem->DirName(fGAliceFileName));
656 geom += "/geometry.root";
657 AliGeomManager::LoadGeometry(geom.Data());
658 if (!gGeoManager) if (fStopOnError) return kFALSE;
661 if (!MisalignGeometry(fLoadAlignData)) if (fStopOnError) return kFALSE;
662 AliSysInfo::AddStamp("LoadGeom");
665 Int_t sameQACycle = kFALSE ;
666 AliQADataMakerSteer qas ;
667 if (fRunQA && fRawReader) {
668 qas.Run(fRunLocalReconstruction, fRawReader) ;
671 // checking the QA of previous steps
675 // local reconstruction
676 if (!fRunLocalReconstruction.IsNull()) {
677 if (!RunLocalReconstruction(fRunLocalReconstruction)) {
678 if (fStopOnError) {CleanUp(); return kFALSE;}
684 if (fRunVertexFinder && !CreateVertexer()) {
690 AliSysInfo::AddStamp("Vertexer");
693 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
699 AliSysInfo::AddStamp("LoadTrackers");
701 // get the possibly already existing ESD file and tree
702 AliESDEvent* esd = new AliESDEvent(); AliESDEvent* hltesd = new AliESDEvent();
703 TFile* fileOld = NULL;
704 TTree* treeOld = NULL; TTree *hlttreeOld = NULL;
705 if (!gSystem->AccessPathName("AliESDs.root")){
706 gSystem->CopyFile("AliESDs.root", "AliESDs.old.root", kTRUE);
707 fileOld = TFile::Open("AliESDs.old.root");
708 if (fileOld && fileOld->IsOpen()) {
709 treeOld = (TTree*) fileOld->Get("esdTree");
710 if (treeOld)esd->ReadFromTree(treeOld);
711 hlttreeOld = (TTree*) fileOld->Get("HLTesdTree");
712 if (hlttreeOld) hltesd->ReadFromTree(hlttreeOld);
716 // create the ESD output file and tree
717 TFile* file = TFile::Open("AliESDs.root", "RECREATE");
718 file->SetCompressionLevel(2);
719 if (!file->IsOpen()) {
720 AliError("opening AliESDs.root failed");
721 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
724 TTree* tree = new TTree("esdTree", "Tree with ESD objects");
725 esd = new AliESDEvent();
726 esd->CreateStdContent();
727 esd->WriteToTree(tree);
729 TTree* hlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
730 hltesd = new AliESDEvent();
731 hltesd->CreateStdContent();
732 hltesd->WriteToTree(hlttree);
735 delete esd; delete hltesd;
736 esd = NULL; hltesd = NULL;
738 // create the branch with ESD additions
742 AliESDfriend *esdf = 0;
743 if (fWriteESDfriend) {
744 esdf = new AliESDfriend();
745 TBranch *br=tree->Branch("ESDfriend.","AliESDfriend", &esdf);
746 br->SetFile("AliESDfriends.root");
747 esd->AddObject(esdf);
751 // Get the GRP CDB entry
752 AliCDBEntry* entryGRP = AliCDBManager::Instance()->Get("GRP/GRP/Data");
755 fGRPList = dynamic_cast<TList*> (entryGRP->GetObject());
757 AliError("No GRP entry found in OCDB!");
760 // Get the diamond profile from OCDB
761 AliCDBEntry* entry = AliCDBManager::Instance()
762 ->Get("GRP/Calib/MeanVertex");
765 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
767 AliError("No diamond profile found in OCDB!");
771 entry = AliCDBManager::Instance()
772 ->Get("GRP/Calib/MeanVertexTPC");
775 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
777 AliError("No diamond profile found in OCDB!");
780 AliVertexerTracks tVertexer(AliTracker::GetBz());
781 if(fDiamondProfile && fMeanVertexConstraint) tVertexer.SetVtxStart(fDiamondProfile);
783 if (fRawReader) fRawReader->RewindEvents();
786 gSystem->GetProcInfo(&ProcInfo);
787 AliInfo(Form("Current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual));
790 //Initialize the QA and start of cycle for out-of-cycle QA
792 TString detStr(fFillESD);
793 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
794 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
795 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
797 AliInfo(Form("Initializing the QA data maker for %s",
798 fgkDetectorName[iDet]));
799 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
800 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
802 // qadm->StartOfCycle(AliQA::kRECPOINTS);
803 // qadm->StartOfCycle(AliQA::kESDS,"same");
808 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
809 AliInfo(Form("Initializing the global QA data maker"));
811 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
812 AliTracker::SetResidualsArray(arr);
813 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
815 qadm->StartOfCycle(AliQA::kRECPOINTS, sameQACycle);
816 qadm->StartOfCycle(AliQA::kESDS, "same");
820 //Initialize the Plane Efficiency framework
821 if (fRunPlaneEff && !InitPlaneEff()) {
822 if(fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
825 //******* The loop over events
826 for (Int_t iEvent = 0; iEvent < fRunLoader->GetNumberOfEvents(); iEvent++) {
827 if (fRawReader) fRawReader->NextEvent();
828 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
829 // copy old ESD to the new one
831 esd->ReadFromTree(treeOld);
832 treeOld->GetEntry(iEvent);
836 esd->ReadFromTree(hlttreeOld);
837 hlttreeOld->GetEntry(iEvent);
843 AliInfo(Form("processing event %d", iEvent));
845 //Start of cycle for the in-loop QA
848 TString detStr(fFillESD);
849 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
850 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
851 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
853 qadm->StartOfCycle(AliQA::kRECPOINTS, sameQACycle);
854 qadm->StartOfCycle(AliQA::kESDS, "same") ;
858 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
859 qadm->StartOfCycle(AliQA::kRECPOINTS, sameQACycle);
860 qadm->StartOfCycle(AliQA::kESDS, "same");
864 fRunLoader->GetEvent(iEvent);
867 sprintf(aFileName, "ESD_%d.%d_final.root",
868 fRunLoader->GetHeader()->GetRun(),
869 fRunLoader->GetHeader()->GetEventNrInRun());
870 if (!gSystem->AccessPathName(aFileName)) continue;
872 // local signle event reconstruction
873 if (!fRunLocalReconstruction.IsNull()) {
874 if (!RunLocalEventReconstruction(fRunLocalReconstruction)) {
875 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
879 esd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
880 hltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
881 esd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
882 hltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
884 // Set magnetic field from the tracker
885 esd->SetMagneticField(AliTracker::GetBz());
886 hltesd->SetMagneticField(AliTracker::GetBz());
890 // Fill raw-data error log into the ESD
891 if (fRawReader) FillRawDataErrorLog(iEvent,esd);
894 if (fRunVertexFinder) {
895 if (!ReadESD(esd, "vertex")) {
896 if (!RunVertexFinder(esd)) {
897 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
899 if (fCheckPointLevel > 0) WriteESD(esd, "vertex");
904 if (!fRunTracking.IsNull()) {
905 if (fRunHLTTracking) {
906 hltesd->SetPrimaryVertexSPD(esd->GetVertex());
907 if (!RunHLTTracking(hltesd)) {
908 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
914 if (!fRunTracking.IsNull()) {
915 if (fRunMuonTracking) {
916 if (!RunMuonTracking(esd)) {
917 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
923 if (!fRunTracking.IsNull()) {
924 if (!ReadESD(esd, "tracking")) {
925 if (!RunTracking(esd)) {
926 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
928 if (fCheckPointLevel > 0) WriteESD(esd, "tracking");
933 if (!fFillESD.IsNull()) {
934 if (!FillESD(esd, fFillESD)) {
935 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
939 // fill Event header information from the RawEventHeader
940 if (fRawReader){FillRawEventHeaderESD(esd);}
943 AliESDpid::MakePID(esd);
944 if (fCheckPointLevel > 1) WriteESD(esd, "PID");
946 if (fFillTriggerESD) {
947 if (!ReadESD(esd, "trigger")) {
948 if (!FillTriggerESD(esd)) {
949 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
951 if (fCheckPointLevel > 1) WriteESD(esd, "trigger");
958 // Propagate track to the beam pipe (if not laready done by ITS)
960 const Int_t ntracks = esd->GetNumberOfTracks();
961 const Double_t kBz = esd->GetMagneticField();
962 const Double_t kRadius = 2.8; //something less than the beam pipe radius
965 UShort_t *selectedIdx=new UShort_t[ntracks];
967 for (Int_t itrack=0; itrack<ntracks; itrack++){
968 const Double_t kMaxStep = 5; //max step over the material
971 AliESDtrack *track = esd->GetTrack(itrack);
972 if (!track) continue;
974 AliExternalTrackParam *tpcTrack =
975 (AliExternalTrackParam *)track->GetTPCInnerParam();
979 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kTRUE);
984 Int_t n=trkArray.GetEntriesFast();
985 selectedIdx[n]=track->GetID();
986 trkArray.AddLast(tpcTrack);
989 if (track->GetX() < kRadius) continue;
992 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kTRUE);
994 track->RelateToVertex(esd->GetPrimaryVertexSPD(), kBz, kRadius);
999 // Improve the reconstructed primary vertex position using the tracks
1001 TObject *obj = fOptions.FindObject("ITS");
1003 TString optITS = obj->GetTitle();
1004 if (optITS.Contains("cosmics") || optITS.Contains("COSMICS"))
1005 fRunVertexFinderTracks=kFALSE;
1007 if (fRunVertexFinderTracks) {
1008 // TPC + ITS primary vertex
1009 tVertexer.SetITSrefitRequired();
1010 if(fDiamondProfile && fMeanVertexConstraint) {
1011 tVertexer.SetVtxStart(fDiamondProfile);
1013 tVertexer.SetConstraintOff();
1015 AliESDVertex *pvtx=tVertexer.FindPrimaryVertex(esd);
1017 if (pvtx->GetStatus()) {
1018 esd->SetPrimaryVertex(pvtx);
1019 for (Int_t i=0; i<ntracks; i++) {
1020 AliESDtrack *t = esd->GetTrack(i);
1021 t->RelateToVertex(pvtx, kBz, kRadius);
1026 // TPC-only primary vertex
1027 tVertexer.SetITSrefitNotRequired();
1028 if(fDiamondProfileTPC && fMeanVertexConstraint) {
1029 tVertexer.SetVtxStart(fDiamondProfileTPC);
1031 tVertexer.SetConstraintOff();
1033 pvtx=tVertexer.FindPrimaryVertex(&trkArray,selectedIdx);
1035 if (pvtx->GetStatus()) {
1036 esd->SetPrimaryVertexTPC(pvtx);
1037 Int_t nsel=trkArray.GetEntriesFast();
1038 for (Int_t i=0; i<nsel; i++) {
1039 AliExternalTrackParam *t =
1040 (AliExternalTrackParam *)trkArray.UncheckedAt(i);
1041 t->PropagateToDCA(pvtx, kBz, kRadius);
1047 delete[] selectedIdx;
1049 if(fDiamondProfile) esd->SetDiamond(fDiamondProfile);
1054 AliV0vertexer vtxer;
1055 vtxer.Tracks2V0vertices(esd);
1057 if (fRunCascadeFinder) {
1059 AliCascadeVertexer cvtxer;
1060 cvtxer.V0sTracks2CascadeVertices(esd);
1065 if (fCleanESD) CleanESD(esd);
1068 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1069 if (qadm) qadm->Exec(AliQA::kESDS, esd);
1072 if (fWriteESDfriend) {
1073 esdf->~AliESDfriend();
1074 new (esdf) AliESDfriend(); // Reset...
1075 esd->GetESDfriend(esdf);
1082 if (fCheckPointLevel > 0) WriteESD(esd, "final");
1085 if (fWriteESDfriend) {
1086 esdf->~AliESDfriend();
1087 new (esdf) AliESDfriend(); // Reset...
1090 gSystem->GetProcInfo(&ProcInfo);
1091 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1094 // End of cycle for the in-loop QA
1097 RunQA(fFillESD.Data(), esd);
1098 TString detStr(fFillESD);
1099 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1100 if (!IsSelected(fgkDetectorName[iDet], detStr))
1102 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1105 qadm->EndOfCycle(AliQA::kRECPOINTS);
1106 qadm->EndOfCycle(AliQA::kESDS);
1111 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1113 qadm->EndOfCycle(AliQA::kRECPOINTS);
1114 qadm->EndOfCycle(AliQA::kESDS);
1120 //******** End of the loop over events
1124 tree->GetUserInfo()->Add(esd);
1125 hlttree->GetUserInfo()->Add(hltesd);
1127 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1128 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1130 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1131 cdbMapCopy->SetOwner(1);
1132 cdbMapCopy->SetName("cdbMap");
1133 TIter iter(cdbMap->GetTable());
1136 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1137 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1138 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1139 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1142 TList *cdbListCopy = new TList();
1143 cdbListCopy->SetOwner(1);
1144 cdbListCopy->SetName("cdbList");
1146 TIter iter2(cdbList);
1149 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1150 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1153 tree->GetUserInfo()->Add(cdbMapCopy);
1154 tree->GetUserInfo()->Add(cdbListCopy);
1157 if(fESDPar.Contains("ESD.par")){
1158 AliInfo("Attaching ESD.par to Tree");
1159 TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1160 tree->GetUserInfo()->Add(fn);
1166 if (fWriteESDfriend)
1167 tree->SetBranchStatus("ESDfriend*",0);
1168 // we want to have only one tree version number
1169 tree->Write(tree->GetName(),TObject::kOverwrite);
1172 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1173 if (fRunPlaneEff && !FinishPlaneEff()) {
1174 AliWarning("Finish PlaneEff evaluation failed");
1178 CleanUp(file, fileOld);
1181 AliWarning("AOD creation not supported anymore during reconstruction. See ANALYSIS/AliAnalysisTaskESDfilter.cxx instead.");
1184 // Create tags for the events in the ESD tree (the ESD tree is always present)
1185 // In case of empty events the tags will contain dummy values
1186 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1187 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPList);
1189 AliWarning("AOD tag creation not supported anymore during reconstruction.");
1192 //Finish QA and end of cycle for out-of-loop QA
1195 qas.Run(fRunLocalReconstruction.Data(), AliQA::kRECPOINTS, sameQACycle);
1197 qas.Run(fRunTracking.Data(), AliQA::kESDS, sameQACycle);
1200 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1202 qadm->EndOfCycle(AliQA::kRECPOINTS);
1203 qadm->EndOfCycle(AliQA::kESDS);
1209 // Cleanup of CDB manager: cache and active storages!
1210 AliCDBManager::Instance()->ClearCache();
1217 //_____________________________________________________________________________
1218 Bool_t AliReconstruction::RunLocalReconstruction(const TString& /*detectors*/)
1220 // run the local reconstruction
1221 static Int_t eventNr=0;
1222 AliCodeTimerAuto("")
1224 // AliCDBManager* man = AliCDBManager::Instance();
1225 // Bool_t origCache = man->GetCacheFlag();
1227 // TString detStr = detectors;
1228 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1229 // if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1230 // AliReconstructor* reconstructor = GetReconstructor(iDet);
1231 // if (!reconstructor) continue;
1232 // if (reconstructor->HasLocalReconstruction()) continue;
1234 // AliCodeTimerStart(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1235 // AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1237 // AliCodeTimerStart(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1238 // AliInfo(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1240 // man->SetCacheFlag(kTRUE);
1241 // TString calibPath = Form("%s/Calib/*", fgkDetectorName[iDet]);
1242 // man->GetAll(calibPath); // entries are cached!
1244 // AliCodeTimerStop(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1246 // if (fRawReader) {
1247 // fRawReader->RewindEvents();
1248 // reconstructor->Reconstruct(fRunLoader, fRawReader);
1250 // reconstructor->Reconstruct(fRunLoader);
1253 // AliCodeTimerStop(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1254 // AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr));
1256 // // unload calibration data
1257 // man->UnloadFromCache(calibPath);
1258 // //man->ClearCache();
1261 // man->SetCacheFlag(origCache);
1263 // if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1264 // AliError(Form("the following detectors were not found: %s",
1266 // if (fStopOnError) return kFALSE;
1273 //_____________________________________________________________________________
1274 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1276 // run the local reconstruction
1278 static Int_t eventNr=0;
1279 AliCodeTimerAuto("")
1281 TString detStr = detectors;
1282 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1283 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1284 AliReconstructor* reconstructor = GetReconstructor(iDet);
1285 if (!reconstructor) continue;
1286 AliLoader* loader = fLoader[iDet];
1288 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1291 // conversion of digits
1292 if (fRawReader && reconstructor->HasDigitConversion()) {
1293 AliInfo(Form("converting raw data digits into root objects for %s",
1294 fgkDetectorName[iDet]));
1295 AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1296 fgkDetectorName[iDet]));
1297 loader->LoadDigits("update");
1298 loader->CleanDigits();
1299 loader->MakeDigitsContainer();
1300 TTree* digitsTree = loader->TreeD();
1301 reconstructor->ConvertDigits(fRawReader, digitsTree);
1302 loader->WriteDigits("OVERWRITE");
1303 loader->UnloadDigits();
1305 // local reconstruction
1306 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1307 AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1308 loader->LoadRecPoints("update");
1309 loader->CleanRecPoints();
1310 loader->MakeRecPointsContainer();
1311 TTree* clustersTree = loader->TreeR();
1312 if (fRawReader && !reconstructor->HasDigitConversion()) {
1313 reconstructor->Reconstruct(fRawReader, clustersTree);
1315 loader->LoadDigits("read");
1316 TTree* digitsTree = loader->TreeD();
1318 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1319 if (fStopOnError) return kFALSE;
1321 reconstructor->Reconstruct(digitsTree, clustersTree);
1323 loader->UnloadDigits();
1326 // In-loop QA for local reconstrucion
1327 if (fRunQA && fInLoopQA) {
1328 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1331 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1333 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1335 qadm->Exec(AliQA::kRECPOINTS, clustersTree) ;
1338 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1342 loader->WriteRecPoints("OVERWRITE");
1343 loader->UnloadRecPoints();
1344 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1347 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1348 AliError(Form("the following detectors were not found: %s",
1350 if (fStopOnError) return kFALSE;
1356 //_____________________________________________________________________________
1357 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1359 // run the barrel tracking
1361 AliCodeTimerAuto("")
1363 AliESDVertex* vertex = NULL;
1364 Double_t vtxPos[3] = {0, 0, 0};
1365 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1366 TArrayF mcVertex(3);
1367 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1368 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1369 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1373 if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1374 AliInfo("running the ITS vertex finder");
1375 if (fLoader[0]) fLoader[0]->LoadRecPoints();
1376 vertex = fVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber());
1377 if (fLoader[0]) fLoader[0]->UnloadRecPoints();
1379 AliWarning("Vertex not found");
1380 vertex = new AliESDVertex();
1381 vertex->SetName("default");
1384 vertex->SetName("reconstructed");
1388 AliInfo("getting the primary vertex from MC");
1389 vertex = new AliESDVertex(vtxPos, vtxErr);
1393 vertex->GetXYZ(vtxPos);
1394 vertex->GetSigmaXYZ(vtxErr);
1396 AliWarning("no vertex reconstructed");
1397 vertex = new AliESDVertex(vtxPos, vtxErr);
1399 esd->SetPrimaryVertexSPD(vertex);
1400 // if SPD multiplicity has been determined, it is stored in the ESD
1401 AliMultiplicity *mult = fVertexer->GetMultiplicity();
1402 if(mult)esd->SetMultiplicity(mult);
1404 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1405 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1412 //_____________________________________________________________________________
1413 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
1415 // run the HLT barrel tracking
1417 AliCodeTimerAuto("")
1420 AliError("Missing runLoader!");
1424 AliInfo("running HLT tracking");
1426 // Get a pointer to the HLT reconstructor
1427 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
1428 if (!reconstructor) return kFALSE;
1431 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1432 TString detName = fgkDetectorName[iDet];
1433 AliDebug(1, Form("%s HLT tracking", detName.Data()));
1434 reconstructor->SetOption(detName.Data());
1435 AliTracker *tracker = reconstructor->CreateTracker();
1437 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
1438 if (fStopOnError) return kFALSE;
1442 Double_t vtxErr[3]={0.005,0.005,0.010};
1443 const AliESDVertex *vertex = esd->GetVertex();
1444 vertex->GetXYZ(vtxPos);
1445 tracker->SetVertex(vtxPos,vtxErr);
1447 fLoader[iDet]->LoadRecPoints("read");
1448 TTree* tree = fLoader[iDet]->TreeR();
1450 AliError(Form("Can't get the %s cluster tree", detName.Data()));
1453 tracker->LoadClusters(tree);
1455 if (tracker->Clusters2Tracks(esd) != 0) {
1456 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
1460 tracker->UnloadClusters();
1468 //_____________________________________________________________________________
1469 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
1471 // run the muon spectrometer tracking
1473 AliCodeTimerAuto("")
1476 AliError("Missing runLoader!");
1479 Int_t iDet = 7; // for MUON
1481 AliInfo("is running...");
1483 // Get a pointer to the MUON reconstructor
1484 AliReconstructor *reconstructor = GetReconstructor(iDet);
1485 if (!reconstructor) return kFALSE;
1488 TString detName = fgkDetectorName[iDet];
1489 AliDebug(1, Form("%s tracking", detName.Data()));
1490 AliTracker *tracker = reconstructor->CreateTracker();
1492 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1497 fLoader[iDet]->LoadRecPoints("read");
1499 tracker->LoadClusters(fLoader[iDet]->TreeR());
1501 Int_t rv = tracker->Clusters2Tracks(esd);
1505 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1509 fLoader[iDet]->UnloadRecPoints();
1511 tracker->UnloadClusters();
1519 //_____________________________________________________________________________
1520 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
1522 // run the barrel tracking
1523 static Int_t eventNr=0;
1524 AliCodeTimerAuto("")
1526 AliInfo("running tracking");
1528 //Fill the ESD with the T0 info (will be used by the TOF)
1529 if (fReconstructor[11] && fLoader[11]) {
1530 fLoader[11]->LoadRecPoints("READ");
1531 TTree *treeR = fLoader[11]->TreeR();
1532 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
1535 // pass 1: TPC + ITS inwards
1536 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1537 if (!fTracker[iDet]) continue;
1538 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
1541 fLoader[iDet]->LoadRecPoints("read");
1542 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
1543 TTree* tree = fLoader[iDet]->TreeR();
1545 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1548 fTracker[iDet]->LoadClusters(tree);
1549 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1551 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
1552 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1555 if (fCheckPointLevel > 1) {
1556 WriteESD(esd, Form("%s.tracking", fgkDetectorName[iDet]));
1558 // preliminary PID in TPC needed by the ITS tracker
1560 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1561 AliESDpid::MakePID(esd);
1563 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
1566 // pass 2: ALL backwards
1568 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1569 if (!fTracker[iDet]) continue;
1570 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
1573 if (iDet > 1) { // all except ITS, TPC
1575 fLoader[iDet]->LoadRecPoints("read");
1576 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
1577 tree = fLoader[iDet]->TreeR();
1579 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1582 fTracker[iDet]->LoadClusters(tree);
1583 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1587 if (iDet>1) // start filling residuals for the "outer" detectors
1588 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1590 if (fTracker[iDet]->PropagateBack(esd) != 0) {
1591 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
1594 if (fCheckPointLevel > 1) {
1595 WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
1599 if (iDet > 2) { // all except ITS, TPC, TRD
1600 fTracker[iDet]->UnloadClusters();
1601 fLoader[iDet]->UnloadRecPoints();
1603 // updated PID in TPC needed by the ITS tracker -MI
1605 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1606 AliESDpid::MakePID(esd);
1608 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1610 //stop filling residuals for the "outer" detectors
1611 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1613 // write space-points to the ESD in case alignment data output
1615 if (fWriteAlignmentData)
1616 WriteAlignmentData(esd);
1618 // pass 3: TRD + TPC + ITS refit inwards
1620 for (Int_t iDet = 2; iDet >= 0; iDet--) {
1621 if (!fTracker[iDet]) continue;
1622 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
1625 if (iDet<2) // start filling residuals for TPC and ITS
1626 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1628 if (fTracker[iDet]->RefitInward(esd) != 0) {
1629 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
1632 // run postprocessing
1633 if (fTracker[iDet]->PostProcess(esd) != 0) {
1634 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
1637 if (fCheckPointLevel > 1) {
1638 WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
1640 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1642 fTracker[iDet]->UnloadClusters();
1643 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
1644 fLoader[iDet]->UnloadRecPoints();
1645 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
1647 // stop filling residuals for TPC and ITS
1648 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1654 //_____________________________________________________________________________
1655 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
1657 // Remove the data which are not needed for the physics analysis.
1660 Int_t nTracks=esd->GetNumberOfTracks();
1661 Int_t nV0s=esd->GetNumberOfV0s();
1663 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
1665 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
1666 Bool_t rc=esd->Clean(cleanPars);
1668 nTracks=esd->GetNumberOfTracks();
1669 nV0s=esd->GetNumberOfV0s();
1671 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
1676 //_____________________________________________________________________________
1677 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
1679 // fill the event summary data
1681 AliCodeTimerAuto("")
1682 static Int_t eventNr=0;
1683 TString detStr = detectors;
1685 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1686 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1687 AliReconstructor* reconstructor = GetReconstructor(iDet);
1688 if (!reconstructor) continue;
1689 if (!ReadESD(esd, fgkDetectorName[iDet])) {
1690 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
1691 TTree* clustersTree = NULL;
1692 if (fLoader[iDet]) {
1693 fLoader[iDet]->LoadRecPoints("read");
1694 clustersTree = fLoader[iDet]->TreeR();
1695 if (!clustersTree) {
1696 AliError(Form("Can't get the %s clusters tree",
1697 fgkDetectorName[iDet]));
1698 if (fStopOnError) return kFALSE;
1701 if (fRawReader && !reconstructor->HasDigitConversion()) {
1702 reconstructor->FillESD(fRawReader, clustersTree, esd);
1704 TTree* digitsTree = NULL;
1705 if (fLoader[iDet]) {
1706 fLoader[iDet]->LoadDigits("read");
1707 digitsTree = fLoader[iDet]->TreeD();
1709 AliError(Form("Can't get the %s digits tree",
1710 fgkDetectorName[iDet]));
1711 if (fStopOnError) return kFALSE;
1714 reconstructor->FillESD(digitsTree, clustersTree, esd);
1715 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
1717 if (fLoader[iDet]) {
1718 fLoader[iDet]->UnloadRecPoints();
1721 if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
1725 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1726 AliError(Form("the following detectors were not found: %s",
1728 if (fStopOnError) return kFALSE;
1730 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
1735 //_____________________________________________________________________________
1736 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
1738 // Reads the trigger decision which is
1739 // stored in Trigger.root file and fills
1740 // the corresponding esd entries
1742 AliCodeTimerAuto("")
1744 AliInfo("Filling trigger information into the ESD");
1746 AliCentralTrigger *aCTP = NULL;
1749 AliCTPRawStream input(fRawReader);
1750 if (!input.Next()) {
1751 AliError("No valid CTP (trigger) DDL raw data is found ! The trigger information is not stored in the ESD !");
1754 esd->SetTriggerMask(input.GetClassMask());
1755 esd->SetTriggerCluster(input.GetClusterMask());
1757 aCTP = new AliCentralTrigger();
1758 TString configstr("");
1759 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
1760 AliError("No trigger configuration found in OCDB! The trigger classes information will no be stored in ESD!");
1765 AliRunLoader *runloader = AliRunLoader::GetRunLoader();
1767 if (!runloader->LoadTrigger()) {
1768 aCTP = runloader->GetTrigger();
1769 esd->SetTriggerMask(aCTP->GetClassMask());
1770 esd->SetTriggerCluster(aCTP->GetClusterMask());
1773 AliWarning("No trigger can be loaded! The trigger information is not stored in the ESD !");
1778 AliError("No run loader is available! The trigger information is not stored in the ESD !");
1783 // Now fill the trigger class names into AliESDRun object
1784 AliTriggerConfiguration *config = aCTP->GetConfiguration();
1786 AliError("No trigger configuration has been found! The trigger classes information will no be stored in ESD!");
1790 const TObjArray& classesArray = config->GetClasses();
1791 Int_t nclasses = classesArray.GetEntriesFast();
1792 for( Int_t j=0; j<nclasses; j++ ) {
1793 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At( j );
1794 Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
1795 esd->SetTriggerClass(trclass->GetName(),trindex);
1805 //_____________________________________________________________________________
1806 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
1809 // Filling information from RawReader Header
1812 AliInfo("Filling information from RawReader Header");
1813 esd->SetBunchCrossNumber(0);
1814 esd->SetOrbitNumber(0);
1815 esd->SetPeriodNumber(0);
1816 esd->SetTimeStamp(0);
1817 esd->SetEventType(0);
1818 const AliRawEventHeaderBase * eventHeader = fRawReader->GetEventHeader();
1821 const UInt_t *id = eventHeader->GetP("Id");
1822 esd->SetBunchCrossNumber((id)[1]&0x00000fff);
1823 esd->SetOrbitNumber((((id)[0]<<20)&0xf00000)|(((id)[1]>>12)&0xfffff));
1824 esd->SetPeriodNumber(((id)[0]>>4)&0x0fffffff);
1826 esd->SetTimeStamp((eventHeader->Get("Timestamp")));
1827 esd->SetEventType((eventHeader->Get("Type")));
1834 //_____________________________________________________________________________
1835 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
1837 // check whether detName is contained in detectors
1838 // if yes, it is removed from detectors
1840 // check if all detectors are selected
1841 if ((detectors.CompareTo("ALL") == 0) ||
1842 detectors.BeginsWith("ALL ") ||
1843 detectors.EndsWith(" ALL") ||
1844 detectors.Contains(" ALL ")) {
1849 // search for the given detector
1850 Bool_t result = kFALSE;
1851 if ((detectors.CompareTo(detName) == 0) ||
1852 detectors.BeginsWith(detName+" ") ||
1853 detectors.EndsWith(" "+detName) ||
1854 detectors.Contains(" "+detName+" ")) {
1855 detectors.ReplaceAll(detName, "");
1859 // clean up the detectors string
1860 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1861 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1862 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1867 //_____________________________________________________________________________
1868 Bool_t AliReconstruction::InitRunLoader()
1870 // get or create the run loader
1872 if (gAlice) delete gAlice;
1875 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
1876 // load all base libraries to get the loader classes
1877 TString libs = gSystem->GetLibraries();
1878 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1879 TString detName = fgkDetectorName[iDet];
1880 if (detName == "HLT") continue;
1881 if (libs.Contains("lib" + detName + "base.so")) continue;
1882 gSystem->Load("lib" + detName + "base.so");
1884 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
1886 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1890 fRunLoader->CdGAFile();
1891 if (gFile->GetKey(AliRunLoader::GetGAliceName())) {
1892 if (fRunLoader->LoadgAlice() == 0) {
1893 gAlice = fRunLoader->GetAliRun();
1894 AliTracker::SetFieldMap(gAlice->Field(),fUniformField);
1897 if (!gAlice && !fRawReader) {
1898 AliError(Form("no gAlice object found in file %s",
1899 fGAliceFileName.Data()));
1904 //PH This is a temporary fix to give access to the kinematics
1905 //PH that is needed for the labels of ITS clusters
1906 fRunLoader->LoadHeader();
1907 fRunLoader->LoadKinematics();
1909 } else { // galice.root does not exist
1911 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
1915 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
1916 AliConfig::GetDefaultEventFolderName(),
1919 AliError(Form("could not create run loader in file %s",
1920 fGAliceFileName.Data()));
1924 fRunLoader->MakeTree("E");
1926 while (fRawReader->NextEvent()) {
1927 fRunLoader->SetEventNumber(iEvent);
1928 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1930 fRunLoader->MakeTree("H");
1931 fRunLoader->TreeE()->Fill();
1934 fRawReader->RewindEvents();
1935 if (fNumberOfEventsPerFile > 0)
1936 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
1938 fRunLoader->SetNumberOfEventsPerFile(iEvent);
1939 fRunLoader->WriteHeader("OVERWRITE");
1940 fRunLoader->CdGAFile();
1941 fRunLoader->Write(0, TObject::kOverwrite);
1942 // AliTracker::SetFieldMap(???);
1948 //_____________________________________________________________________________
1949 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
1951 // get the reconstructor object and the loader for a detector
1953 if (fReconstructor[iDet]) return fReconstructor[iDet];
1955 // load the reconstructor object
1956 TPluginManager* pluginManager = gROOT->GetPluginManager();
1957 TString detName = fgkDetectorName[iDet];
1958 TString recName = "Ali" + detName + "Reconstructor";
1959 if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT")) return NULL;
1961 AliReconstructor* reconstructor = NULL;
1962 // first check if a plugin is defined for the reconstructor
1963 TPluginHandler* pluginHandler =
1964 pluginManager->FindHandler("AliReconstructor", detName);
1965 // if not, add a plugin for it
1966 if (!pluginHandler) {
1967 AliDebug(1, Form("defining plugin for %s", recName.Data()));
1968 TString libs = gSystem->GetLibraries();
1969 if (libs.Contains("lib" + detName + "base.so") ||
1970 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
1971 pluginManager->AddHandler("AliReconstructor", detName,
1972 recName, detName + "rec", recName + "()");
1974 pluginManager->AddHandler("AliReconstructor", detName,
1975 recName, detName, recName + "()");
1977 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
1979 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
1980 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
1982 if (reconstructor) {
1983 TObject* obj = fOptions.FindObject(detName.Data());
1984 if (obj) reconstructor->SetOption(obj->GetTitle());
1985 reconstructor->Init();
1986 fReconstructor[iDet] = reconstructor;
1989 // get or create the loader
1990 if (detName != "HLT") {
1991 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
1992 if (!fLoader[iDet]) {
1993 AliConfig::Instance()
1994 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
1996 // first check if a plugin is defined for the loader
1998 pluginManager->FindHandler("AliLoader", detName);
1999 // if not, add a plugin for it
2000 if (!pluginHandler) {
2001 TString loaderName = "Ali" + detName + "Loader";
2002 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2003 pluginManager->AddHandler("AliLoader", detName,
2004 loaderName, detName + "base",
2005 loaderName + "(const char*, TFolder*)");
2006 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2008 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2010 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2011 fRunLoader->GetEventFolder());
2013 if (!fLoader[iDet]) { // use default loader
2014 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2016 if (!fLoader[iDet]) {
2017 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2018 if (fStopOnError) return NULL;
2020 fRunLoader->AddLoader(fLoader[iDet]);
2021 fRunLoader->CdGAFile();
2022 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2023 fRunLoader->Write(0, TObject::kOverwrite);
2028 return reconstructor;
2031 //_____________________________________________________________________________
2032 Bool_t AliReconstruction::CreateVertexer()
2034 // create the vertexer
2037 AliReconstructor* itsReconstructor = GetReconstructor(0);
2038 if (itsReconstructor) {
2039 fVertexer = itsReconstructor->CreateVertexer();
2042 AliWarning("couldn't create a vertexer for ITS");
2043 if (fStopOnError) return kFALSE;
2049 //_____________________________________________________________________________
2050 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2052 // create the trackers
2054 TString detStr = detectors;
2055 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2056 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2057 AliReconstructor* reconstructor = GetReconstructor(iDet);
2058 if (!reconstructor) continue;
2059 TString detName = fgkDetectorName[iDet];
2060 if (detName == "HLT") {
2061 fRunHLTTracking = kTRUE;
2064 if (detName == "MUON") {
2065 fRunMuonTracking = kTRUE;
2070 fTracker[iDet] = reconstructor->CreateTracker();
2071 if (!fTracker[iDet] && (iDet < 7)) {
2072 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2073 if (fStopOnError) return kFALSE;
2075 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2081 //_____________________________________________________________________________
2082 void AliReconstruction::CleanUp(TFile* file, TFile* fileOld)
2084 // delete trackers and the run loader and close and delete the file
2086 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2087 delete fReconstructor[iDet];
2088 fReconstructor[iDet] = NULL;
2089 fLoader[iDet] = NULL;
2090 delete fTracker[iDet];
2091 fTracker[iDet] = NULL;
2092 // delete fQADataMaker[iDet];
2093 // fQADataMaker[iDet] = NULL;
2098 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2099 delete fDiamondProfile;
2100 fDiamondProfile = NULL;
2101 delete fDiamondProfileTPC;
2102 fDiamondProfileTPC = NULL;
2112 if (fParentRawReader) delete fParentRawReader;
2113 fParentRawReader=NULL;
2123 gSystem->Unlink("AliESDs.old.root");
2127 //_____________________________________________________________________________
2129 Bool_t AliReconstruction::ReadESD(AliESDEvent*& esd, const char* recStep) const
2131 // read the ESD event from a file
2133 if (!esd) return kFALSE;
2135 sprintf(fileName, "ESD_%d.%d_%s.root",
2136 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2137 if (gSystem->AccessPathName(fileName)) return kFALSE;
2139 AliInfo(Form("reading ESD from file %s", fileName));
2140 AliDebug(1, Form("reading ESD from file %s", fileName));
2141 TFile* file = TFile::Open(fileName);
2142 if (!file || !file->IsOpen()) {
2143 AliError(Form("opening %s failed", fileName));
2150 esd = (AliESDEvent*) file->Get("ESD");
2159 //_____________________________________________________________________________
2160 void AliReconstruction::WriteESD(AliESDEvent* esd, const char* recStep) const
2162 // write the ESD event to a file
2166 sprintf(fileName, "ESD_%d.%d_%s.root",
2167 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2169 AliDebug(1, Form("writing ESD to file %s", fileName));
2170 TFile* file = TFile::Open(fileName, "recreate");
2171 if (!file || !file->IsOpen()) {
2172 AliError(Form("opening %s failed", fileName));
2181 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2183 // Write space-points which are then used in the alignment procedures
2184 // For the moment only ITS, TRD and TPC
2186 // Load TOF clusters
2188 fLoader[3]->LoadRecPoints("read");
2189 TTree* tree = fLoader[3]->TreeR();
2191 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[3]));
2194 fTracker[3]->LoadClusters(tree);
2196 Int_t ntracks = esd->GetNumberOfTracks();
2197 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2199 AliESDtrack *track = esd->GetTrack(itrack);
2202 for (Int_t iDet = 3; iDet >= 0; iDet--)
2203 nsp += track->GetNcls(iDet);
2205 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2206 track->SetTrackPointArray(sp);
2208 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2209 AliTracker *tracker = fTracker[iDet];
2210 if (!tracker) continue;
2211 Int_t nspdet = track->GetNcls(iDet);
2212 if (nspdet <= 0) continue;
2213 track->GetClusters(iDet,idx);
2217 while (isp2 < nspdet) {
2219 TString dets = fgkDetectorName[iDet];
2220 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2221 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2222 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2223 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2224 isvalid = tracker->GetTrackPointTrackingError(idx[isp2],p,track);
2226 isvalid = tracker->GetTrackPoint(idx[isp2],p);
2229 const Int_t kNTPCmax = 159;
2230 if (iDet==1 && isp2>kNTPCmax) break; // to be fixed
2231 if (!isvalid) continue;
2232 sp->AddPoint(isptrack,&p); isptrack++; isp++;
2238 fTracker[3]->UnloadClusters();
2239 fLoader[3]->UnloadRecPoints();
2243 //_____________________________________________________________________________
2244 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2246 // The method reads the raw-data error log
2247 // accumulated within the rawReader.
2248 // It extracts the raw-data errors related to
2249 // the current event and stores them into
2250 // a TClonesArray inside the esd object.
2252 if (!fRawReader) return;
2254 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2256 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2258 if (iEvent != log->GetEventNumber()) continue;
2260 esd->AddRawDataErrorLog(log);
2265 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString fName){
2266 // Dump a file content into a char in TNamed
2268 in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2269 Int_t kBytes = (Int_t)in.tellg();
2270 printf("Size: %d \n",kBytes);
2273 char* memblock = new char [kBytes];
2274 in.seekg (0, ios::beg);
2275 in.read (memblock, kBytes);
2277 TString fData(memblock,kBytes);
2278 fn = new TNamed(fName,fData);
2279 printf("fData Size: %d \n",fData.Sizeof());
2280 printf("fName Size: %d \n",fName.Sizeof());
2281 printf("fn Size: %d \n",fn->Sizeof());
2285 AliInfo(Form("Could not Open %s\n",fPath.Data()));
2291 void AliReconstruction::TNamedToFile(TTree* fTree, TString fName){
2292 // This is not really needed in AliReconstruction at the moment
2293 // but can serve as a template
2295 TList *fList = fTree->GetUserInfo();
2296 TNamed *fn = (TNamed*)fList->FindObject(fName.Data());
2297 printf("fn Size: %d \n",fn->Sizeof());
2299 TString fTmp(fn->GetName()); // to be 100% sure in principle fName also works
2300 const char* cdata = fn->GetTitle();
2301 printf("fTmp Size %d\n",fTmp.Sizeof());
2303 int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2304 printf("calculated size %d\n",size);
2305 ofstream out(fName.Data(),ios::out | ios::binary);
2306 out.write(cdata,size);
2311 //_____________________________________________________________________________
2312 AliQADataMakerRec * AliReconstruction::GetQADataMaker(Int_t iDet)
2314 // get the quality assurance data maker object and the loader for a detector
2316 if (fQADataMaker[iDet])
2317 return fQADataMaker[iDet];
2319 AliQADataMakerRec * qadm = NULL;
2320 if (iDet == fgkNDetectors) { //Global QA
2321 qadm = new AliGlobalQADataMaker();
2322 fQADataMaker[iDet] = qadm;
2326 // load the QA data maker object
2327 TPluginManager* pluginManager = gROOT->GetPluginManager();
2328 TString detName = fgkDetectorName[iDet];
2329 TString qadmName = "Ali" + detName + "QADataMakerRec";
2330 if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT"))
2333 // first check if a plugin is defined for the quality assurance data maker
2334 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2335 // if not, add a plugin for it
2336 if (!pluginHandler) {
2337 AliDebug(1, Form("defining plugin for %s", qadmName.Data()));
2338 TString libs = gSystem->GetLibraries();
2339 if (libs.Contains("lib" + detName + "base.so") ||
2340 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2341 pluginManager->AddHandler("AliQADataMakerRec", detName,
2342 qadmName, detName + "qadm", qadmName + "()");
2344 pluginManager->AddHandler("AliQADataMakerRec", detName,
2345 qadmName, detName, qadmName + "()");
2347 pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2349 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2350 qadm = (AliQADataMakerRec *) pluginHandler->ExecPlugin(0);
2353 fQADataMaker[iDet] = qadm;
2358 //_____________________________________________________________________________
2359 Bool_t AliReconstruction::RunQA(const char* detectors, AliESDEvent *& esd)
2361 // run the Quality Assurance data producer
2363 AliCodeTimerAuto("")
2364 TString detStr = detectors;
2365 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2366 if (!IsSelected(fgkDetectorName[iDet], detStr))
2368 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
2371 AliCodeTimerStart(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2372 AliInfo(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2374 qadm->Exec(AliQA::kESDS, esd) ;
2377 AliCodeTimerStop(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2379 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2380 AliError(Form("the following detectors were not found: %s",
2390 //_____________________________________________________________________________
2391 void AliReconstruction::CheckQA()
2393 // check the QA of SIM for this run and remove the detectors
2394 // with status Fatal
2396 TString newRunLocalReconstruction ;
2397 TString newRunTracking ;
2398 TString newFillESD ;
2400 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
2401 TString detName(AliQA::GetDetName(iDet)) ;
2402 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX_t(iDet)) ;
2403 if ( qa->IsSet(AliQA::DETECTORINDEX_t(iDet), AliQA::kSIM, AliQA::kFATAL)) {
2404 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
2406 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
2407 fRunLocalReconstruction.Contains("ALL") ) {
2408 newRunLocalReconstruction += detName ;
2409 newRunLocalReconstruction += " " ;
2411 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
2412 fRunTracking.Contains("ALL") ) {
2413 newRunTracking += detName ;
2414 newRunTracking += " " ;
2416 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
2417 fFillESD.Contains("ALL") ) {
2418 newFillESD += detName ;
2423 fRunLocalReconstruction = newRunLocalReconstruction ;
2424 fRunTracking = newRunTracking ;
2425 fFillESD = newFillESD ;
2428 //_____________________________________________________________________________
2429 Int_t AliReconstruction::GetDetIndex(const char* detector)
2431 // return the detector index corresponding to detector
2433 for (index = 0; index < fgkNDetectors ; index++) {
2434 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2439 //_____________________________________________________________________________
2440 Bool_t AliReconstruction::FinishPlaneEff() {
2442 // Here execute all the necessary operationis, at the end of the tracking phase,
2443 // in case that evaluation of PlaneEfficiencies was required for some detector.
2444 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
2446 // This Preliminary version works only FOR ITS !!!!!
2447 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2450 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2453 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2454 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
2455 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2456 if(fTracker[iDet]) {
2457 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
2458 ret=planeeff->WriteIntoCDB();
2459 if(planeeff->GetCreateHistos()) {
2460 TString name="PlaneEffHisto";
2461 name+=fgkDetectorName[iDet];
2463 ret*=planeeff->WriteHistosToFile(name,"RECREATE");
2469 //_____________________________________________________________________________
2470 Bool_t AliReconstruction::InitPlaneEff() {
2472 // Here execute all the necessary operations, before of the tracking phase,
2473 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
2474 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
2475 // which should be updated/recalculated.
2477 // This Preliminary version will work only FOR ITS !!!!!
2478 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2481 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2483 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));