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>
132 #include "AliReconstruction.h"
133 #include "AliCodeTimer.h"
134 #include "AliReconstructor.h"
136 #include "AliRunLoader.h"
138 #include "AliRawReaderFile.h"
139 #include "AliRawReaderDate.h"
140 #include "AliRawReaderRoot.h"
141 #include "AliRawEventHeaderBase.h"
142 #include "AliESDEvent.h"
143 #include "AliESDMuonTrack.h"
144 #include "AliESDfriend.h"
145 #include "AliESDVertex.h"
146 #include "AliESDcascade.h"
147 #include "AliESDkink.h"
148 #include "AliESDtrack.h"
149 #include "AliESDCaloCluster.h"
150 #include "AliESDCaloCells.h"
151 #include "AliMultiplicity.h"
152 #include "AliTracker.h"
153 #include "AliVertexer.h"
154 #include "AliVertexerTracks.h"
155 #include "AliV0vertexer.h"
156 #include "AliCascadeVertexer.h"
157 #include "AliHeader.h"
158 #include "AliGenEventHeader.h"
160 #include "AliESDpid.h"
161 #include "AliESDtrack.h"
162 #include "AliESDPmdTrack.h"
164 #include "AliESDTagCreator.h"
165 #include "AliAODTagCreator.h"
167 #include "AliGeomManager.h"
168 #include "AliTrackPointArray.h"
169 #include "AliCDBManager.h"
170 #include "AliCDBStorage.h"
171 #include "AliCDBEntry.h"
172 #include "AliAlignObj.h"
174 #include "AliCentralTrigger.h"
175 #include "AliTriggerConfiguration.h"
176 #include "AliTriggerClass.h"
177 #include "AliCTPRawStream.h"
179 #include "AliAODEvent.h"
180 #include "AliAODHeader.h"
181 #include "AliAODTrack.h"
182 #include "AliAODVertex.h"
183 #include "AliAODv0.h"
184 #include "AliAODJet.h"
185 #include "AliAODCaloCells.h"
186 #include "AliAODCaloCluster.h"
187 #include "AliAODPmdCluster.h"
188 #include "AliAODFmdCluster.h"
189 #include "AliAODTracklets.h"
191 #include "AliQADataMakerRec.h"
192 #include "AliGlobalQADataMaker.h"
194 #include "AliQADataMakerSteer.h"
196 #include "AliPlaneEff.h"
198 #include "AliSysInfo.h" // memory snapshots
199 #include "AliRawHLTManager.h"
202 ClassImp(AliReconstruction)
205 //_____________________________________________________________________________
206 const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
208 //_____________________________________________________________________________
209 AliReconstruction::AliReconstruction(const char* gAliceFilename,
210 const char* name, const char* title) :
213 fUniformField(kTRUE),
214 fRunVertexFinder(kTRUE),
215 fRunVertexFinderTracks(kTRUE),
216 fRunHLTTracking(kFALSE),
217 fRunMuonTracking(kFALSE),
219 fRunCascadeFinder(kTRUE),
220 fStopOnError(kFALSE),
221 fWriteAlignmentData(kFALSE),
222 fWriteESDfriend(kFALSE),
224 fFillTriggerESD(kTRUE),
232 fRunLocalReconstruction("ALL"),
235 fUseTrackingErrorsForAlignment(""),
236 fGAliceFileName(gAliceFilename),
241 fNumberOfEventsPerFile(1),
244 fLoadAlignFromCDB(kTRUE),
245 fLoadAlignData("ALL"),
251 fParentRawReader(NULL),
254 fDiamondProfile(NULL),
255 fMeanVertexConstraint(kTRUE),
259 fAlignObjArray(NULL),
262 fInitCDBCalled(kFALSE),
263 fSetRunNumberFromDataCalled(kFALSE),
265 fRunGlobalQA(kFALSE),
270 // create reconstruction object with default parameters
272 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
273 fReconstructor[iDet] = NULL;
274 fLoader[iDet] = NULL;
275 fTracker[iDet] = NULL;
276 fQADataMaker[iDet] = NULL;
277 fQACycles[iDet] = 999999;
279 fQADataMaker[fgkNDetectors]=NULL; //Global QA
283 //_____________________________________________________________________________
284 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
287 fUniformField(rec.fUniformField),
288 fRunVertexFinder(rec.fRunVertexFinder),
289 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
290 fRunHLTTracking(rec.fRunHLTTracking),
291 fRunMuonTracking(rec.fRunMuonTracking),
292 fRunV0Finder(rec.fRunV0Finder),
293 fRunCascadeFinder(rec.fRunCascadeFinder),
294 fStopOnError(rec.fStopOnError),
295 fWriteAlignmentData(rec.fWriteAlignmentData),
296 fWriteESDfriend(rec.fWriteESDfriend),
297 fWriteAOD(rec.fWriteAOD),
298 fFillTriggerESD(rec.fFillTriggerESD),
300 fCleanESD(rec.fCleanESD),
301 fV0DCAmax(rec.fV0DCAmax),
302 fV0CsPmin(rec.fV0CsPmin),
306 fRunLocalReconstruction(rec.fRunLocalReconstruction),
307 fRunTracking(rec.fRunTracking),
308 fFillESD(rec.fFillESD),
309 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
310 fGAliceFileName(rec.fGAliceFileName),
312 fEquipIdMap(rec.fEquipIdMap),
313 fFirstEvent(rec.fFirstEvent),
314 fLastEvent(rec.fLastEvent),
315 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
318 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
319 fLoadAlignData(rec.fLoadAlignData),
320 fESDPar(rec.fESDPar),
326 fDiamondProfile(NULL),
327 fMeanVertexConstraint(rec.fMeanVertexConstraint),
331 fAlignObjArray(rec.fAlignObjArray),
332 fCDBUri(rec.fCDBUri),
334 fInitCDBCalled(rec.fInitCDBCalled),
335 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
337 fRunGlobalQA(rec.fRunGlobalQA),
338 fInLoopQA(rec.fInLoopQA),
339 fRunPlaneEff(rec.fRunPlaneEff)
343 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
344 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
346 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
347 fReconstructor[iDet] = NULL;
348 fLoader[iDet] = NULL;
349 fTracker[iDet] = NULL;
350 fQADataMaker[iDet] = NULL;
351 fQACycles[iDet] = rec.fQACycles[iDet];
353 fQADataMaker[fgkNDetectors]=NULL; //Global QA
354 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
355 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
359 //_____________________________________________________________________________
360 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
362 // assignment operator
364 this->~AliReconstruction();
365 new(this) AliReconstruction(rec);
369 //_____________________________________________________________________________
370 AliReconstruction::~AliReconstruction()
376 fSpecCDBUri.Delete();
378 AliCodeTimer::Instance()->Print();
381 //_____________________________________________________________________________
382 void AliReconstruction::InitCDB()
384 // activate a default CDB storage
385 // First check if we have any CDB storage set, because it is used
386 // to retrieve the calibration and alignment constants
388 if (fInitCDBCalled) return;
389 fInitCDBCalled = kTRUE;
391 AliCDBManager* man = AliCDBManager::Instance();
392 if (man->IsDefaultStorageSet())
394 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
395 AliWarning("Default CDB storage has been already set !");
396 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
397 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
398 fCDBUri = man->GetDefaultStorage()->GetURI();
401 if (fCDBUri.Length() > 0)
403 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
404 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
405 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
407 fCDBUri="local://$ALICE_ROOT";
408 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
409 AliWarning("Default CDB storage not yet set !!!!");
410 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
411 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
414 man->SetDefaultStorage(fCDBUri);
417 // Now activate the detector specific CDB storage locations
418 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
419 TObject* obj = fSpecCDBUri[i];
421 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
422 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
423 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
424 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
429 //_____________________________________________________________________________
430 void AliReconstruction::SetDefaultStorage(const char* uri) {
431 // Store the desired default CDB storage location
432 // Activate it later within the Run() method
438 //_____________________________________________________________________________
439 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
440 // Store a detector-specific CDB storage location
441 // Activate it later within the Run() method
443 AliCDBPath aPath(calibType);
444 if(!aPath.IsValid()){
445 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
446 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
447 if(!strcmp(calibType, fgkDetectorName[iDet])) {
448 aPath.SetPath(Form("%s/*", calibType));
449 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
453 if(!aPath.IsValid()){
454 AliError(Form("Not a valid path or detector: %s", calibType));
459 // // check that calibType refers to a "valid" detector name
460 // Bool_t isDetector = kFALSE;
461 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
462 // TString detName = fgkDetectorName[iDet];
463 // if(aPath.GetLevel0() == detName) {
464 // isDetector = kTRUE;
470 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
474 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
475 if (obj) fSpecCDBUri.Remove(obj);
476 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
480 //_____________________________________________________________________________
481 Bool_t AliReconstruction::SetRunNumberFromData()
483 // The method is called in Run() in order
484 // to set a correct run number.
485 // In case of raw data reconstruction the
486 // run number is taken from the raw data header
488 if (fSetRunNumberFromDataCalled) return kTRUE;
489 fSetRunNumberFromDataCalled = kTRUE;
491 AliCDBManager* man = AliCDBManager::Instance();
493 if(man->GetRun() > 0) {
494 AliWarning("Run number is taken from event header! Ignoring settings in AliCDBManager!");
498 AliError("No run loader is found !");
501 // read run number from gAlice
502 if(fRunLoader->GetAliRun())
503 AliCDBManager::Instance()->SetRun(fRunLoader->GetHeader()->GetRun());
506 if(fRawReader->NextEvent()) {
507 AliCDBManager::Instance()->SetRun(fRawReader->GetRunNumber());
508 fRawReader->RewindEvents();
511 AliError("No raw-data events found !");
516 AliError("Neither gAlice nor RawReader objects are found !");
526 //_____________________________________________________________________________
527 void AliReconstruction::SetCDBLock() {
528 // Set CDB lock: from now on it is forbidden to reset the run number
529 // or the default storage or to activate any further storage!
531 AliCDBManager::Instance()->SetLock(1);
534 //_____________________________________________________________________________
535 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
537 // Read the alignment objects from CDB.
538 // Each detector is supposed to have the
539 // alignment objects in DET/Align/Data CDB path.
540 // All the detector objects are then collected,
541 // sorted by geometry level (starting from ALIC) and
542 // then applied to the TGeo geometry.
543 // Finally an overlaps check is performed.
545 // Load alignment data from CDB and fill fAlignObjArray
546 if(fLoadAlignFromCDB){
548 TString detStr = detectors;
549 TString loadAlObjsListOfDets = "";
551 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
552 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
553 loadAlObjsListOfDets += fgkDetectorName[iDet];
554 loadAlObjsListOfDets += " ";
555 } // end loop over detectors
556 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
557 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
559 // Check if the array with alignment objects was
560 // provided by the user. If yes, apply the objects
561 // to the present TGeo geometry
562 if (fAlignObjArray) {
563 if (gGeoManager && gGeoManager->IsClosed()) {
564 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
565 AliError("The misalignment of one or more volumes failed!"
566 "Compare the list of simulated detectors and the list of detector alignment data!");
571 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
577 delete fAlignObjArray; fAlignObjArray=0;
582 //_____________________________________________________________________________
583 void AliReconstruction::SetGAliceFile(const char* fileName)
585 // set the name of the galice file
587 fGAliceFileName = fileName;
590 //_____________________________________________________________________________
591 void AliReconstruction::SetOption(const char* detector, const char* option)
593 // set options for the reconstruction of a detector
595 TObject* obj = fOptions.FindObject(detector);
596 if (obj) fOptions.Remove(obj);
597 fOptions.Add(new TNamed(detector, option));
601 //_____________________________________________________________________________
602 Bool_t AliReconstruction::Run(const char* input, Bool_t IsOnline)
604 // run the reconstruction
610 if (!input) input = fInput.Data();
611 TString fileName(input);
612 if (fileName.EndsWith("/")) {
613 fRawReader = new AliRawReaderFile(fileName);
614 } else if (fileName.EndsWith(".root")) {
615 fRawReader = new AliRawReaderRoot(fileName);
616 } else if (!fileName.IsNull()) {
617 fRawReader = new AliRawReaderDate(fileName);
622 AliError("Null pointer to the event structure!");
625 fRawReader = new AliRawReaderDate((void *)input);
628 if (!fEquipIdMap.IsNull() && fRawReader)
629 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
631 if (!fUseHLTData.IsNull()) {
632 // create the RawReaderHLT which performs redirection of HLT input data for
633 // the specified detectors
634 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
636 fParentRawReader=fRawReader;
637 fRawReader=pRawReader;
639 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
643 AliSysInfo::AddStamp("Start");
644 // get the run loader
645 if (!InitRunLoader()) return kFALSE;
646 AliSysInfo::AddStamp("LoadLoader");
648 // Initialize the CDB storage
651 AliSysInfo::AddStamp("LoadCDB");
653 // Set run number in CDBManager (if it is not already set by the user)
654 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
656 // Set CDB lock: from now on it is forbidden to reset the run number
657 // or the default storage or to activate any further storage!
660 // Import ideal TGeo geometry and apply misalignment
662 TString geom(gSystem->DirName(fGAliceFileName));
663 geom += "/geometry.root";
664 AliGeomManager::LoadGeometry(geom.Data());
665 if (!gGeoManager) if (fStopOnError) return kFALSE;
668 if (!MisalignGeometry(fLoadAlignData)) if (fStopOnError) return kFALSE;
669 AliSysInfo::AddStamp("LoadGeom");
672 AliQADataMakerSteer qas ;
673 if (fRunQA && fRawReader) qas.Run(fRunLocalReconstruction, fRawReader) ;
674 // checking the QA of previous steps
678 // local reconstruction
679 if (!fRunLocalReconstruction.IsNull()) {
680 if (!RunLocalReconstruction(fRunLocalReconstruction)) {
681 if (fStopOnError) {CleanUp(); return kFALSE;}
687 if (fRunVertexFinder && !CreateVertexer()) {
693 AliSysInfo::AddStamp("Vertexer");
696 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
702 AliSysInfo::AddStamp("LoadTrackers");
704 // get the possibly already existing ESD file and tree
705 AliESDEvent* esd = new AliESDEvent(); AliESDEvent* hltesd = new AliESDEvent();
706 TFile* fileOld = NULL;
707 TTree* treeOld = NULL; TTree *hlttreeOld = NULL;
708 if (!gSystem->AccessPathName("AliESDs.root")){
709 gSystem->CopyFile("AliESDs.root", "AliESDs.old.root", kTRUE);
710 fileOld = TFile::Open("AliESDs.old.root");
711 if (fileOld && fileOld->IsOpen()) {
712 treeOld = (TTree*) fileOld->Get("esdTree");
713 if (treeOld)esd->ReadFromTree(treeOld);
714 hlttreeOld = (TTree*) fileOld->Get("HLTesdTree");
715 if (hlttreeOld) hltesd->ReadFromTree(hlttreeOld);
719 // create the ESD output file and tree
720 TFile* file = TFile::Open("AliESDs.root", "RECREATE");
721 file->SetCompressionLevel(2);
722 if (!file->IsOpen()) {
723 AliError("opening AliESDs.root failed");
724 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
727 TTree* tree = new TTree("esdTree", "Tree with ESD objects");
728 esd = new AliESDEvent();
729 esd->CreateStdContent();
730 esd->WriteToTree(tree);
732 TTree* hlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
733 hltesd = new AliESDEvent();
734 hltesd->CreateStdContent();
735 hltesd->WriteToTree(hlttree);
738 delete esd; delete hltesd;
739 esd = NULL; hltesd = NULL;
741 // create the branch with ESD additions
745 AliESDfriend *esdf = 0;
746 if (fWriteESDfriend) {
747 esdf = new AliESDfriend();
748 TBranch *br=tree->Branch("ESDfriend.","AliESDfriend", &esdf);
749 br->SetFile("AliESDfriends.root");
750 esd->AddObject(esdf);
754 // Get the GRP CDB entry
755 AliCDBEntry* entryGRP = AliCDBManager::Instance()->Get("GRP/GRP/Data");
758 fGRPList = dynamic_cast<TList*> (entryGRP->GetObject());
760 AliError("No GRP entry found in OCDB!");
763 // Get the diamond profile from OCDB
764 AliCDBEntry* entry = AliCDBManager::Instance()
765 ->Get("GRP/Calib/MeanVertex");
768 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
770 AliError("No diamond profile found in OCDB!");
773 AliVertexerTracks tVertexer(AliTracker::GetBz());
774 if(fDiamondProfile && fMeanVertexConstraint) tVertexer.SetVtxStart(fDiamondProfile);
776 if (fRawReader) fRawReader->RewindEvents();
779 gSystem->GetProcInfo(&ProcInfo);
780 AliInfo(Form("Current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual));
783 //Initialize the QA and start of cycle for out-of-cycle QA
785 TString detStr(fFillESD);
786 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
787 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
788 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
790 AliInfo(Form("Initializing the QA data maker for %s",
791 fgkDetectorName[iDet]));
792 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
793 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
795 qadm->StartOfCycle(AliQA::kRECPOINTS);
796 qadm->StartOfCycle(AliQA::kESDS,"same");
800 AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
801 AliInfo(Form("Initializing the global QA data maker"));
803 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
804 AliTracker::SetResidualsArray(arr);
805 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
807 qadm->StartOfCycle(AliQA::kRECPOINTS);
808 qadm->StartOfCycle(AliQA::kESDS,"same");
813 //Initialize the Plane Efficiency framework
814 if (fRunPlaneEff && !InitPlaneEff()) {
815 if(fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
818 //******* The loop over events
819 for (Int_t iEvent = 0; iEvent < fRunLoader->GetNumberOfEvents(); iEvent++) {
820 if (fRawReader) fRawReader->NextEvent();
821 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
822 // copy old ESD to the new one
824 esd->ReadFromTree(treeOld);
825 treeOld->GetEntry(iEvent);
829 esd->ReadFromTree(hlttreeOld);
830 hlttreeOld->GetEntry(iEvent);
836 AliInfo(Form("processing event %d", iEvent));
838 //Start of cycle for the in-loop QA
839 if (fRunQA && fInLoopQA) {
840 TString detStr(fFillESD);
841 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
842 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
843 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
845 qadm->StartOfCycle(AliQA::kRECPOINTS);
846 qadm->StartOfCycle(AliQA::kESDS, "same") ;
849 AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
850 qadm->StartOfCycle(AliQA::kRECPOINTS);
851 qadm->StartOfCycle(AliQA::kESDS,"same");
855 fRunLoader->GetEvent(iEvent);
858 sprintf(aFileName, "ESD_%d.%d_final.root",
859 fRunLoader->GetHeader()->GetRun(),
860 fRunLoader->GetHeader()->GetEventNrInRun());
861 if (!gSystem->AccessPathName(aFileName)) continue;
863 // local signle event reconstruction
864 if (!fRunLocalReconstruction.IsNull()) {
865 if (!RunLocalEventReconstruction(fRunLocalReconstruction)) {
866 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
870 esd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
871 hltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
872 esd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
873 hltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
875 // Set magnetic field from the tracker
876 esd->SetMagneticField(AliTracker::GetBz());
877 hltesd->SetMagneticField(AliTracker::GetBz());
881 // Fill raw-data error log into the ESD
882 if (fRawReader) FillRawDataErrorLog(iEvent,esd);
885 if (fRunVertexFinder) {
886 if (!ReadESD(esd, "vertex")) {
887 if (!RunVertexFinder(esd)) {
888 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
890 if (fCheckPointLevel > 0) WriteESD(esd, "vertex");
895 if (!fRunTracking.IsNull()) {
896 if (fRunHLTTracking) {
897 hltesd->SetPrimaryVertexSPD(esd->GetVertex());
898 if (!RunHLTTracking(hltesd)) {
899 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
905 if (!fRunTracking.IsNull()) {
906 if (fRunMuonTracking) {
907 if (!RunMuonTracking(esd)) {
908 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
914 if (!fRunTracking.IsNull()) {
915 if (!ReadESD(esd, "tracking")) {
916 if (!RunTracking(esd)) {
917 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
919 if (fCheckPointLevel > 0) WriteESD(esd, "tracking");
924 if (!fFillESD.IsNull()) {
925 if (!FillESD(esd, fFillESD)) {
926 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
930 // fill Event header information from the RawEventHeader
931 if (fRawReader){FillRawEventHeaderESD(esd);}
934 AliESDpid::MakePID(esd);
935 if (fCheckPointLevel > 1) WriteESD(esd, "PID");
937 if (fFillTriggerESD) {
938 if (!ReadESD(esd, "trigger")) {
939 if (!FillTriggerESD(esd)) {
940 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
942 if (fCheckPointLevel > 1) WriteESD(esd, "trigger");
948 //Try to improve the reconstructed primary vertex position using the tracks
949 AliESDVertex *pvtx=0;
950 TObject* obj = fOptions.FindObject("ITS");
952 TString optITS = obj->GetTitle();
953 if (optITS.Contains("cosmics") || optITS.Contains("COSMICS"))
954 fRunVertexFinderTracks=kFALSE;
956 if(fRunVertexFinderTracks) pvtx=tVertexer.FindPrimaryVertex(esd);
957 if(fDiamondProfile) esd->SetDiamond(fDiamondProfile);
960 if (pvtx->GetStatus()) {
961 // Store the improved primary vertex
962 esd->SetPrimaryVertex(pvtx);
963 // Propagate the tracks to the DCA to the improved primary vertex
964 Double_t somethingbig = 777.;
965 Double_t bz = esd->GetMagneticField();
966 Int_t nt=esd->GetNumberOfTracks();
968 AliESDtrack *t = esd->GetTrack(nt);
969 t->RelateToVertex(pvtx, bz, somethingbig);
976 vtxer.Tracks2V0vertices(esd);
978 if (fRunCascadeFinder) {
980 AliCascadeVertexer cvtxer;
981 cvtxer.V0sTracks2CascadeVertices(esd);
986 if (fCleanESD) CleanESD(esd);
988 if (fRunQA && fRunGlobalQA) {
989 AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
990 if (qadm) qadm->Exec(AliQA::kESDS, esd);
993 if (fWriteESDfriend) {
994 esdf->~AliESDfriend();
995 new (esdf) AliESDfriend(); // Reset...
996 esd->GetESDfriend(esdf);
1003 if (fCheckPointLevel > 0) WriteESD(esd, "final");
1006 if (fWriteESDfriend) {
1007 esdf->~AliESDfriend();
1008 new (esdf) AliESDfriend(); // Reset...
1011 gSystem->GetProcInfo(&ProcInfo);
1012 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1015 // End of cycle for the in-loop QA
1016 if (fRunQA && fInLoopQA) {
1017 RunQA(fFillESD.Data(), esd);
1018 TString detStr(fFillESD);
1019 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1020 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1021 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1022 if (!qadm) continue;
1023 qadm->EndOfCycle(AliQA::kRECPOINTS);
1024 qadm->EndOfCycle(AliQA::kESDS);
1028 AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
1030 qadm->EndOfCycle(AliQA::kRECPOINTS);
1031 qadm->EndOfCycle(AliQA::kESDS);
1037 //******** End of the loop over events
1041 tree->GetUserInfo()->Add(esd);
1042 hlttree->GetUserInfo()->Add(hltesd);
1044 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1045 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1047 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1048 cdbMapCopy->SetOwner(1);
1049 cdbMapCopy->SetName("cdbMap");
1050 TIter iter(cdbMap->GetTable());
1053 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1054 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1055 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1056 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1059 TList *cdbListCopy = new TList();
1060 cdbListCopy->SetOwner(1);
1061 cdbListCopy->SetName("cdbList");
1063 TIter iter2(cdbList);
1066 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1067 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1070 tree->GetUserInfo()->Add(cdbMapCopy);
1071 tree->GetUserInfo()->Add(cdbListCopy);
1074 if(fESDPar.Contains("ESD.par")){
1075 AliInfo("Attaching ESD.par to Tree");
1076 TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1077 tree->GetUserInfo()->Add(fn);
1083 if (fWriteESDfriend)
1084 tree->SetBranchStatus("ESDfriend*",0);
1085 // we want to have only one tree version number
1086 tree->Write(tree->GetName(),TObject::kOverwrite);
1089 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1090 if (fRunPlaneEff && !FinishPlaneEff()) {
1091 AliWarning("Finish PlaneEff evaluation failed");
1095 CleanUp(file, fileOld);
1098 TFile *esdFile = TFile::Open("AliESDs.root", "READONLY");
1099 TFile *aodFile = TFile::Open("AliAOD.root", "RECREATE");
1100 ESDFile2AODFile(esdFile, aodFile);
1105 // Create tags for the events in the ESD tree (the ESD tree is always present)
1106 // In case of empty events the tags will contain dummy values
1107 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1108 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPList);
1110 AliAODTagCreator *aodtagCreator = new AliAODTagCreator();
1111 aodtagCreator->CreateAODTags(fFirstEvent,fLastEvent,fGRPList);
1114 //Finish QA and end of cycle for out-of-loop QA
1115 if (fRunQA && !fInLoopQA) {
1116 qas.Run(fRunLocalReconstruction.Data(), AliQA::kRECPOINTS);
1118 qas.Run(fRunTracking.Data(), AliQA::kESDS);
1121 AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
1123 qadm->EndOfCycle(AliQA::kRECPOINTS);
1124 qadm->EndOfCycle(AliQA::kESDS);
1130 // Cleanup of CDB manager: cache and active storages!
1131 AliCDBManager::Instance()->ClearCache();
1138 //_____________________________________________________________________________
1139 Bool_t AliReconstruction::RunLocalReconstruction(const TString& /*detectors*/)
1141 // run the local reconstruction
1142 static Int_t eventNr=0;
1143 AliCodeTimerAuto("")
1145 // AliCDBManager* man = AliCDBManager::Instance();
1146 // Bool_t origCache = man->GetCacheFlag();
1148 // TString detStr = detectors;
1149 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1150 // if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1151 // AliReconstructor* reconstructor = GetReconstructor(iDet);
1152 // if (!reconstructor) continue;
1153 // if (reconstructor->HasLocalReconstruction()) continue;
1155 // AliCodeTimerStart(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1156 // AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1158 // AliCodeTimerStart(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1159 // AliInfo(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1161 // man->SetCacheFlag(kTRUE);
1162 // TString calibPath = Form("%s/Calib/*", fgkDetectorName[iDet]);
1163 // man->GetAll(calibPath); // entries are cached!
1165 // AliCodeTimerStop(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1167 // if (fRawReader) {
1168 // fRawReader->RewindEvents();
1169 // reconstructor->Reconstruct(fRunLoader, fRawReader);
1171 // reconstructor->Reconstruct(fRunLoader);
1174 // AliCodeTimerStop(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1175 // AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr));
1177 // // unload calibration data
1178 // man->UnloadFromCache(calibPath);
1179 // //man->ClearCache();
1182 // man->SetCacheFlag(origCache);
1184 // if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1185 // AliError(Form("the following detectors were not found: %s",
1187 // if (fStopOnError) return kFALSE;
1194 //_____________________________________________________________________________
1195 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1197 // run the local reconstruction
1199 static Int_t eventNr=0;
1200 AliCodeTimerAuto("")
1202 TString detStr = detectors;
1203 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1204 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1205 AliReconstructor* reconstructor = GetReconstructor(iDet);
1206 if (!reconstructor) continue;
1207 AliLoader* loader = fLoader[iDet];
1209 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1212 // conversion of digits
1213 if (fRawReader && reconstructor->HasDigitConversion()) {
1214 AliInfo(Form("converting raw data digits into root objects for %s",
1215 fgkDetectorName[iDet]));
1216 AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1217 fgkDetectorName[iDet]));
1218 loader->LoadDigits("update");
1219 loader->CleanDigits();
1220 loader->MakeDigitsContainer();
1221 TTree* digitsTree = loader->TreeD();
1222 reconstructor->ConvertDigits(fRawReader, digitsTree);
1223 loader->WriteDigits("OVERWRITE");
1224 loader->UnloadDigits();
1226 // local reconstruction
1227 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1228 AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1229 loader->LoadRecPoints("update");
1230 loader->CleanRecPoints();
1231 loader->MakeRecPointsContainer();
1232 TTree* clustersTree = loader->TreeR();
1233 if (fRawReader && !reconstructor->HasDigitConversion()) {
1234 reconstructor->Reconstruct(fRawReader, clustersTree);
1236 loader->LoadDigits("read");
1237 TTree* digitsTree = loader->TreeD();
1239 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1240 if (fStopOnError) return kFALSE;
1242 reconstructor->Reconstruct(digitsTree, clustersTree);
1244 loader->UnloadDigits();
1247 // In-loop QA for local reconstrucion
1248 if (fRunQA && fInLoopQA) {
1249 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1252 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1254 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1256 qadm->Exec(AliQA::kRECPOINTS, clustersTree) ;
1259 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1263 loader->WriteRecPoints("OVERWRITE");
1264 loader->UnloadRecPoints();
1265 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1268 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1269 AliError(Form("the following detectors were not found: %s",
1271 if (fStopOnError) return kFALSE;
1277 //_____________________________________________________________________________
1278 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1280 // run the barrel tracking
1282 AliCodeTimerAuto("")
1284 AliESDVertex* vertex = NULL;
1285 Double_t vtxPos[3] = {0, 0, 0};
1286 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1287 TArrayF mcVertex(3);
1288 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1289 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1290 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1294 if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1295 AliInfo("running the ITS vertex finder");
1296 if (fLoader[0]) fLoader[0]->LoadRecPoints();
1297 vertex = fVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber());
1298 if (fLoader[0]) fLoader[0]->UnloadRecPoints();
1300 AliWarning("Vertex not found");
1301 vertex = new AliESDVertex();
1302 vertex->SetName("default");
1305 vertex->SetName("reconstructed");
1309 AliInfo("getting the primary vertex from MC");
1310 vertex = new AliESDVertex(vtxPos, vtxErr);
1314 vertex->GetXYZ(vtxPos);
1315 vertex->GetSigmaXYZ(vtxErr);
1317 AliWarning("no vertex reconstructed");
1318 vertex = new AliESDVertex(vtxPos, vtxErr);
1320 esd->SetPrimaryVertexSPD(vertex);
1321 // if SPD multiplicity has been determined, it is stored in the ESD
1322 AliMultiplicity *mult = fVertexer->GetMultiplicity();
1323 if(mult)esd->SetMultiplicity(mult);
1325 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1326 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1333 //_____________________________________________________________________________
1334 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
1336 // run the HLT barrel tracking
1338 AliCodeTimerAuto("")
1341 AliError("Missing runLoader!");
1345 AliInfo("running HLT tracking");
1347 // Get a pointer to the HLT reconstructor
1348 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
1349 if (!reconstructor) return kFALSE;
1352 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1353 TString detName = fgkDetectorName[iDet];
1354 AliDebug(1, Form("%s HLT tracking", detName.Data()));
1355 reconstructor->SetOption(detName.Data());
1356 AliTracker *tracker = reconstructor->CreateTracker();
1358 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
1359 if (fStopOnError) return kFALSE;
1363 Double_t vtxErr[3]={0.005,0.005,0.010};
1364 const AliESDVertex *vertex = esd->GetVertex();
1365 vertex->GetXYZ(vtxPos);
1366 tracker->SetVertex(vtxPos,vtxErr);
1368 fLoader[iDet]->LoadRecPoints("read");
1369 TTree* tree = fLoader[iDet]->TreeR();
1371 AliError(Form("Can't get the %s cluster tree", detName.Data()));
1374 tracker->LoadClusters(tree);
1376 if (tracker->Clusters2Tracks(esd) != 0) {
1377 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
1381 tracker->UnloadClusters();
1389 //_____________________________________________________________________________
1390 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
1392 // run the muon spectrometer tracking
1394 AliCodeTimerAuto("")
1397 AliError("Missing runLoader!");
1400 Int_t iDet = 7; // for MUON
1402 AliInfo("is running...");
1404 // Get a pointer to the MUON reconstructor
1405 AliReconstructor *reconstructor = GetReconstructor(iDet);
1406 if (!reconstructor) return kFALSE;
1409 TString detName = fgkDetectorName[iDet];
1410 AliDebug(1, Form("%s tracking", detName.Data()));
1411 AliTracker *tracker = reconstructor->CreateTracker();
1413 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1418 fLoader[iDet]->LoadRecPoints("read");
1420 tracker->LoadClusters(fLoader[iDet]->TreeR());
1422 Int_t rv = tracker->Clusters2Tracks(esd);
1426 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1430 fLoader[iDet]->UnloadRecPoints();
1432 tracker->UnloadClusters();
1440 //_____________________________________________________________________________
1441 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
1443 // run the barrel tracking
1444 static Int_t eventNr=0;
1445 AliCodeTimerAuto("")
1447 AliInfo("running tracking");
1449 //Fill the ESD with the T0 info (will be used by the TOF)
1450 if (fReconstructor[11] && fLoader[11]) {
1451 fLoader[11]->LoadRecPoints("READ");
1452 TTree *treeR = fLoader[11]->TreeR();
1453 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
1456 // pass 1: TPC + ITS inwards
1457 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1458 if (!fTracker[iDet]) continue;
1459 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
1462 fLoader[iDet]->LoadRecPoints("read");
1463 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
1464 TTree* tree = fLoader[iDet]->TreeR();
1466 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1469 fTracker[iDet]->LoadClusters(tree);
1470 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1472 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
1473 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1476 if (fCheckPointLevel > 1) {
1477 WriteESD(esd, Form("%s.tracking", fgkDetectorName[iDet]));
1479 // preliminary PID in TPC needed by the ITS tracker
1481 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1482 AliESDpid::MakePID(esd);
1484 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
1487 // pass 2: ALL backwards
1489 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1490 if (!fTracker[iDet]) continue;
1491 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
1494 if (iDet > 1) { // all except ITS, TPC
1496 fLoader[iDet]->LoadRecPoints("read");
1497 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
1498 tree = fLoader[iDet]->TreeR();
1500 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1503 fTracker[iDet]->LoadClusters(tree);
1504 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1508 if (fTracker[iDet]->PropagateBack(esd) != 0) {
1509 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
1512 if (fCheckPointLevel > 1) {
1513 WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
1517 if (iDet > 2) { // all except ITS, TPC, TRD
1518 fTracker[iDet]->UnloadClusters();
1519 fLoader[iDet]->UnloadRecPoints();
1521 // updated PID in TPC needed by the ITS tracker -MI
1523 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1524 AliESDpid::MakePID(esd);
1526 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1529 // write space-points to the ESD in case alignment data output
1531 if (fWriteAlignmentData)
1532 WriteAlignmentData(esd);
1534 // pass 3: TRD + TPC + ITS refit inwards
1536 if (fRunQA && fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1538 for (Int_t iDet = 2; iDet >= 0; iDet--) {
1539 if (!fTracker[iDet]) continue;
1540 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
1543 if (fTracker[iDet]->RefitInward(esd) != 0) {
1544 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
1547 // run postprocessing
1548 if (fTracker[iDet]->PostProcess(esd) != 0) {
1549 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
1552 if (fCheckPointLevel > 1) {
1553 WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
1555 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1557 fTracker[iDet]->UnloadClusters();
1558 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
1559 fLoader[iDet]->UnloadRecPoints();
1560 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
1563 if (fRunQA && fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1566 // Propagate track to the vertex - if not done by ITS
1568 Int_t ntracks = esd->GetNumberOfTracks();
1569 for (Int_t itrack=0; itrack<ntracks; itrack++){
1570 const Double_t kRadius = 3; // beam pipe radius
1571 const Double_t kMaxStep = 5; // max step
1572 const Double_t kMaxD = 123456; // max distance to prim vertex
1573 Double_t fieldZ = AliTracker::GetBz(); //
1574 AliESDtrack * track = esd->GetTrack(itrack);
1575 if (!track) continue;
1576 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1577 AliTracker::PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kTRUE);
1578 track->RelateToVertex(esd->GetVertex(),fieldZ, kMaxD);
1584 //_____________________________________________________________________________
1585 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
1587 // Remove the data which are not needed for the physics analysis.
1590 Int_t nTracks=esd->GetNumberOfTracks();
1591 Int_t nV0s=esd->GetNumberOfV0s();
1593 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
1595 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
1596 Bool_t rc=esd->Clean(cleanPars);
1598 nTracks=esd->GetNumberOfTracks();
1599 nV0s=esd->GetNumberOfV0s();
1601 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
1606 //_____________________________________________________________________________
1607 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
1609 // fill the event summary data
1611 AliCodeTimerAuto("")
1612 static Int_t eventNr=0;
1613 TString detStr = detectors;
1615 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1616 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1617 AliReconstructor* reconstructor = GetReconstructor(iDet);
1618 if (!reconstructor) continue;
1619 if (!ReadESD(esd, fgkDetectorName[iDet])) {
1620 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
1621 TTree* clustersTree = NULL;
1622 if (fLoader[iDet]) {
1623 fLoader[iDet]->LoadRecPoints("read");
1624 clustersTree = fLoader[iDet]->TreeR();
1625 if (!clustersTree) {
1626 AliError(Form("Can't get the %s clusters tree",
1627 fgkDetectorName[iDet]));
1628 if (fStopOnError) return kFALSE;
1631 if (fRawReader && !reconstructor->HasDigitConversion()) {
1632 reconstructor->FillESD(fRawReader, clustersTree, esd);
1634 TTree* digitsTree = NULL;
1635 if (fLoader[iDet]) {
1636 fLoader[iDet]->LoadDigits("read");
1637 digitsTree = fLoader[iDet]->TreeD();
1639 AliError(Form("Can't get the %s digits tree",
1640 fgkDetectorName[iDet]));
1641 if (fStopOnError) return kFALSE;
1644 reconstructor->FillESD(digitsTree, clustersTree, esd);
1645 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
1647 if (fLoader[iDet]) {
1648 fLoader[iDet]->UnloadRecPoints();
1651 if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
1655 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1656 AliError(Form("the following detectors were not found: %s",
1658 if (fStopOnError) return kFALSE;
1660 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
1665 //_____________________________________________________________________________
1666 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
1668 // Reads the trigger decision which is
1669 // stored in Trigger.root file and fills
1670 // the corresponding esd entries
1672 AliCodeTimerAuto("")
1674 AliInfo("Filling trigger information into the ESD");
1676 AliCentralTrigger *aCTP = NULL;
1679 AliCTPRawStream input(fRawReader);
1680 if (!input.Next()) {
1681 AliError("No valid CTP (trigger) DDL raw data is found ! The trigger information is not stored in the ESD !");
1684 esd->SetTriggerMask(input.GetClassMask());
1685 esd->SetTriggerCluster(input.GetClusterMask());
1687 aCTP = new AliCentralTrigger();
1688 TString configstr("");
1689 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
1690 AliError("No trigger configuration found in OCDB! The trigger classes information will no be stored in ESD!");
1695 AliRunLoader *runloader = AliRunLoader::GetRunLoader();
1697 if (!runloader->LoadTrigger()) {
1698 aCTP = runloader->GetTrigger();
1699 esd->SetTriggerMask(aCTP->GetClassMask());
1700 esd->SetTriggerCluster(aCTP->GetClusterMask());
1703 AliWarning("No trigger can be loaded! The trigger information is not stored in the ESD !");
1708 AliError("No run loader is available! The trigger information is not stored in the ESD !");
1713 // Now fill the trigger class names into AliESDRun object
1714 AliTriggerConfiguration *config = aCTP->GetConfiguration();
1716 AliError("No trigger configuration has been found! The trigger classes information will no be stored in ESD!");
1720 const TObjArray& classesArray = config->GetClasses();
1721 Int_t nclasses = classesArray.GetEntriesFast();
1722 for( Int_t j=0; j<nclasses; j++ ) {
1723 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At( j );
1724 Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
1725 esd->SetTriggerClass(trclass->GetName(),trindex);
1735 //_____________________________________________________________________________
1736 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
1739 // Filling information from RawReader Header
1742 AliInfo("Filling information from RawReader Header");
1743 esd->SetBunchCrossNumber(0);
1744 esd->SetOrbitNumber(0);
1745 esd->SetPeriodNumber(0);
1746 esd->SetTimeStamp(0);
1747 esd->SetEventType(0);
1748 const AliRawEventHeaderBase * eventHeader = fRawReader->GetEventHeader();
1751 const UInt_t *id = eventHeader->GetP("Id");
1752 esd->SetBunchCrossNumber((id)[1]&0x00000fff);
1753 esd->SetOrbitNumber((((id)[0]<<20)&0xf00000)|(((id)[1]>>12)&0xfffff));
1754 esd->SetPeriodNumber(((id)[0]>>4)&0x0fffffff);
1756 esd->SetTimeStamp((eventHeader->Get("Timestamp")));
1757 esd->SetEventType((eventHeader->Get("Type")));
1764 //_____________________________________________________________________________
1765 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
1767 // check whether detName is contained in detectors
1768 // if yes, it is removed from detectors
1770 // check if all detectors are selected
1771 if ((detectors.CompareTo("ALL") == 0) ||
1772 detectors.BeginsWith("ALL ") ||
1773 detectors.EndsWith(" ALL") ||
1774 detectors.Contains(" ALL ")) {
1779 // search for the given detector
1780 Bool_t result = kFALSE;
1781 if ((detectors.CompareTo(detName) == 0) ||
1782 detectors.BeginsWith(detName+" ") ||
1783 detectors.EndsWith(" "+detName) ||
1784 detectors.Contains(" "+detName+" ")) {
1785 detectors.ReplaceAll(detName, "");
1789 // clean up the detectors string
1790 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1791 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1792 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1797 //_____________________________________________________________________________
1798 Bool_t AliReconstruction::InitRunLoader()
1800 // get or create the run loader
1802 if (gAlice) delete gAlice;
1805 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
1806 // load all base libraries to get the loader classes
1807 TString libs = gSystem->GetLibraries();
1808 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1809 TString detName = fgkDetectorName[iDet];
1810 if (detName == "HLT") continue;
1811 if (libs.Contains("lib" + detName + "base.so")) continue;
1812 gSystem->Load("lib" + detName + "base.so");
1814 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
1816 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1820 fRunLoader->CdGAFile();
1821 if (gFile->GetKey(AliRunLoader::GetGAliceName())) {
1822 if (fRunLoader->LoadgAlice() == 0) {
1823 gAlice = fRunLoader->GetAliRun();
1824 AliTracker::SetFieldMap(gAlice->Field(),fUniformField);
1827 if (!gAlice && !fRawReader) {
1828 AliError(Form("no gAlice object found in file %s",
1829 fGAliceFileName.Data()));
1834 //PH This is a temporary fix to give access to the kinematics
1835 //PH that is needed for the labels of ITS clusters
1836 fRunLoader->LoadHeader();
1837 fRunLoader->LoadKinematics();
1839 } else { // galice.root does not exist
1841 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
1845 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
1846 AliConfig::GetDefaultEventFolderName(),
1849 AliError(Form("could not create run loader in file %s",
1850 fGAliceFileName.Data()));
1854 fRunLoader->MakeTree("E");
1856 while (fRawReader->NextEvent()) {
1857 fRunLoader->SetEventNumber(iEvent);
1858 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1860 fRunLoader->MakeTree("H");
1861 fRunLoader->TreeE()->Fill();
1864 fRawReader->RewindEvents();
1865 if (fNumberOfEventsPerFile > 0)
1866 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
1868 fRunLoader->SetNumberOfEventsPerFile(iEvent);
1869 fRunLoader->WriteHeader("OVERWRITE");
1870 fRunLoader->CdGAFile();
1871 fRunLoader->Write(0, TObject::kOverwrite);
1872 // AliTracker::SetFieldMap(???);
1878 //_____________________________________________________________________________
1879 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
1881 // get the reconstructor object and the loader for a detector
1883 if (fReconstructor[iDet]) return fReconstructor[iDet];
1885 // load the reconstructor object
1886 TPluginManager* pluginManager = gROOT->GetPluginManager();
1887 TString detName = fgkDetectorName[iDet];
1888 TString recName = "Ali" + detName + "Reconstructor";
1889 if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT")) return NULL;
1891 AliReconstructor* reconstructor = NULL;
1892 // first check if a plugin is defined for the reconstructor
1893 TPluginHandler* pluginHandler =
1894 pluginManager->FindHandler("AliReconstructor", detName);
1895 // if not, add a plugin for it
1896 if (!pluginHandler) {
1897 AliDebug(1, Form("defining plugin for %s", recName.Data()));
1898 TString libs = gSystem->GetLibraries();
1899 if (libs.Contains("lib" + detName + "base.so") ||
1900 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
1901 pluginManager->AddHandler("AliReconstructor", detName,
1902 recName, detName + "rec", recName + "()");
1904 pluginManager->AddHandler("AliReconstructor", detName,
1905 recName, detName, recName + "()");
1907 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
1909 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
1910 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
1912 if (reconstructor) {
1913 TObject* obj = fOptions.FindObject(detName.Data());
1914 if (obj) reconstructor->SetOption(obj->GetTitle());
1915 reconstructor->Init();
1916 fReconstructor[iDet] = reconstructor;
1919 // get or create the loader
1920 if (detName != "HLT") {
1921 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
1922 if (!fLoader[iDet]) {
1923 AliConfig::Instance()
1924 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
1926 // first check if a plugin is defined for the loader
1928 pluginManager->FindHandler("AliLoader", detName);
1929 // if not, add a plugin for it
1930 if (!pluginHandler) {
1931 TString loaderName = "Ali" + detName + "Loader";
1932 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
1933 pluginManager->AddHandler("AliLoader", detName,
1934 loaderName, detName + "base",
1935 loaderName + "(const char*, TFolder*)");
1936 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
1938 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
1940 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
1941 fRunLoader->GetEventFolder());
1943 if (!fLoader[iDet]) { // use default loader
1944 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
1946 if (!fLoader[iDet]) {
1947 AliWarning(Form("couldn't get loader for %s", detName.Data()));
1948 if (fStopOnError) return NULL;
1950 fRunLoader->AddLoader(fLoader[iDet]);
1951 fRunLoader->CdGAFile();
1952 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
1953 fRunLoader->Write(0, TObject::kOverwrite);
1958 return reconstructor;
1961 //_____________________________________________________________________________
1962 Bool_t AliReconstruction::CreateVertexer()
1964 // create the vertexer
1967 AliReconstructor* itsReconstructor = GetReconstructor(0);
1968 if (itsReconstructor) {
1969 fVertexer = itsReconstructor->CreateVertexer();
1972 AliWarning("couldn't create a vertexer for ITS");
1973 if (fStopOnError) return kFALSE;
1979 //_____________________________________________________________________________
1980 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
1982 // create the trackers
1984 TString detStr = detectors;
1985 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1986 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1987 AliReconstructor* reconstructor = GetReconstructor(iDet);
1988 if (!reconstructor) continue;
1989 TString detName = fgkDetectorName[iDet];
1990 if (detName == "HLT") {
1991 fRunHLTTracking = kTRUE;
1994 if (detName == "MUON") {
1995 fRunMuonTracking = kTRUE;
2000 fTracker[iDet] = reconstructor->CreateTracker();
2001 if (!fTracker[iDet] && (iDet < 7)) {
2002 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2003 if (fStopOnError) return kFALSE;
2005 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2011 //_____________________________________________________________________________
2012 void AliReconstruction::CleanUp(TFile* file, TFile* fileOld)
2014 // delete trackers and the run loader and close and delete the file
2016 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2017 delete fReconstructor[iDet];
2018 fReconstructor[iDet] = NULL;
2019 fLoader[iDet] = NULL;
2020 delete fTracker[iDet];
2021 fTracker[iDet] = NULL;
2022 // delete fQADataMaker[iDet];
2023 // fQADataMaker[iDet] = NULL;
2028 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2029 delete fDiamondProfile;
2030 fDiamondProfile = NULL;
2040 if (fParentRawReader) delete fParentRawReader;
2041 fParentRawReader=NULL;
2051 gSystem->Unlink("AliESDs.old.root");
2055 //_____________________________________________________________________________
2057 Bool_t AliReconstruction::ReadESD(AliESDEvent*& esd, const char* recStep) const
2059 // read the ESD event from a file
2061 if (!esd) return kFALSE;
2063 sprintf(fileName, "ESD_%d.%d_%s.root",
2064 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2065 if (gSystem->AccessPathName(fileName)) return kFALSE;
2067 AliInfo(Form("reading ESD from file %s", fileName));
2068 AliDebug(1, Form("reading ESD from file %s", fileName));
2069 TFile* file = TFile::Open(fileName);
2070 if (!file || !file->IsOpen()) {
2071 AliError(Form("opening %s failed", fileName));
2078 esd = (AliESDEvent*) file->Get("ESD");
2087 //_____________________________________________________________________________
2088 void AliReconstruction::WriteESD(AliESDEvent* esd, const char* recStep) const
2090 // write the ESD event to a file
2094 sprintf(fileName, "ESD_%d.%d_%s.root",
2095 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2097 AliDebug(1, Form("writing ESD to file %s", fileName));
2098 TFile* file = TFile::Open(fileName, "recreate");
2099 if (!file || !file->IsOpen()) {
2100 AliError(Form("opening %s failed", fileName));
2112 //_____________________________________________________________________________
2113 void AliReconstruction::ESDFile2AODFile(TFile* esdFile, TFile* aodFile)
2115 // write all files from the given esd file to an aod file
2117 // create an AliAOD object
2118 AliAODEvent *aod = new AliAODEvent();
2119 aod->CreateStdContent();
2125 TTree *aodTree = new TTree("aodTree", "AliAOD tree");
2126 aodTree->Branch(aod->GetList());
2129 TTree *t = (TTree*) esdFile->Get("esdTree");
2130 AliESDEvent *esd = new AliESDEvent();
2131 esd->ReadFromTree(t);
2133 Int_t nEvents = t->GetEntries();
2135 // set arrays and pointers
2145 // loop over events and fill them
2146 for (Int_t iEvent = 0; iEvent < nEvents; ++iEvent) {
2147 //cout << "event: " << iEvent << endl;
2148 t->GetEntry(iEvent);
2150 // Multiplicity information needed by the header (to be revised!)
2151 Int_t nTracks = esd->GetNumberOfTracks();
2152 Int_t nPosTracks = 0;
2153 for (Int_t iTrack=0; iTrack<nTracks; ++iTrack)
2154 if (esd->GetTrack(iTrack)->Charge()> 0) nPosTracks++;
2156 // Access the header
2157 AliAODHeader *header = aod->GetHeader();
2160 header->SetRunNumber (esd->GetRunNumber() );
2161 header->SetBunchCrossNumber(esd->GetBunchCrossNumber());
2162 header->SetOrbitNumber (esd->GetOrbitNumber() );
2163 header->SetPeriodNumber (esd->GetPeriodNumber() );
2164 header->SetTriggerMask (esd->GetTriggerMask() );
2165 header->SetTriggerCluster (esd->GetTriggerCluster() );
2166 header->SetEventType (esd->GetEventType() );
2167 header->SetMagneticField (esd->GetMagneticField() );
2168 header->SetZDCN1Energy (esd->GetZDCN1Energy() );
2169 header->SetZDCP1Energy (esd->GetZDCP1Energy() );
2170 header->SetZDCN2Energy (esd->GetZDCN2Energy() );
2171 header->SetZDCP2Energy (esd->GetZDCP2Energy() );
2172 header->SetZDCEMEnergy (esd->GetZDCEMEnergy(0),esd->GetZDCEMEnergy(1));
2173 header->SetRefMultiplicity (nTracks);
2174 header->SetRefMultiplicityPos(nPosTracks);
2175 header->SetRefMultiplicityNeg(nTracks - nPosTracks);
2176 header->SetMuonMagFieldScale(-999.); // FIXME
2177 header->SetCentrality(-999.); // FIXME
2179 Int_t nV0s = esd->GetNumberOfV0s();
2180 Int_t nCascades = esd->GetNumberOfCascades();
2181 Int_t nKinks = esd->GetNumberOfKinks();
2182 Int_t nVertices = nV0s + 2*nCascades /*could lead to two vertices, one V0 and the Xi */+ nKinks + 1 /* = prim. vtx*/;
2184 Int_t nCaloClus = esd->GetNumberOfCaloClusters();
2186 Int_t nPmdClus = esd->GetNumberOfPmdTracks();
2188 aod->ResetStd(nTracks, nVertices, nV0s+nCascades, nJets, nCaloClus, nFmdClus, nPmdClus);
2190 // Array to take into account the tracks already added to the AOD
2191 Bool_t * usedTrack = NULL;
2193 usedTrack = new Bool_t[nTracks];
2194 for (Int_t iTrack=0; iTrack<nTracks; ++iTrack) usedTrack[iTrack]=kFALSE;
2196 // Array to take into account the V0s already added to the AOD
2197 Bool_t * usedV0 = NULL;
2199 usedV0 = new Bool_t[nV0s];
2200 for (Int_t iV0=0; iV0<nV0s; ++iV0) usedV0[iV0]=kFALSE;
2202 // Array to take into account the kinks already added to the AOD
2203 Bool_t * usedKink = NULL;
2205 usedKink = new Bool_t[nKinks];
2206 for (Int_t iKink=0; iKink<nKinks; ++iKink) usedKink[iKink]=kFALSE;
2209 // Access to the AOD container of vertices
2210 TClonesArray &vertices = *(aod->GetVertices());
2213 // Access to the AOD container of tracks
2214 TClonesArray &tracks = *(aod->GetTracks());
2217 // Access to the AOD container of V0s
2218 TClonesArray &V0s = *(aod->GetV0s());
2221 // Add primary vertex. The primary tracks will be defined
2222 // after the loops on the composite objects (V0, cascades, kinks)
2223 const AliESDVertex *vtx = esd->GetPrimaryVertex();
2225 vtx->GetXYZ(pos); // position
2226 vtx->GetCovMatrix(covVtx); //covariance matrix
2228 AliAODVertex * primary = new(vertices[jVertices++])
2229 AliAODVertex(pos, covVtx, vtx->GetChi2toNDF(), NULL, -1, AliAODVertex::kPrimary);
2232 AliAODTrack *aodTrack = 0x0;
2234 // Create vertices starting from the most complex objects
2237 for (Int_t nCascade = 0; nCascade < nCascades; ++nCascade) {
2238 AliESDcascade *cascade = esd->GetCascade(nCascade);
2240 cascade->GetXYZ(pos[0], pos[1], pos[2]);
2241 cascade->GetPosCovXi(covVtx);
2243 // Add the cascade vertex
2244 AliAODVertex * vcascade = new(vertices[jVertices++]) AliAODVertex(pos,
2246 cascade->GetChi2Xi(), // = chi2/NDF since NDF = 2*2-3
2249 AliAODVertex::kCascade);
2251 primary->AddDaughter(vcascade); // the cascade 'particle' (represented by a vertex) is added as a daughter to the primary vertex
2253 // Add the V0 from the cascade. The ESD class have to be optimized...
2254 // Now we have to search for the corresponding V0 in the list of V0s
2255 // using the indeces of the positive and negative tracks
2257 Int_t posFromV0 = cascade->GetPindex();
2258 Int_t negFromV0 = cascade->GetNindex();
2261 AliESDv0 * v0 = 0x0;
2264 for (Int_t iV0=0; iV0<nV0s; ++iV0) {
2266 v0 = esd->GetV0(iV0);
2267 Int_t posV0 = v0->GetPindex();
2268 Int_t negV0 = v0->GetNindex();
2270 if (posV0==posFromV0 && negV0==negFromV0) {
2276 AliAODVertex * vV0FromCascade = 0x0;
2278 if (indV0>-1 && !usedV0[indV0]) {
2280 // the V0 exists in the array of V0s and is not used
2282 usedV0[indV0] = kTRUE;
2284 v0->GetXYZ(pos[0], pos[1], pos[2]);
2285 v0->GetPosCov(covVtx);
2287 vV0FromCascade = new(vertices[jVertices++]) AliAODVertex(pos,
2289 v0->GetChi2V0(), // = chi2/NDF since NDF = 2*2-3
2295 // the V0 doesn't exist in the array of V0s or was used
2296 cerr << "Error: event " << iEvent << " cascade " << nCascade
2297 << " The V0 " << indV0
2298 << " doesn't exist in the array of V0s or was used!" << endl;
2300 cascade->GetXYZ(pos[0], pos[1], pos[2]);
2301 cascade->GetPosCov(covVtx);
2303 vV0FromCascade = new(vertices[jVertices++]) AliAODVertex(pos,
2305 v0->GetChi2V0(), // = chi2/NDF since NDF = 2*2-3
2309 vcascade->AddDaughter(vV0FromCascade);
2313 // Add the positive tracks from the V0
2315 if (! usedTrack[posFromV0]) {
2317 usedTrack[posFromV0] = kTRUE;
2319 AliESDtrack *esdTrack = esd->GetTrack(posFromV0);
2320 esdTrack->GetPxPyPz(p_pos);
2321 esdTrack->GetXYZ(pos);
2322 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2323 esdTrack->GetESDpid(pid);
2325 vV0FromCascade->AddDaughter(aodTrack =
2326 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2327 esdTrack->GetLabel(),
2333 (Short_t)esdTrack->Charge(),
2334 esdTrack->GetITSClusterMap(),
2337 kTRUE, // check if this is right
2338 kFALSE, // check if this is right
2339 AliAODTrack::kSecondary)
2341 aodTrack->ConvertAliPIDtoAODPID();
2344 cerr << "Error: event " << iEvent << " cascade " << nCascade
2345 << " track " << posFromV0 << " has already been used!" << endl;
2348 // Add the negative tracks from the V0
2350 if (!usedTrack[negFromV0]) {
2352 usedTrack[negFromV0] = kTRUE;
2354 AliESDtrack *esdTrack = esd->GetTrack(negFromV0);
2355 esdTrack->GetPxPyPz(p_neg);
2356 esdTrack->GetXYZ(pos);
2357 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2358 esdTrack->GetESDpid(pid);
2360 vV0FromCascade->AddDaughter(aodTrack =
2361 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2362 esdTrack->GetLabel(),
2368 (Short_t)esdTrack->Charge(),
2369 esdTrack->GetITSClusterMap(),
2372 kTRUE, // check if this is right
2373 kFALSE, // check if this is right
2374 AliAODTrack::kSecondary)
2376 aodTrack->ConvertAliPIDtoAODPID();
2379 cerr << "Error: event " << iEvent << " cascade " << nCascade
2380 << " track " << negFromV0 << " has already been used!" << endl;
2383 // add it to the V0 array as well
2384 Double_t d0[2] = { -999., -99.};
2385 // counting is probably wrong
2386 new(V0s[jV0s++]) AliAODv0(vV0FromCascade, -999., -99., p_pos, p_neg, d0); // to be refined
2388 // Add the bachelor track from the cascade
2390 Int_t bachelor = cascade->GetBindex();
2392 if(!usedTrack[bachelor]) {
2394 usedTrack[bachelor] = kTRUE;
2396 AliESDtrack *esdTrack = esd->GetTrack(bachelor);
2397 esdTrack->GetPxPyPz(p);
2398 esdTrack->GetXYZ(pos);
2399 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2400 esdTrack->GetESDpid(pid);
2402 vcascade->AddDaughter(aodTrack =
2403 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2404 esdTrack->GetLabel(),
2410 (Short_t)esdTrack->Charge(),
2411 esdTrack->GetITSClusterMap(),
2414 kTRUE, // check if this is right
2415 kFALSE, // check if this is right
2416 AliAODTrack::kSecondary)
2418 aodTrack->ConvertAliPIDtoAODPID();
2421 cerr << "Error: event " << iEvent << " cascade " << nCascade
2422 << " track " << bachelor << " has already been used!" << endl;
2425 // Add the primary track of the cascade (if any)
2427 } // end of the loop on cascades
2431 for (Int_t nV0 = 0; nV0 < nV0s; ++nV0) {
2433 if (usedV0[nV0]) continue; // skip if aready added to the AOD
2435 AliESDv0 *v0 = esd->GetV0(nV0);
2437 v0->GetXYZ(pos[0], pos[1], pos[2]);
2438 v0->GetPosCov(covVtx);
2440 AliAODVertex * vV0 =
2441 new(vertices[jVertices++]) AliAODVertex(pos,
2443 v0->GetChi2V0(), // = chi2/NDF since NDF = 2*2-3
2447 primary->AddDaughter(vV0);
2449 Int_t posFromV0 = v0->GetPindex();
2450 Int_t negFromV0 = v0->GetNindex();
2452 // Add the positive tracks from the V0
2454 if (!usedTrack[posFromV0]) {
2456 usedTrack[posFromV0] = kTRUE;
2458 AliESDtrack *esdTrack = esd->GetTrack(posFromV0);
2459 esdTrack->GetPxPyPz(p_pos);
2460 esdTrack->GetXYZ(pos);
2461 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2462 esdTrack->GetESDpid(pid);
2464 vV0->AddDaughter(aodTrack =
2465 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2466 esdTrack->GetLabel(),
2472 (Short_t)esdTrack->Charge(),
2473 esdTrack->GetITSClusterMap(),
2476 kTRUE, // check if this is right
2477 kFALSE, // check if this is right
2478 AliAODTrack::kSecondary)
2480 aodTrack->ConvertAliPIDtoAODPID();
2483 cerr << "Error: event " << iEvent << " V0 " << nV0
2484 << " track " << posFromV0 << " has already been used!" << endl;
2487 // Add the negative tracks from the V0
2489 if (!usedTrack[negFromV0]) {
2491 usedTrack[negFromV0] = kTRUE;
2493 AliESDtrack *esdTrack = esd->GetTrack(negFromV0);
2494 esdTrack->GetPxPyPz(p_neg);
2495 esdTrack->GetXYZ(pos);
2496 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2497 esdTrack->GetESDpid(pid);
2499 vV0->AddDaughter(aodTrack =
2500 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2501 esdTrack->GetLabel(),
2507 (Short_t)esdTrack->Charge(),
2508 esdTrack->GetITSClusterMap(),
2511 kTRUE, // check if this is right
2512 kFALSE, // check if this is right
2513 AliAODTrack::kSecondary)
2515 aodTrack->ConvertAliPIDtoAODPID();
2518 cerr << "Error: event " << iEvent << " V0 " << nV0
2519 << " track " << negFromV0 << " has already been used!" << endl;
2522 // add it to the V0 array as well
2523 Double_t d0[2] = { 999., 99.};
2524 new(V0s[jV0s++]) AliAODv0(vV0, 999., 99., p_pos, p_neg, d0); // to be refined
2527 // end of the loop on V0s
2529 // Kinks: it is a big mess the access to the information in the kinks
2530 // The loop is on the tracks in order to find the mother and daugther of each kink
2533 for (Int_t iTrack=0; iTrack<nTracks; ++iTrack) {
2535 AliESDtrack * esdTrack = esd->GetTrack(iTrack);
2537 Int_t ikink = esdTrack->GetKinkIndex(0);
2540 // Negative kink index: mother, positive: daughter
2542 // Search for the second track of the kink
2544 for (Int_t jTrack = iTrack+1; jTrack<nTracks; ++jTrack) {
2546 AliESDtrack * esdTrack1 = esd->GetTrack(jTrack);
2548 Int_t jkink = esdTrack1->GetKinkIndex(0);
2550 if ( TMath::Abs(ikink)==TMath::Abs(jkink) ) {
2552 // The two tracks are from the same kink
2554 if (usedKink[TMath::Abs(ikink)-1]) continue; // skip used kinks
2557 Int_t idaughter = -1;
2559 if (ikink<0 && jkink>0) {
2564 else if (ikink>0 && jkink<0) {
2570 cerr << "Error: Wrong combination of kink indexes: "
2571 << ikink << " " << jkink << endl;
2575 // Add the mother track
2577 AliAODTrack * mother = NULL;
2579 if (!usedTrack[imother]) {
2581 usedTrack[imother] = kTRUE;
2583 AliESDtrack *esdTrack = esd->GetTrack(imother);
2584 esdTrack->GetPxPyPz(p);
2585 esdTrack->GetXYZ(pos);
2586 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2587 esdTrack->GetESDpid(pid);
2590 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2591 esdTrack->GetLabel(),
2597 (Short_t)esdTrack->Charge(),
2598 esdTrack->GetITSClusterMap(),
2601 kTRUE, // check if this is right
2602 kTRUE, // check if this is right
2603 AliAODTrack::kPrimary);
2604 primary->AddDaughter(mother);
2605 mother->ConvertAliPIDtoAODPID();
2608 cerr << "Error: event " << iEvent << " kink " << TMath::Abs(ikink)-1
2609 << " track " << imother << " has already been used!" << endl;
2612 // Add the kink vertex
2613 AliESDkink * kink = esd->GetKink(TMath::Abs(ikink)-1);
2615 AliAODVertex * vkink =
2616 new(vertices[jVertices++]) AliAODVertex(kink->GetPosition(),
2620 esdTrack->GetID(), // This is the track ID of the mother's track!
2621 AliAODVertex::kKink);
2622 // Add the daughter track
2624 AliAODTrack * daughter = NULL;
2626 if (!usedTrack[idaughter]) {
2628 usedTrack[idaughter] = kTRUE;
2630 AliESDtrack *esdTrack = esd->GetTrack(idaughter);
2631 esdTrack->GetPxPyPz(p);
2632 esdTrack->GetXYZ(pos);
2633 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2634 esdTrack->GetESDpid(pid);
2637 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2638 esdTrack->GetLabel(),
2644 (Short_t)esdTrack->Charge(),
2645 esdTrack->GetITSClusterMap(),
2648 kTRUE, // check if this is right
2649 kTRUE, // check if this is right
2650 AliAODTrack::kPrimary);
2651 vkink->AddDaughter(daughter);
2652 daughter->ConvertAliPIDtoAODPID();
2655 cerr << "Error: event " << iEvent << " kink " << TMath::Abs(ikink)-1
2656 << " track " << idaughter << " has already been used!" << endl;
2662 vertices.Expand(jVertices);
2664 // Tracks (primary and orphan)
2665 for (Int_t nTrack = 0; nTrack < nTracks; ++nTrack) {
2667 if (usedTrack[nTrack]) continue;
2669 AliESDtrack *esdTrack = esd->GetTrack(nTrack);
2670 esdTrack->GetPxPyPz(p);
2671 esdTrack->GetXYZ(pos);
2672 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2673 esdTrack->GetESDpid(pid);
2675 Float_t impactXY, impactZ;
2677 esdTrack->GetImpactParameters(impactXY,impactZ);
2680 // track inside the beam pipe
2682 primary->AddDaughter(aodTrack =
2683 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2684 esdTrack->GetLabel(),
2690 (Short_t)esdTrack->Charge(),
2691 esdTrack->GetITSClusterMap(),
2694 kTRUE, // check if this is right
2695 kTRUE, // check if this is right
2696 AliAODTrack::kPrimary)
2698 aodTrack->ConvertAliPIDtoAODPID();
2701 // outside the beam pipe: orphan track
2702 // Don't write them anymore!
2705 } // end of loop on tracks
2708 Int_t nMuTracks = esd->GetNumberOfMuonTracks();
2709 for (Int_t nMuTrack = 0; nMuTrack < nMuTracks; ++nMuTrack) {
2711 AliESDMuonTrack *esdMuTrack = esd->GetMuonTrack(nMuTrack);
2712 p[0] = esdMuTrack->Px();
2713 p[1] = esdMuTrack->Py();
2714 p[2] = esdMuTrack->Pz();
2715 pos[0] = primary->GetX();
2716 pos[1] = primary->GetY();
2717 pos[2] = primary->GetZ();
2719 // has to be changed once the muon pid is provided by the ESD
2720 for (Int_t i = 0; i < 10; pid[i++] = 0.); pid[AliAODTrack::kMuon]=1.;
2722 primary->AddDaughter(aodTrack =
2723 new(tracks[jTracks++]) AliAODTrack(0, // no ID provided
2724 0, // no label provided
2729 NULL, // no covariance matrix provided
2730 esdMuTrack->Charge(),
2731 0, // ITSClusterMap is set below
2734 kFALSE, // muon tracks are not used to fit the primary vtx
2735 kFALSE, // not used for vertex fit
2736 AliAODTrack::kPrimary)
2739 aodTrack->SetHitsPatternInTrigCh(esdMuTrack->GetHitsPatternInTrigCh());
2740 Int_t track2Trigger = esdMuTrack->GetMatchTrigger();
2741 aodTrack->SetMatchTrigger(track2Trigger);
2743 aodTrack->SetChi2MatchTrigger(esdMuTrack->GetChi2MatchTrigger());
2745 aodTrack->SetChi2MatchTrigger(0.);
2747 tracks.Expand(jTracks); // remove 'empty slots' due to unwritten tracks
2749 // Access to the AOD container of PMD clusters
2750 TClonesArray &pmdClusters = *(aod->GetPmdClusters());
2751 Int_t jPmdClusters=0;
2753 for (Int_t iPmd = 0; iPmd < nPmdClus; ++iPmd) {
2754 // file pmd clusters, to be revised!
2755 AliESDPmdTrack *pmdTrack = esd->GetPmdTrack(iPmd);
2758 Double_t pos[3] = { pmdTrack->GetClusterX(), pmdTrack->GetClusterY(), pmdTrack->GetClusterZ() };
2759 Double_t pid[9] = { 0., 0., 0., 0., 0., 0., 0., 0., 0. }; // to be revised!
2761 // assoc cluster not set
2762 new(pmdClusters[jPmdClusters++]) AliAODPmdCluster(iPmd, nLabel, label, pmdTrack->GetClusterADC(), pos, pid);
2765 // Access to the AOD container of clusters
2766 TClonesArray &caloClusters = *(aod->GetCaloClusters());
2769 for (Int_t iClust=0; iClust<nCaloClus; ++iClust) {
2771 AliESDCaloCluster * cluster = esd->GetCaloCluster(iClust);
2773 Int_t id = cluster->GetID();
2776 Float_t energy = cluster->E();
2777 cluster->GetPosition(posF);
2778 Char_t ttype=AliAODCluster::kUndef;
2780 if (cluster->GetClusterType() == AliESDCaloCluster::kPHOSCluster) {
2781 ttype=AliAODCluster::kPHOSNeutral;
2783 else if (cluster->GetClusterType() == AliESDCaloCluster::kEMCALClusterv1) {
2784 ttype = AliAODCluster::kEMCALClusterv1;
2788 AliAODCaloCluster *caloCluster = new(caloClusters[jClusters++]) AliAODCaloCluster(id,
2796 caloCluster->SetCaloCluster(); // to be refined!
2799 caloClusters.Expand(jClusters); // resize TObjArray to 'remove' slots for pseudo clusters
2800 // end of loop on calo clusters
2802 // fill EMCAL cell info
2803 if (esd->GetEMCALCells()) { // protection against missing ESD information
2804 AliESDCaloCells &esdEMcells = *(esd->GetEMCALCells());
2805 Int_t nEMcell = esdEMcells.GetNumberOfCells() ;
2807 AliAODCaloCells &aodEMcells = *(aod->GetEMCALCells());
2808 aodEMcells.CreateContainer(nEMcell);
2809 aodEMcells.SetType(AliAODCaloCells::kEMCAL);
2810 for (Int_t iCell = 0; iCell < nEMcell; iCell++) {
2811 aodEMcells.SetCell(iCell,esdEMcells.GetCellNumber(iCell),esdEMcells.GetAmplitude(iCell));
2816 // fill PHOS cell info
2817 if (esd->GetPHOSCells()) { // protection against missing ESD information
2818 AliESDCaloCells &esdPHcells = *(esd->GetPHOSCells());
2819 Int_t nPHcell = esdPHcells.GetNumberOfCells() ;
2821 AliAODCaloCells &aodPHcells = *(aod->GetPHOSCells());
2822 aodPHcells.CreateContainer(nPHcell);
2823 aodPHcells.SetType(AliAODCaloCells::kPHOS);
2824 for (Int_t iCell = 0; iCell < nPHcell; iCell++) {
2825 aodPHcells.SetCell(iCell,esdPHcells.GetCellNumber(iCell),esdPHcells.GetAmplitude(iCell));
2831 AliAODTracklets &SPDTracklets = *(aod->GetTracklets());
2832 const AliMultiplicity *mult = esd->GetMultiplicity();
2834 if (mult->GetNumberOfTracklets()>0) {
2835 SPDTracklets.CreateContainer(mult->GetNumberOfTracklets());
2837 for (Int_t n=0; n<mult->GetNumberOfTracklets(); n++) {
2838 SPDTracklets.SetTracklet(n, mult->GetTheta(n), mult->GetPhi(n), mult->GetDeltaPhi(n), mult->GetLabel(n));
2842 Printf("ERROR: AliMultiplicity could not be retrieved from ESD");
2845 delete [] usedTrack;
2849 // fill the tree for this event
2851 } // end of event loop
2853 aodTree->GetUserInfo()->Add(aod);
2855 // write the tree to the specified file
2856 aodFile = aodTree->GetCurrentFile();
2863 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2865 // Write space-points which are then used in the alignment procedures
2866 // For the moment only ITS, TRD and TPC
2868 // Load TOF clusters
2870 fLoader[3]->LoadRecPoints("read");
2871 TTree* tree = fLoader[3]->TreeR();
2873 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[3]));
2876 fTracker[3]->LoadClusters(tree);
2878 Int_t ntracks = esd->GetNumberOfTracks();
2879 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2881 AliESDtrack *track = esd->GetTrack(itrack);
2884 for (Int_t iDet = 3; iDet >= 0; iDet--)
2885 nsp += track->GetNcls(iDet);
2887 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2888 track->SetTrackPointArray(sp);
2890 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2891 AliTracker *tracker = fTracker[iDet];
2892 if (!tracker) continue;
2893 Int_t nspdet = track->GetNcls(iDet);
2894 if (nspdet <= 0) continue;
2895 track->GetClusters(iDet,idx);
2899 while (isp2 < nspdet) {
2901 TString dets = fgkDetectorName[iDet];
2902 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2903 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2904 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2905 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2906 isvalid = tracker->GetTrackPointTrackingError(idx[isp2],p,track);
2908 isvalid = tracker->GetTrackPoint(idx[isp2],p);
2911 const Int_t kNTPCmax = 159;
2912 if (iDet==1 && isp2>kNTPCmax) break; // to be fixed
2913 if (!isvalid) continue;
2914 sp->AddPoint(isptrack,&p); isptrack++; isp++;
2920 fTracker[3]->UnloadClusters();
2921 fLoader[3]->UnloadRecPoints();
2925 //_____________________________________________________________________________
2926 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2928 // The method reads the raw-data error log
2929 // accumulated within the rawReader.
2930 // It extracts the raw-data errors related to
2931 // the current event and stores them into
2932 // a TClonesArray inside the esd object.
2934 if (!fRawReader) return;
2936 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2938 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2940 if (iEvent != log->GetEventNumber()) continue;
2942 esd->AddRawDataErrorLog(log);
2947 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString fName){
2948 // Dump a file content into a char in TNamed
2950 in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2951 Int_t kBytes = (Int_t)in.tellg();
2952 printf("Size: %d \n",kBytes);
2955 char* memblock = new char [kBytes];
2956 in.seekg (0, ios::beg);
2957 in.read (memblock, kBytes);
2959 TString fData(memblock,kBytes);
2960 fn = new TNamed(fName,fData);
2961 printf("fData Size: %d \n",fData.Sizeof());
2962 printf("fName Size: %d \n",fName.Sizeof());
2963 printf("fn Size: %d \n",fn->Sizeof());
2967 AliInfo(Form("Could not Open %s\n",fPath.Data()));
2973 void AliReconstruction::TNamedToFile(TTree* fTree, TString fName){
2974 // This is not really needed in AliReconstruction at the moment
2975 // but can serve as a template
2977 TList *fList = fTree->GetUserInfo();
2978 TNamed *fn = (TNamed*)fList->FindObject(fName.Data());
2979 printf("fn Size: %d \n",fn->Sizeof());
2981 TString fTmp(fn->GetName()); // to be 100% sure in principle fName also works
2982 const char* cdata = fn->GetTitle();
2983 printf("fTmp Size %d\n",fTmp.Sizeof());
2985 int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2986 printf("calculated size %d\n",size);
2987 ofstream out(fName.Data(),ios::out | ios::binary);
2988 out.write(cdata,size);
2993 //_____________________________________________________________________________
2994 AliQADataMakerRec * AliReconstruction::GetQADataMaker(Int_t iDet)
2996 // get the quality assurance data maker object and the loader for a detector
2998 if (fQADataMaker[iDet])
2999 return fQADataMaker[iDet];
3001 AliQADataMakerRec * qadm = NULL;
3002 if (iDet == fgkNDetectors) { //Global QA
3003 qadm = new AliGlobalQADataMaker();
3004 fQADataMaker[iDet] = qadm;
3008 // load the QA data maker object
3009 TPluginManager* pluginManager = gROOT->GetPluginManager();
3010 TString detName = fgkDetectorName[iDet];
3011 TString qadmName = "Ali" + detName + "QADataMakerRec";
3012 if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT"))
3015 // first check if a plugin is defined for the quality assurance data maker
3016 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
3017 // if not, add a plugin for it
3018 if (!pluginHandler) {
3019 AliDebug(1, Form("defining plugin for %s", qadmName.Data()));
3020 TString libs = gSystem->GetLibraries();
3021 if (libs.Contains("lib" + detName + "base.so") ||
3022 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
3023 pluginManager->AddHandler("AliQADataMakerRec", detName,
3024 qadmName, detName + "qadm", qadmName + "()");
3026 pluginManager->AddHandler("AliQADataMakerRec", detName,
3027 qadmName, detName, qadmName + "()");
3029 pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
3031 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
3032 qadm = (AliQADataMakerRec *) pluginHandler->ExecPlugin(0);
3035 fQADataMaker[iDet] = qadm;
3040 //_____________________________________________________________________________
3041 Bool_t AliReconstruction::RunQA(const char* detectors, AliESDEvent *& esd)
3043 // run the Quality Assurance data producer
3045 AliCodeTimerAuto("")
3046 TString detStr = detectors;
3047 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3048 if (!IsSelected(fgkDetectorName[iDet], detStr))
3050 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
3053 AliCodeTimerStart(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
3054 AliInfo(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
3056 qadm->Exec(AliQA::kESDS, esd) ;
3059 AliCodeTimerStop(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
3061 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
3062 AliError(Form("the following detectors were not found: %s",
3072 //_____________________________________________________________________________
3073 void AliReconstruction::CheckQA()
3075 // check the QA of SIM for this run and remove the detectors
3076 // with status Fatal
3078 TString newRunLocalReconstruction ;
3079 TString newRunTracking ;
3080 TString newFillESD ;
3082 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
3083 TString detName(AliQA::GetDetName(iDet)) ;
3084 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX(iDet)) ;
3085 if ( qa->IsSet(AliQA::DETECTORINDEX(iDet), AliQA::kSIM, AliQA::kFATAL)) {
3086 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
3088 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
3089 fRunLocalReconstruction.Contains("ALL") ) {
3090 newRunLocalReconstruction += detName ;
3091 newRunLocalReconstruction += " " ;
3093 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
3094 fRunTracking.Contains("ALL") ) {
3095 newRunTracking += detName ;
3096 newRunTracking += " " ;
3098 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
3099 fFillESD.Contains("ALL") ) {
3100 newFillESD += detName ;
3105 fRunLocalReconstruction = newRunLocalReconstruction ;
3106 fRunTracking = newRunTracking ;
3107 fFillESD = newFillESD ;
3110 //_____________________________________________________________________________
3111 Int_t AliReconstruction::GetDetIndex(const char* detector)
3113 // return the detector index corresponding to detector
3115 for (index = 0; index < fgkNDetectors ; index++) {
3116 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
3121 //_____________________________________________________________________________
3122 Bool_t AliReconstruction::FinishPlaneEff() {
3124 // Here execute all the necessary operationis, at the end of the tracking phase,
3125 // in case that evaluation of PlaneEfficiencies was required for some detector.
3126 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
3128 // This Preliminary version works only FOR ITS !!!!!
3129 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3132 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3135 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3136 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
3137 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3138 if(fTracker[iDet]) ret=fTracker[iDet]->GetPlaneEff()->WriteIntoCDB();
3142 //_____________________________________________________________________________
3143 Bool_t AliReconstruction::InitPlaneEff() {
3145 // Here execute all the necessary operations, before of the tracking phase,
3146 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
3147 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
3148 // which should be updated/recalculated.
3150 // This Preliminary version will work only FOR ITS !!!!!
3151 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3154 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3156 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));