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);
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);
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->SetVertex(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);
987 if (fWriteESDfriend) {
988 esdf->~AliESDfriend();
989 new (esdf) AliESDfriend(); // Reset...
990 esd->GetESDfriend(esdf);
997 if (fCheckPointLevel > 0) WriteESD(esd, "final");
1000 if (fWriteESDfriend) {
1001 esdf->~AliESDfriend();
1002 new (esdf) AliESDfriend(); // Reset...
1005 gSystem->GetProcInfo(&ProcInfo);
1006 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1009 // End of cycle for the in-loop QA
1010 if (fRunQA && fInLoopQA) {
1011 RunQA(fFillESD.Data(), esd);
1012 TString detStr(fFillESD);
1013 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1014 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1015 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1016 if (!qadm) continue;
1017 qadm->EndOfCycle(AliQA::kRECPOINTS);
1018 qadm->EndOfCycle(AliQA::kESDS);
1022 AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
1024 qadm->EndOfCycle(AliQA::kRECPOINTS);
1025 qadm->EndOfCycle(AliQA::kESDS);
1031 //******** End of the loop over events
1035 tree->GetUserInfo()->Add(esd);
1036 hlttree->GetUserInfo()->Add(hltesd);
1038 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1039 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1041 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1042 cdbMapCopy->SetOwner(1);
1043 cdbMapCopy->SetName("cdbMap");
1044 TIter iter(cdbMap->GetTable());
1047 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1048 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1049 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1050 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1053 TList *cdbListCopy = new TList();
1054 cdbListCopy->SetOwner(1);
1055 cdbListCopy->SetName("cdbList");
1057 TIter iter2(cdbList);
1060 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1061 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1064 tree->GetUserInfo()->Add(cdbMapCopy);
1065 tree->GetUserInfo()->Add(cdbListCopy);
1068 if(fESDPar.Contains("ESD.par")){
1069 AliInfo("Attaching ESD.par to Tree");
1070 TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1071 tree->GetUserInfo()->Add(fn);
1077 if (fWriteESDfriend)
1078 tree->SetBranchStatus("ESDfriend*",0);
1079 // we want to have only one tree version number
1080 tree->Write(tree->GetName(),TObject::kOverwrite);
1083 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1084 if (fRunPlaneEff && !FinishPlaneEff()) {
1085 AliWarning("Finish PlaneEff evaluation failed");
1089 CleanUp(file, fileOld);
1092 TFile *esdFile = TFile::Open("AliESDs.root", "READONLY");
1093 TFile *aodFile = TFile::Open("AliAOD.root", "RECREATE");
1094 ESDFile2AODFile(esdFile, aodFile);
1099 // Create tags for the events in the ESD tree (the ESD tree is always present)
1100 // In case of empty events the tags will contain dummy values
1101 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1102 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPList);
1104 AliAODTagCreator *aodtagCreator = new AliAODTagCreator();
1105 aodtagCreator->CreateAODTags(fFirstEvent,fLastEvent,fGRPList);
1108 //Finish QA and end of cycle for out-of-loop QA
1109 if (fRunQA && !fInLoopQA) {
1110 qas.Run(fRunLocalReconstruction.Data(), AliQA::kRECPOINTS);
1112 qas.Run(fRunTracking.Data(), AliQA::kESDS);
1115 AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
1117 qadm->EndOfCycle(AliQA::kRECPOINTS);
1118 qadm->EndOfCycle(AliQA::kESDS);
1124 // Cleanup of CDB manager: cache and active storages!
1125 AliCDBManager::Instance()->ClearCache();
1132 //_____________________________________________________________________________
1133 Bool_t AliReconstruction::RunLocalReconstruction(const TString& /*detectors*/)
1135 // run the local reconstruction
1136 static Int_t eventNr=0;
1137 AliCodeTimerAuto("")
1139 // AliCDBManager* man = AliCDBManager::Instance();
1140 // Bool_t origCache = man->GetCacheFlag();
1142 // TString detStr = detectors;
1143 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1144 // if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1145 // AliReconstructor* reconstructor = GetReconstructor(iDet);
1146 // if (!reconstructor) continue;
1147 // if (reconstructor->HasLocalReconstruction()) continue;
1149 // AliCodeTimerStart(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1150 // AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1152 // AliCodeTimerStart(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1153 // AliInfo(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1155 // man->SetCacheFlag(kTRUE);
1156 // TString calibPath = Form("%s/Calib/*", fgkDetectorName[iDet]);
1157 // man->GetAll(calibPath); // entries are cached!
1159 // AliCodeTimerStop(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1161 // if (fRawReader) {
1162 // fRawReader->RewindEvents();
1163 // reconstructor->Reconstruct(fRunLoader, fRawReader);
1165 // reconstructor->Reconstruct(fRunLoader);
1168 // AliCodeTimerStop(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1169 // AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr));
1171 // // unload calibration data
1172 // man->UnloadFromCache(calibPath);
1173 // //man->ClearCache();
1176 // man->SetCacheFlag(origCache);
1178 // if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1179 // AliError(Form("the following detectors were not found: %s",
1181 // if (fStopOnError) return kFALSE;
1188 //_____________________________________________________________________________
1189 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1191 // run the local reconstruction
1193 static Int_t eventNr=0;
1194 AliCodeTimerAuto("")
1196 TString detStr = detectors;
1197 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1198 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1199 AliReconstructor* reconstructor = GetReconstructor(iDet);
1200 if (!reconstructor) continue;
1201 AliLoader* loader = fLoader[iDet];
1203 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1206 // conversion of digits
1207 if (fRawReader && reconstructor->HasDigitConversion()) {
1208 AliInfo(Form("converting raw data digits into root objects for %s",
1209 fgkDetectorName[iDet]));
1210 AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1211 fgkDetectorName[iDet]));
1212 loader->LoadDigits("update");
1213 loader->CleanDigits();
1214 loader->MakeDigitsContainer();
1215 TTree* digitsTree = loader->TreeD();
1216 reconstructor->ConvertDigits(fRawReader, digitsTree);
1217 loader->WriteDigits("OVERWRITE");
1218 loader->UnloadDigits();
1220 // local reconstruction
1221 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1222 AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1223 loader->LoadRecPoints("update");
1224 loader->CleanRecPoints();
1225 loader->MakeRecPointsContainer();
1226 TTree* clustersTree = loader->TreeR();
1227 if (fRawReader && !reconstructor->HasDigitConversion()) {
1228 reconstructor->Reconstruct(fRawReader, clustersTree);
1230 loader->LoadDigits("read");
1231 TTree* digitsTree = loader->TreeD();
1233 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1234 if (fStopOnError) return kFALSE;
1236 reconstructor->Reconstruct(digitsTree, clustersTree);
1238 loader->UnloadDigits();
1241 // In-loop QA for local reconstrucion
1242 if (fRunQA && fInLoopQA) {
1243 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1246 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1248 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1250 qadm->Exec(AliQA::kRECPOINTS, clustersTree) ;
1253 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1257 loader->WriteRecPoints("OVERWRITE");
1258 loader->UnloadRecPoints();
1259 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1262 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1263 AliError(Form("the following detectors were not found: %s",
1265 if (fStopOnError) return kFALSE;
1271 //_____________________________________________________________________________
1272 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1274 // run the barrel tracking
1276 AliCodeTimerAuto("")
1278 AliESDVertex* vertex = NULL;
1279 Double_t vtxPos[3] = {0, 0, 0};
1280 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1281 TArrayF mcVertex(3);
1282 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1283 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1284 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1288 if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1289 AliInfo("running the ITS vertex finder");
1290 if (fLoader[0]) fLoader[0]->LoadRecPoints();
1291 vertex = fVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber());
1292 if (fLoader[0]) fLoader[0]->UnloadRecPoints();
1294 AliWarning("Vertex not found");
1295 vertex = new AliESDVertex();
1296 vertex->SetName("default");
1299 vertex->SetName("reconstructed");
1303 AliInfo("getting the primary vertex from MC");
1304 vertex = new AliESDVertex(vtxPos, vtxErr);
1308 vertex->GetXYZ(vtxPos);
1309 vertex->GetSigmaXYZ(vtxErr);
1311 AliWarning("no vertex reconstructed");
1312 vertex = new AliESDVertex(vtxPos, vtxErr);
1314 esd->SetVertex(vertex);
1315 // if SPD multiplicity has been determined, it is stored in the ESD
1316 AliMultiplicity *mult = fVertexer->GetMultiplicity();
1317 if(mult)esd->SetMultiplicity(mult);
1319 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1320 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1327 //_____________________________________________________________________________
1328 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
1330 // run the HLT barrel tracking
1332 AliCodeTimerAuto("")
1335 AliError("Missing runLoader!");
1339 AliInfo("running HLT tracking");
1341 // Get a pointer to the HLT reconstructor
1342 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
1343 if (!reconstructor) return kFALSE;
1346 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1347 TString detName = fgkDetectorName[iDet];
1348 AliDebug(1, Form("%s HLT tracking", detName.Data()));
1349 reconstructor->SetOption(detName.Data());
1350 AliTracker *tracker = reconstructor->CreateTracker();
1352 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
1353 if (fStopOnError) return kFALSE;
1357 Double_t vtxErr[3]={0.005,0.005,0.010};
1358 const AliESDVertex *vertex = esd->GetVertex();
1359 vertex->GetXYZ(vtxPos);
1360 tracker->SetVertex(vtxPos,vtxErr);
1362 fLoader[iDet]->LoadRecPoints("read");
1363 TTree* tree = fLoader[iDet]->TreeR();
1365 AliError(Form("Can't get the %s cluster tree", detName.Data()));
1368 tracker->LoadClusters(tree);
1370 if (tracker->Clusters2Tracks(esd) != 0) {
1371 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
1375 tracker->UnloadClusters();
1383 //_____________________________________________________________________________
1384 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
1386 // run the muon spectrometer tracking
1388 AliCodeTimerAuto("")
1391 AliError("Missing runLoader!");
1394 Int_t iDet = 7; // for MUON
1396 AliInfo("is running...");
1398 // Get a pointer to the MUON reconstructor
1399 AliReconstructor *reconstructor = GetReconstructor(iDet);
1400 if (!reconstructor) return kFALSE;
1403 TString detName = fgkDetectorName[iDet];
1404 AliDebug(1, Form("%s tracking", detName.Data()));
1405 AliTracker *tracker = reconstructor->CreateTracker();
1407 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1412 fLoader[iDet]->LoadRecPoints("read");
1414 tracker->LoadClusters(fLoader[iDet]->TreeR());
1416 Int_t rv = tracker->Clusters2Tracks(esd);
1420 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1424 fLoader[iDet]->UnloadRecPoints();
1426 tracker->UnloadClusters();
1434 //_____________________________________________________________________________
1435 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
1437 // run the barrel tracking
1438 static Int_t eventNr=0;
1439 AliCodeTimerAuto("")
1441 AliInfo("running tracking");
1443 //Fill the ESD with the T0 info (will be used by the TOF)
1444 if (fReconstructor[11] && fLoader[11]) {
1445 fLoader[11]->LoadRecPoints("READ");
1446 TTree *treeR = fLoader[11]->TreeR();
1447 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
1450 // pass 1: TPC + ITS inwards
1451 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1452 if (!fTracker[iDet]) continue;
1453 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
1456 fLoader[iDet]->LoadRecPoints("read");
1457 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
1458 TTree* tree = fLoader[iDet]->TreeR();
1460 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1463 fTracker[iDet]->LoadClusters(tree);
1464 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1466 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
1467 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1470 if (fCheckPointLevel > 1) {
1471 WriteESD(esd, Form("%s.tracking", fgkDetectorName[iDet]));
1473 // preliminary PID in TPC needed by the ITS tracker
1475 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1476 AliESDpid::MakePID(esd);
1478 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
1481 // pass 2: ALL backwards
1483 if (fRunQA && fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1485 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1486 if (!fTracker[iDet]) continue;
1487 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
1490 if (iDet > 1) { // all except ITS, TPC
1492 fLoader[iDet]->LoadRecPoints("read");
1493 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
1494 tree = fLoader[iDet]->TreeR();
1496 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1499 fTracker[iDet]->LoadClusters(tree);
1500 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1504 if (fTracker[iDet]->PropagateBack(esd) != 0) {
1505 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
1508 if (fCheckPointLevel > 1) {
1509 WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
1513 if (iDet > 2) { // all except ITS, TPC, TRD
1514 fTracker[iDet]->UnloadClusters();
1515 fLoader[iDet]->UnloadRecPoints();
1517 // updated PID in TPC needed by the ITS tracker -MI
1519 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1520 AliESDpid::MakePID(esd);
1522 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1525 if (fRunQA && fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1527 // write space-points to the ESD in case alignment data output
1529 if (fWriteAlignmentData)
1530 WriteAlignmentData(esd);
1532 // pass 3: TRD + TPC + ITS refit inwards
1535 for (Int_t iDet = 2; iDet >= 0; iDet--) {
1536 if (!fTracker[iDet]) continue;
1537 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
1540 if (fTracker[iDet]->RefitInward(esd) != 0) {
1541 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
1544 // run postprocessing
1545 if (fTracker[iDet]->PostProcess(esd) != 0) {
1546 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
1549 if (fCheckPointLevel > 1) {
1550 WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
1552 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1554 fTracker[iDet]->UnloadClusters();
1555 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
1556 fLoader[iDet]->UnloadRecPoints();
1557 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
1561 // Propagate track to the vertex - if not done by ITS
1563 Int_t ntracks = esd->GetNumberOfTracks();
1564 for (Int_t itrack=0; itrack<ntracks; itrack++){
1565 const Double_t kRadius = 3; // beam pipe radius
1566 const Double_t kMaxStep = 5; // max step
1567 const Double_t kMaxD = 123456; // max distance to prim vertex
1568 Double_t fieldZ = AliTracker::GetBz(); //
1569 AliESDtrack * track = esd->GetTrack(itrack);
1570 if (!track) continue;
1571 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1572 AliTracker::PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kTRUE);
1573 track->RelateToVertex(esd->GetVertex(),fieldZ, kMaxD);
1579 //_____________________________________________________________________________
1580 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
1582 // Remove the data which are not needed for the physics analysis.
1585 Int_t nTracks=esd->GetNumberOfTracks();
1586 Int_t nV0s=esd->GetNumberOfV0s();
1588 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
1590 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
1591 Bool_t rc=esd->Clean(cleanPars);
1593 nTracks=esd->GetNumberOfTracks();
1594 nV0s=esd->GetNumberOfV0s();
1596 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
1601 //_____________________________________________________________________________
1602 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
1604 // fill the event summary data
1606 AliCodeTimerAuto("")
1607 static Int_t eventNr=0;
1608 TString detStr = detectors;
1610 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1611 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1612 AliReconstructor* reconstructor = GetReconstructor(iDet);
1613 if (!reconstructor) continue;
1614 if (!ReadESD(esd, fgkDetectorName[iDet])) {
1615 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
1616 TTree* clustersTree = NULL;
1617 if (fLoader[iDet]) {
1618 fLoader[iDet]->LoadRecPoints("read");
1619 clustersTree = fLoader[iDet]->TreeR();
1620 if (!clustersTree) {
1621 AliError(Form("Can't get the %s clusters tree",
1622 fgkDetectorName[iDet]));
1623 if (fStopOnError) return kFALSE;
1626 if (fRawReader && !reconstructor->HasDigitConversion()) {
1627 reconstructor->FillESD(fRawReader, clustersTree, esd);
1629 TTree* digitsTree = NULL;
1630 if (fLoader[iDet]) {
1631 fLoader[iDet]->LoadDigits("read");
1632 digitsTree = fLoader[iDet]->TreeD();
1634 AliError(Form("Can't get the %s digits tree",
1635 fgkDetectorName[iDet]));
1636 if (fStopOnError) return kFALSE;
1639 reconstructor->FillESD(digitsTree, clustersTree, esd);
1640 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
1642 if (fLoader[iDet]) {
1643 fLoader[iDet]->UnloadRecPoints();
1646 if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
1650 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1651 AliError(Form("the following detectors were not found: %s",
1653 if (fStopOnError) return kFALSE;
1655 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
1660 //_____________________________________________________________________________
1661 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
1663 // Reads the trigger decision which is
1664 // stored in Trigger.root file and fills
1665 // the corresponding esd entries
1667 AliCodeTimerAuto("")
1669 AliInfo("Filling trigger information into the ESD");
1671 AliCentralTrigger *aCTP = NULL;
1674 AliCTPRawStream input(fRawReader);
1675 if (!input.Next()) {
1676 AliError("No valid CTP (trigger) DDL raw data is found ! The trigger information is not stored in the ESD !");
1679 esd->SetTriggerMask(input.GetClassMask());
1680 esd->SetTriggerCluster(input.GetClusterMask());
1682 aCTP = new AliCentralTrigger();
1683 TString configstr("");
1684 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
1685 AliError("No trigger configuration found in OCDB! The trigger classes information will no be stored in ESD!");
1690 AliRunLoader *runloader = AliRunLoader::GetRunLoader();
1692 if (!runloader->LoadTrigger()) {
1693 aCTP = runloader->GetTrigger();
1694 esd->SetTriggerMask(aCTP->GetClassMask());
1695 esd->SetTriggerCluster(aCTP->GetClusterMask());
1698 AliWarning("No trigger can be loaded! The trigger information is not stored in the ESD !");
1703 AliError("No run loader is available! The trigger information is not stored in the ESD !");
1708 // Now fill the trigger class names into AliESDRun object
1709 AliTriggerConfiguration *config = aCTP->GetConfiguration();
1711 AliError("No trigger configuration has been found! The trigger classes information will no be stored in ESD!");
1715 const TObjArray& classesArray = config->GetClasses();
1716 Int_t nclasses = classesArray.GetEntriesFast();
1717 for( Int_t j=0; j<nclasses; j++ ) {
1718 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At( j );
1719 Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
1720 esd->SetTriggerClass(trclass->GetName(),trindex);
1730 //_____________________________________________________________________________
1731 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
1734 // Filling information from RawReader Header
1737 AliInfo("Filling information from RawReader Header");
1738 esd->SetBunchCrossNumber(0);
1739 esd->SetOrbitNumber(0);
1740 esd->SetPeriodNumber(0);
1741 esd->SetTimeStamp(0);
1742 esd->SetEventType(0);
1743 const AliRawEventHeaderBase * eventHeader = fRawReader->GetEventHeader();
1746 const UInt_t *id = eventHeader->GetP("Id");
1747 esd->SetBunchCrossNumber((id)[1]&0x00000fff);
1748 esd->SetOrbitNumber((((id)[0]<<20)&0xf00000)|(((id)[1]>>12)&0xfffff));
1749 esd->SetPeriodNumber(((id)[0]>>4)&0x0fffffff);
1751 esd->SetTimeStamp((eventHeader->Get("Timestamp")));
1752 esd->SetEventType((eventHeader->Get("Type")));
1759 //_____________________________________________________________________________
1760 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
1762 // check whether detName is contained in detectors
1763 // if yes, it is removed from detectors
1765 // check if all detectors are selected
1766 if ((detectors.CompareTo("ALL") == 0) ||
1767 detectors.BeginsWith("ALL ") ||
1768 detectors.EndsWith(" ALL") ||
1769 detectors.Contains(" ALL ")) {
1774 // search for the given detector
1775 Bool_t result = kFALSE;
1776 if ((detectors.CompareTo(detName) == 0) ||
1777 detectors.BeginsWith(detName+" ") ||
1778 detectors.EndsWith(" "+detName) ||
1779 detectors.Contains(" "+detName+" ")) {
1780 detectors.ReplaceAll(detName, "");
1784 // clean up the detectors string
1785 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1786 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1787 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1792 //_____________________________________________________________________________
1793 Bool_t AliReconstruction::InitRunLoader()
1795 // get or create the run loader
1797 if (gAlice) delete gAlice;
1800 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
1801 // load all base libraries to get the loader classes
1802 TString libs = gSystem->GetLibraries();
1803 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1804 TString detName = fgkDetectorName[iDet];
1805 if (detName == "HLT") continue;
1806 if (libs.Contains("lib" + detName + "base.so")) continue;
1807 gSystem->Load("lib" + detName + "base.so");
1809 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
1811 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1815 fRunLoader->CdGAFile();
1816 if (gFile->GetKey(AliRunLoader::GetGAliceName())) {
1817 if (fRunLoader->LoadgAlice() == 0) {
1818 gAlice = fRunLoader->GetAliRun();
1819 AliTracker::SetFieldMap(gAlice->Field(),fUniformField);
1822 if (!gAlice && !fRawReader) {
1823 AliError(Form("no gAlice object found in file %s",
1824 fGAliceFileName.Data()));
1829 //PH This is a temporary fix to give access to the kinematics
1830 //PH that is needed for the labels of ITS clusters
1831 fRunLoader->LoadHeader();
1832 fRunLoader->LoadKinematics();
1834 } else { // galice.root does not exist
1836 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
1840 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
1841 AliConfig::GetDefaultEventFolderName(),
1844 AliError(Form("could not create run loader in file %s",
1845 fGAliceFileName.Data()));
1849 fRunLoader->MakeTree("E");
1851 while (fRawReader->NextEvent()) {
1852 fRunLoader->SetEventNumber(iEvent);
1853 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1855 fRunLoader->MakeTree("H");
1856 fRunLoader->TreeE()->Fill();
1859 fRawReader->RewindEvents();
1860 if (fNumberOfEventsPerFile > 0)
1861 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
1863 fRunLoader->SetNumberOfEventsPerFile(iEvent);
1864 fRunLoader->WriteHeader("OVERWRITE");
1865 fRunLoader->CdGAFile();
1866 fRunLoader->Write(0, TObject::kOverwrite);
1867 // AliTracker::SetFieldMap(???);
1873 //_____________________________________________________________________________
1874 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
1876 // get the reconstructor object and the loader for a detector
1878 if (fReconstructor[iDet]) return fReconstructor[iDet];
1880 // load the reconstructor object
1881 TPluginManager* pluginManager = gROOT->GetPluginManager();
1882 TString detName = fgkDetectorName[iDet];
1883 TString recName = "Ali" + detName + "Reconstructor";
1884 if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT")) return NULL;
1886 AliReconstructor* reconstructor = NULL;
1887 // first check if a plugin is defined for the reconstructor
1888 TPluginHandler* pluginHandler =
1889 pluginManager->FindHandler("AliReconstructor", detName);
1890 // if not, add a plugin for it
1891 if (!pluginHandler) {
1892 AliDebug(1, Form("defining plugin for %s", recName.Data()));
1893 TString libs = gSystem->GetLibraries();
1894 if (libs.Contains("lib" + detName + "base.so") ||
1895 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
1896 pluginManager->AddHandler("AliReconstructor", detName,
1897 recName, detName + "rec", recName + "()");
1899 pluginManager->AddHandler("AliReconstructor", detName,
1900 recName, detName, recName + "()");
1902 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
1904 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
1905 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
1907 if (reconstructor) {
1908 TObject* obj = fOptions.FindObject(detName.Data());
1909 if (obj) reconstructor->SetOption(obj->GetTitle());
1910 reconstructor->Init();
1911 fReconstructor[iDet] = reconstructor;
1914 // get or create the loader
1915 if (detName != "HLT") {
1916 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
1917 if (!fLoader[iDet]) {
1918 AliConfig::Instance()
1919 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
1921 // first check if a plugin is defined for the loader
1923 pluginManager->FindHandler("AliLoader", detName);
1924 // if not, add a plugin for it
1925 if (!pluginHandler) {
1926 TString loaderName = "Ali" + detName + "Loader";
1927 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
1928 pluginManager->AddHandler("AliLoader", detName,
1929 loaderName, detName + "base",
1930 loaderName + "(const char*, TFolder*)");
1931 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
1933 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
1935 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
1936 fRunLoader->GetEventFolder());
1938 if (!fLoader[iDet]) { // use default loader
1939 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
1941 if (!fLoader[iDet]) {
1942 AliWarning(Form("couldn't get loader for %s", detName.Data()));
1943 if (fStopOnError) return NULL;
1945 fRunLoader->AddLoader(fLoader[iDet]);
1946 fRunLoader->CdGAFile();
1947 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
1948 fRunLoader->Write(0, TObject::kOverwrite);
1953 return reconstructor;
1956 //_____________________________________________________________________________
1957 Bool_t AliReconstruction::CreateVertexer()
1959 // create the vertexer
1962 AliReconstructor* itsReconstructor = GetReconstructor(0);
1963 if (itsReconstructor) {
1964 fVertexer = itsReconstructor->CreateVertexer();
1967 AliWarning("couldn't create a vertexer for ITS");
1968 if (fStopOnError) return kFALSE;
1974 //_____________________________________________________________________________
1975 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
1977 // create the trackers
1979 TString detStr = detectors;
1980 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1981 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1982 AliReconstructor* reconstructor = GetReconstructor(iDet);
1983 if (!reconstructor) continue;
1984 TString detName = fgkDetectorName[iDet];
1985 if (detName == "HLT") {
1986 fRunHLTTracking = kTRUE;
1989 if (detName == "MUON") {
1990 fRunMuonTracking = kTRUE;
1995 fTracker[iDet] = reconstructor->CreateTracker();
1996 if (!fTracker[iDet] && (iDet < 7)) {
1997 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1998 if (fStopOnError) return kFALSE;
2000 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2006 //_____________________________________________________________________________
2007 void AliReconstruction::CleanUp(TFile* file, TFile* fileOld)
2009 // delete trackers and the run loader and close and delete the file
2011 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2012 delete fReconstructor[iDet];
2013 fReconstructor[iDet] = NULL;
2014 fLoader[iDet] = NULL;
2015 delete fTracker[iDet];
2016 fTracker[iDet] = NULL;
2017 // delete fQADataMaker[iDet];
2018 // fQADataMaker[iDet] = NULL;
2023 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2024 delete fDiamondProfile;
2025 fDiamondProfile = NULL;
2035 if (fParentRawReader) delete fParentRawReader;
2036 fParentRawReader=NULL;
2046 gSystem->Unlink("AliESDs.old.root");
2050 //_____________________________________________________________________________
2052 Bool_t AliReconstruction::ReadESD(AliESDEvent*& esd, const char* recStep) const
2054 // read the ESD event from a file
2056 if (!esd) return kFALSE;
2058 sprintf(fileName, "ESD_%d.%d_%s.root",
2059 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2060 if (gSystem->AccessPathName(fileName)) return kFALSE;
2062 AliInfo(Form("reading ESD from file %s", fileName));
2063 AliDebug(1, Form("reading ESD from file %s", fileName));
2064 TFile* file = TFile::Open(fileName);
2065 if (!file || !file->IsOpen()) {
2066 AliError(Form("opening %s failed", fileName));
2073 esd = (AliESDEvent*) file->Get("ESD");
2082 //_____________________________________________________________________________
2083 void AliReconstruction::WriteESD(AliESDEvent* esd, const char* recStep) const
2085 // write the ESD event to a file
2089 sprintf(fileName, "ESD_%d.%d_%s.root",
2090 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2092 AliDebug(1, Form("writing ESD to file %s", fileName));
2093 TFile* file = TFile::Open(fileName, "recreate");
2094 if (!file || !file->IsOpen()) {
2095 AliError(Form("opening %s failed", fileName));
2107 //_____________________________________________________________________________
2108 void AliReconstruction::ESDFile2AODFile(TFile* esdFile, TFile* aodFile)
2110 // write all files from the given esd file to an aod file
2112 // create an AliAOD object
2113 AliAODEvent *aod = new AliAODEvent();
2114 aod->CreateStdContent();
2120 TTree *aodTree = new TTree("aodTree", "AliAOD tree");
2121 aodTree->Branch(aod->GetList());
2124 TTree *t = (TTree*) esdFile->Get("esdTree");
2125 AliESDEvent *esd = new AliESDEvent();
2126 esd->ReadFromTree(t);
2128 Int_t nEvents = t->GetEntries();
2130 // set arrays and pointers
2140 // loop over events and fill them
2141 for (Int_t iEvent = 0; iEvent < nEvents; ++iEvent) {
2142 //cout << "event: " << iEvent << endl;
2143 t->GetEntry(iEvent);
2145 // Multiplicity information needed by the header (to be revised!)
2146 Int_t nTracks = esd->GetNumberOfTracks();
2147 Int_t nPosTracks = 0;
2148 for (Int_t iTrack=0; iTrack<nTracks; ++iTrack)
2149 if (esd->GetTrack(iTrack)->Charge()> 0) nPosTracks++;
2151 // Access the header
2152 AliAODHeader *header = aod->GetHeader();
2155 header->SetRunNumber (esd->GetRunNumber() );
2156 header->SetBunchCrossNumber(esd->GetBunchCrossNumber());
2157 header->SetOrbitNumber (esd->GetOrbitNumber() );
2158 header->SetPeriodNumber (esd->GetPeriodNumber() );
2159 header->SetTriggerMask (esd->GetTriggerMask() );
2160 header->SetTriggerCluster (esd->GetTriggerCluster() );
2161 header->SetEventType (esd->GetEventType() );
2162 header->SetMagneticField (esd->GetMagneticField() );
2163 header->SetZDCN1Energy (esd->GetZDCN1Energy() );
2164 header->SetZDCP1Energy (esd->GetZDCP1Energy() );
2165 header->SetZDCN2Energy (esd->GetZDCN2Energy() );
2166 header->SetZDCP2Energy (esd->GetZDCP2Energy() );
2167 header->SetZDCEMEnergy (esd->GetZDCEMEnergy(0),esd->GetZDCEMEnergy(1));
2168 header->SetRefMultiplicity (nTracks);
2169 header->SetRefMultiplicityPos(nPosTracks);
2170 header->SetRefMultiplicityNeg(nTracks - nPosTracks);
2171 header->SetMuonMagFieldScale(-999.); // FIXME
2172 header->SetCentrality(-999.); // FIXME
2174 Int_t nV0s = esd->GetNumberOfV0s();
2175 Int_t nCascades = esd->GetNumberOfCascades();
2176 Int_t nKinks = esd->GetNumberOfKinks();
2177 Int_t nVertices = nV0s + 2*nCascades /*could lead to two vertices, one V0 and the Xi */+ nKinks + 1 /* = prim. vtx*/;
2179 Int_t nCaloClus = esd->GetNumberOfCaloClusters();
2181 Int_t nPmdClus = esd->GetNumberOfPmdTracks();
2183 aod->ResetStd(nTracks, nVertices, nV0s+nCascades, nJets, nCaloClus, nFmdClus, nPmdClus);
2185 // Array to take into account the tracks already added to the AOD
2186 Bool_t * usedTrack = NULL;
2188 usedTrack = new Bool_t[nTracks];
2189 for (Int_t iTrack=0; iTrack<nTracks; ++iTrack) usedTrack[iTrack]=kFALSE;
2191 // Array to take into account the V0s already added to the AOD
2192 Bool_t * usedV0 = NULL;
2194 usedV0 = new Bool_t[nV0s];
2195 for (Int_t iV0=0; iV0<nV0s; ++iV0) usedV0[iV0]=kFALSE;
2197 // Array to take into account the kinks already added to the AOD
2198 Bool_t * usedKink = NULL;
2200 usedKink = new Bool_t[nKinks];
2201 for (Int_t iKink=0; iKink<nKinks; ++iKink) usedKink[iKink]=kFALSE;
2204 // Access to the AOD container of vertices
2205 TClonesArray &vertices = *(aod->GetVertices());
2208 // Access to the AOD container of tracks
2209 TClonesArray &tracks = *(aod->GetTracks());
2212 // Access to the AOD container of V0s
2213 TClonesArray &V0s = *(aod->GetV0s());
2216 // Add primary vertex. The primary tracks will be defined
2217 // after the loops on the composite objects (V0, cascades, kinks)
2218 const AliESDVertex *vtx = esd->GetPrimaryVertex();
2220 vtx->GetXYZ(pos); // position
2221 vtx->GetCovMatrix(covVtx); //covariance matrix
2223 AliAODVertex * primary = new(vertices[jVertices++])
2224 AliAODVertex(pos, covVtx, vtx->GetChi2toNDF(), NULL, -1, AliAODVertex::kPrimary);
2227 AliAODTrack *aodTrack = 0x0;
2229 // Create vertices starting from the most complex objects
2232 for (Int_t nCascade = 0; nCascade < nCascades; ++nCascade) {
2233 AliESDcascade *cascade = esd->GetCascade(nCascade);
2235 cascade->GetXYZ(pos[0], pos[1], pos[2]);
2236 cascade->GetPosCovXi(covVtx);
2238 // Add the cascade vertex
2239 AliAODVertex * vcascade = new(vertices[jVertices++]) AliAODVertex(pos,
2241 cascade->GetChi2Xi(), // = chi2/NDF since NDF = 2*2-3
2244 AliAODVertex::kCascade);
2246 primary->AddDaughter(vcascade); // the cascade 'particle' (represented by a vertex) is added as a daughter to the primary vertex
2248 // Add the V0 from the cascade. The ESD class have to be optimized...
2249 // Now we have to search for the corresponding V0 in the list of V0s
2250 // using the indeces of the positive and negative tracks
2252 Int_t posFromV0 = cascade->GetPindex();
2253 Int_t negFromV0 = cascade->GetNindex();
2256 AliESDv0 * v0 = 0x0;
2259 for (Int_t iV0=0; iV0<nV0s; ++iV0) {
2261 v0 = esd->GetV0(iV0);
2262 Int_t posV0 = v0->GetPindex();
2263 Int_t negV0 = v0->GetNindex();
2265 if (posV0==posFromV0 && negV0==negFromV0) {
2271 AliAODVertex * vV0FromCascade = 0x0;
2273 if (indV0>-1 && !usedV0[indV0]) {
2275 // the V0 exists in the array of V0s and is not used
2277 usedV0[indV0] = kTRUE;
2279 v0->GetXYZ(pos[0], pos[1], pos[2]);
2280 v0->GetPosCov(covVtx);
2282 vV0FromCascade = new(vertices[jVertices++]) AliAODVertex(pos,
2284 v0->GetChi2V0(), // = chi2/NDF since NDF = 2*2-3
2290 // the V0 doesn't exist in the array of V0s or was used
2291 cerr << "Error: event " << iEvent << " cascade " << nCascade
2292 << " The V0 " << indV0
2293 << " doesn't exist in the array of V0s or was used!" << endl;
2295 cascade->GetXYZ(pos[0], pos[1], pos[2]);
2296 cascade->GetPosCov(covVtx);
2298 vV0FromCascade = new(vertices[jVertices++]) AliAODVertex(pos,
2300 v0->GetChi2V0(), // = chi2/NDF since NDF = 2*2-3
2304 vcascade->AddDaughter(vV0FromCascade);
2308 // Add the positive tracks from the V0
2310 if (! usedTrack[posFromV0]) {
2312 usedTrack[posFromV0] = kTRUE;
2314 AliESDtrack *esdTrack = esd->GetTrack(posFromV0);
2315 esdTrack->GetPxPyPz(p_pos);
2316 esdTrack->GetXYZ(pos);
2317 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2318 esdTrack->GetESDpid(pid);
2320 vV0FromCascade->AddDaughter(aodTrack =
2321 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2322 esdTrack->GetLabel(),
2328 (Short_t)esdTrack->Charge(),
2329 esdTrack->GetITSClusterMap(),
2332 kTRUE, // check if this is right
2333 kFALSE, // check if this is right
2334 AliAODTrack::kSecondary)
2336 aodTrack->ConvertAliPIDtoAODPID();
2339 cerr << "Error: event " << iEvent << " cascade " << nCascade
2340 << " track " << posFromV0 << " has already been used!" << endl;
2343 // Add the negative tracks from the V0
2345 if (!usedTrack[negFromV0]) {
2347 usedTrack[negFromV0] = kTRUE;
2349 AliESDtrack *esdTrack = esd->GetTrack(negFromV0);
2350 esdTrack->GetPxPyPz(p_neg);
2351 esdTrack->GetXYZ(pos);
2352 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2353 esdTrack->GetESDpid(pid);
2355 vV0FromCascade->AddDaughter(aodTrack =
2356 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2357 esdTrack->GetLabel(),
2363 (Short_t)esdTrack->Charge(),
2364 esdTrack->GetITSClusterMap(),
2367 kTRUE, // check if this is right
2368 kFALSE, // check if this is right
2369 AliAODTrack::kSecondary)
2371 aodTrack->ConvertAliPIDtoAODPID();
2374 cerr << "Error: event " << iEvent << " cascade " << nCascade
2375 << " track " << negFromV0 << " has already been used!" << endl;
2378 // add it to the V0 array as well
2379 Double_t d0[2] = { -999., -99.};
2380 // counting is probably wrong
2381 new(V0s[jV0s++]) AliAODv0(vV0FromCascade, -999., -99., p_pos, p_neg, d0); // to be refined
2383 // Add the bachelor track from the cascade
2385 Int_t bachelor = cascade->GetBindex();
2387 if(!usedTrack[bachelor]) {
2389 usedTrack[bachelor] = kTRUE;
2391 AliESDtrack *esdTrack = esd->GetTrack(bachelor);
2392 esdTrack->GetPxPyPz(p);
2393 esdTrack->GetXYZ(pos);
2394 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2395 esdTrack->GetESDpid(pid);
2397 vcascade->AddDaughter(aodTrack =
2398 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2399 esdTrack->GetLabel(),
2405 (Short_t)esdTrack->Charge(),
2406 esdTrack->GetITSClusterMap(),
2409 kTRUE, // check if this is right
2410 kFALSE, // check if this is right
2411 AliAODTrack::kSecondary)
2413 aodTrack->ConvertAliPIDtoAODPID();
2416 cerr << "Error: event " << iEvent << " cascade " << nCascade
2417 << " track " << bachelor << " has already been used!" << endl;
2420 // Add the primary track of the cascade (if any)
2422 } // end of the loop on cascades
2426 for (Int_t nV0 = 0; nV0 < nV0s; ++nV0) {
2428 if (usedV0[nV0]) continue; // skip if aready added to the AOD
2430 AliESDv0 *v0 = esd->GetV0(nV0);
2432 v0->GetXYZ(pos[0], pos[1], pos[2]);
2433 v0->GetPosCov(covVtx);
2435 AliAODVertex * vV0 =
2436 new(vertices[jVertices++]) AliAODVertex(pos,
2438 v0->GetChi2V0(), // = chi2/NDF since NDF = 2*2-3
2442 primary->AddDaughter(vV0);
2444 Int_t posFromV0 = v0->GetPindex();
2445 Int_t negFromV0 = v0->GetNindex();
2447 // Add the positive tracks from the V0
2449 if (!usedTrack[posFromV0]) {
2451 usedTrack[posFromV0] = kTRUE;
2453 AliESDtrack *esdTrack = esd->GetTrack(posFromV0);
2454 esdTrack->GetPxPyPz(p_pos);
2455 esdTrack->GetXYZ(pos);
2456 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2457 esdTrack->GetESDpid(pid);
2459 vV0->AddDaughter(aodTrack =
2460 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2461 esdTrack->GetLabel(),
2467 (Short_t)esdTrack->Charge(),
2468 esdTrack->GetITSClusterMap(),
2471 kTRUE, // check if this is right
2472 kFALSE, // check if this is right
2473 AliAODTrack::kSecondary)
2475 aodTrack->ConvertAliPIDtoAODPID();
2478 cerr << "Error: event " << iEvent << " V0 " << nV0
2479 << " track " << posFromV0 << " has already been used!" << endl;
2482 // Add the negative tracks from the V0
2484 if (!usedTrack[negFromV0]) {
2486 usedTrack[negFromV0] = kTRUE;
2488 AliESDtrack *esdTrack = esd->GetTrack(negFromV0);
2489 esdTrack->GetPxPyPz(p_neg);
2490 esdTrack->GetXYZ(pos);
2491 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2492 esdTrack->GetESDpid(pid);
2494 vV0->AddDaughter(aodTrack =
2495 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2496 esdTrack->GetLabel(),
2502 (Short_t)esdTrack->Charge(),
2503 esdTrack->GetITSClusterMap(),
2506 kTRUE, // check if this is right
2507 kFALSE, // check if this is right
2508 AliAODTrack::kSecondary)
2510 aodTrack->ConvertAliPIDtoAODPID();
2513 cerr << "Error: event " << iEvent << " V0 " << nV0
2514 << " track " << negFromV0 << " has already been used!" << endl;
2517 // add it to the V0 array as well
2518 Double_t d0[2] = { 999., 99.};
2519 new(V0s[jV0s++]) AliAODv0(vV0, 999., 99., p_pos, p_neg, d0); // to be refined
2522 // end of the loop on V0s
2524 // Kinks: it is a big mess the access to the information in the kinks
2525 // The loop is on the tracks in order to find the mother and daugther of each kink
2528 for (Int_t iTrack=0; iTrack<nTracks; ++iTrack) {
2530 AliESDtrack * esdTrack = esd->GetTrack(iTrack);
2532 Int_t ikink = esdTrack->GetKinkIndex(0);
2535 // Negative kink index: mother, positive: daughter
2537 // Search for the second track of the kink
2539 for (Int_t jTrack = iTrack+1; jTrack<nTracks; ++jTrack) {
2541 AliESDtrack * esdTrack1 = esd->GetTrack(jTrack);
2543 Int_t jkink = esdTrack1->GetKinkIndex(0);
2545 if ( TMath::Abs(ikink)==TMath::Abs(jkink) ) {
2547 // The two tracks are from the same kink
2549 if (usedKink[TMath::Abs(ikink)-1]) continue; // skip used kinks
2552 Int_t idaughter = -1;
2554 if (ikink<0 && jkink>0) {
2559 else if (ikink>0 && jkink<0) {
2565 cerr << "Error: Wrong combination of kink indexes: "
2566 << ikink << " " << jkink << endl;
2570 // Add the mother track
2572 AliAODTrack * mother = NULL;
2574 if (!usedTrack[imother]) {
2576 usedTrack[imother] = kTRUE;
2578 AliESDtrack *esdTrack = esd->GetTrack(imother);
2579 esdTrack->GetPxPyPz(p);
2580 esdTrack->GetXYZ(pos);
2581 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2582 esdTrack->GetESDpid(pid);
2585 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2586 esdTrack->GetLabel(),
2592 (Short_t)esdTrack->Charge(),
2593 esdTrack->GetITSClusterMap(),
2596 kTRUE, // check if this is right
2597 kTRUE, // check if this is right
2598 AliAODTrack::kPrimary);
2599 primary->AddDaughter(mother);
2600 mother->ConvertAliPIDtoAODPID();
2603 cerr << "Error: event " << iEvent << " kink " << TMath::Abs(ikink)-1
2604 << " track " << imother << " has already been used!" << endl;
2607 // Add the kink vertex
2608 AliESDkink * kink = esd->GetKink(TMath::Abs(ikink)-1);
2610 AliAODVertex * vkink =
2611 new(vertices[jVertices++]) AliAODVertex(kink->GetPosition(),
2615 esdTrack->GetID(), // This is the track ID of the mother's track!
2616 AliAODVertex::kKink);
2617 // Add the daughter track
2619 AliAODTrack * daughter = NULL;
2621 if (!usedTrack[idaughter]) {
2623 usedTrack[idaughter] = kTRUE;
2625 AliESDtrack *esdTrack = esd->GetTrack(idaughter);
2626 esdTrack->GetPxPyPz(p);
2627 esdTrack->GetXYZ(pos);
2628 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2629 esdTrack->GetESDpid(pid);
2632 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2633 esdTrack->GetLabel(),
2639 (Short_t)esdTrack->Charge(),
2640 esdTrack->GetITSClusterMap(),
2643 kTRUE, // check if this is right
2644 kTRUE, // check if this is right
2645 AliAODTrack::kPrimary);
2646 vkink->AddDaughter(daughter);
2647 daughter->ConvertAliPIDtoAODPID();
2650 cerr << "Error: event " << iEvent << " kink " << TMath::Abs(ikink)-1
2651 << " track " << idaughter << " has already been used!" << endl;
2657 vertices.Expand(jVertices);
2659 // Tracks (primary and orphan)
2660 for (Int_t nTrack = 0; nTrack < nTracks; ++nTrack) {
2662 if (usedTrack[nTrack]) continue;
2664 AliESDtrack *esdTrack = esd->GetTrack(nTrack);
2665 esdTrack->GetPxPyPz(p);
2666 esdTrack->GetXYZ(pos);
2667 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2668 esdTrack->GetESDpid(pid);
2670 Float_t impactXY, impactZ;
2672 esdTrack->GetImpactParameters(impactXY,impactZ);
2675 // track inside the beam pipe
2677 primary->AddDaughter(aodTrack =
2678 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2679 esdTrack->GetLabel(),
2685 (Short_t)esdTrack->Charge(),
2686 esdTrack->GetITSClusterMap(),
2689 kTRUE, // check if this is right
2690 kTRUE, // check if this is right
2691 AliAODTrack::kPrimary)
2693 aodTrack->ConvertAliPIDtoAODPID();
2696 // outside the beam pipe: orphan track
2697 // Don't write them anymore!
2700 } // end of loop on tracks
2703 Int_t nMuTracks = esd->GetNumberOfMuonTracks();
2704 for (Int_t nMuTrack = 0; nMuTrack < nMuTracks; ++nMuTrack) {
2706 AliESDMuonTrack *esdMuTrack = esd->GetMuonTrack(nMuTrack);
2707 p[0] = esdMuTrack->Px();
2708 p[1] = esdMuTrack->Py();
2709 p[2] = esdMuTrack->Pz();
2710 pos[0] = primary->GetX();
2711 pos[1] = primary->GetY();
2712 pos[2] = primary->GetZ();
2714 // has to be changed once the muon pid is provided by the ESD
2715 for (Int_t i = 0; i < 10; pid[i++] = 0.); pid[AliAODTrack::kMuon]=1.;
2717 primary->AddDaughter(aodTrack =
2718 new(tracks[jTracks++]) AliAODTrack(0, // no ID provided
2719 0, // no label provided
2724 NULL, // no covariance matrix provided
2725 esdMuTrack->Charge(),
2726 0, // ITSClusterMap is set below
2729 kFALSE, // muon tracks are not used to fit the primary vtx
2730 kFALSE, // not used for vertex fit
2731 AliAODTrack::kPrimary)
2734 aodTrack->SetHitsPatternInTrigCh(esdMuTrack->GetHitsPatternInTrigCh());
2735 Int_t track2Trigger = esdMuTrack->GetMatchTrigger();
2736 aodTrack->SetMatchTrigger(track2Trigger);
2738 aodTrack->SetChi2MatchTrigger(esdMuTrack->GetChi2MatchTrigger());
2740 aodTrack->SetChi2MatchTrigger(0.);
2742 tracks.Expand(jTracks); // remove 'empty slots' due to unwritten tracks
2744 // Access to the AOD container of PMD clusters
2745 TClonesArray &pmdClusters = *(aod->GetPmdClusters());
2746 Int_t jPmdClusters=0;
2748 for (Int_t iPmd = 0; iPmd < nPmdClus; ++iPmd) {
2749 // file pmd clusters, to be revised!
2750 AliESDPmdTrack *pmdTrack = esd->GetPmdTrack(iPmd);
2753 Double_t pos[3] = { pmdTrack->GetClusterX(), pmdTrack->GetClusterY(), pmdTrack->GetClusterZ() };
2754 Double_t pid[9] = { 0., 0., 0., 0., 0., 0., 0., 0., 0. }; // to be revised!
2756 // assoc cluster not set
2757 new(pmdClusters[jPmdClusters++]) AliAODPmdCluster(iPmd, nLabel, label, pmdTrack->GetClusterADC(), pos, pid);
2760 // Access to the AOD container of clusters
2761 TClonesArray &caloClusters = *(aod->GetCaloClusters());
2764 for (Int_t iClust=0; iClust<nCaloClus; ++iClust) {
2766 AliESDCaloCluster * cluster = esd->GetCaloCluster(iClust);
2768 Int_t id = cluster->GetID();
2771 Float_t energy = cluster->E();
2772 cluster->GetPosition(posF);
2773 Char_t ttype=AliAODCluster::kUndef;
2775 if (cluster->GetClusterType() == AliESDCaloCluster::kPHOSCluster) {
2776 ttype=AliAODCluster::kPHOSNeutral;
2778 else if (cluster->GetClusterType() == AliESDCaloCluster::kEMCALClusterv1) {
2779 ttype = AliAODCluster::kEMCALClusterv1;
2783 AliAODCaloCluster *caloCluster = new(caloClusters[jClusters++]) AliAODCaloCluster(id,
2791 caloCluster->SetCaloCluster(); // to be refined!
2794 caloClusters.Expand(jClusters); // resize TObjArray to 'remove' slots for pseudo clusters
2795 // end of loop on calo clusters
2797 // fill EMCAL cell info
2798 if (esd->GetEMCALCells()) { // protection against missing ESD information
2799 AliESDCaloCells &esdEMcells = *(esd->GetEMCALCells());
2800 Int_t nEMcell = esdEMcells.GetNumberOfCells() ;
2802 AliAODCaloCells &aodEMcells = *(aod->GetEMCALCells());
2803 aodEMcells.CreateContainer(nEMcell);
2804 aodEMcells.SetType(AliAODCaloCells::kEMCAL);
2805 for (Int_t iCell = 0; iCell < nEMcell; iCell++) {
2806 aodEMcells.SetCell(iCell,esdEMcells.GetCellNumber(iCell),esdEMcells.GetAmplitude(iCell));
2811 // fill PHOS cell info
2812 if (esd->GetPHOSCells()) { // protection against missing ESD information
2813 AliESDCaloCells &esdPHcells = *(esd->GetPHOSCells());
2814 Int_t nPHcell = esdPHcells.GetNumberOfCells() ;
2816 AliAODCaloCells &aodPHcells = *(aod->GetPHOSCells());
2817 aodPHcells.CreateContainer(nPHcell);
2818 aodPHcells.SetType(AliAODCaloCells::kPHOS);
2819 for (Int_t iCell = 0; iCell < nPHcell; iCell++) {
2820 aodPHcells.SetCell(iCell,esdPHcells.GetCellNumber(iCell),esdPHcells.GetAmplitude(iCell));
2826 AliAODTracklets &SPDTracklets = *(aod->GetTracklets());
2827 const AliMultiplicity *mult = esd->GetMultiplicity();
2829 if (mult->GetNumberOfTracklets()>0) {
2830 SPDTracklets.CreateContainer(mult->GetNumberOfTracklets());
2832 for (Int_t n=0; n<mult->GetNumberOfTracklets(); n++) {
2833 SPDTracklets.SetTracklet(n, mult->GetTheta(n), mult->GetPhi(n), mult->GetDeltaPhi(n), mult->GetLabel(n));
2837 Printf("ERROR: AliMultiplicity could not be retrieved from ESD");
2840 delete [] usedTrack;
2844 // fill the tree for this event
2846 } // end of event loop
2848 aodTree->GetUserInfo()->Add(aod);
2850 // write the tree to the specified file
2851 aodFile = aodTree->GetCurrentFile();
2858 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2860 // Write space-points which are then used in the alignment procedures
2861 // For the moment only ITS, TRD and TPC
2863 // Load TOF clusters
2865 fLoader[3]->LoadRecPoints("read");
2866 TTree* tree = fLoader[3]->TreeR();
2868 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[3]));
2871 fTracker[3]->LoadClusters(tree);
2873 Int_t ntracks = esd->GetNumberOfTracks();
2874 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2876 AliESDtrack *track = esd->GetTrack(itrack);
2879 for (Int_t iDet = 3; iDet >= 0; iDet--)
2880 nsp += track->GetNcls(iDet);
2882 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2883 track->SetTrackPointArray(sp);
2885 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2886 AliTracker *tracker = fTracker[iDet];
2887 if (!tracker) continue;
2888 Int_t nspdet = track->GetNcls(iDet);
2889 if (nspdet <= 0) continue;
2890 track->GetClusters(iDet,idx);
2894 while (isp2 < nspdet) {
2896 TString dets = fgkDetectorName[iDet];
2897 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2898 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2899 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2900 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2901 isvalid = tracker->GetTrackPointTrackingError(idx[isp2],p,track);
2903 isvalid = tracker->GetTrackPoint(idx[isp2],p);
2906 const Int_t kNTPCmax = 159;
2907 if (iDet==1 && isp2>kNTPCmax) break; // to be fixed
2908 if (!isvalid) continue;
2909 sp->AddPoint(isptrack,&p); isptrack++; isp++;
2915 fTracker[3]->UnloadClusters();
2916 fLoader[3]->UnloadRecPoints();
2920 //_____________________________________________________________________________
2921 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2923 // The method reads the raw-data error log
2924 // accumulated within the rawReader.
2925 // It extracts the raw-data errors related to
2926 // the current event and stores them into
2927 // a TClonesArray inside the esd object.
2929 if (!fRawReader) return;
2931 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2933 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2935 if (iEvent != log->GetEventNumber()) continue;
2937 esd->AddRawDataErrorLog(log);
2942 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString fName){
2943 // Dump a file content into a char in TNamed
2945 in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2946 Int_t kBytes = (Int_t)in.tellg();
2947 printf("Size: %d \n",kBytes);
2950 char* memblock = new char [kBytes];
2951 in.seekg (0, ios::beg);
2952 in.read (memblock, kBytes);
2954 TString fData(memblock,kBytes);
2955 fn = new TNamed(fName,fData);
2956 printf("fData Size: %d \n",fData.Sizeof());
2957 printf("fName Size: %d \n",fName.Sizeof());
2958 printf("fn Size: %d \n",fn->Sizeof());
2962 AliInfo(Form("Could not Open %s\n",fPath.Data()));
2968 void AliReconstruction::TNamedToFile(TTree* fTree, TString fName){
2969 // This is not really needed in AliReconstruction at the moment
2970 // but can serve as a template
2972 TList *fList = fTree->GetUserInfo();
2973 TNamed *fn = (TNamed*)fList->FindObject(fName.Data());
2974 printf("fn Size: %d \n",fn->Sizeof());
2976 TString fTmp(fn->GetName()); // to be 100% sure in principle fName also works
2977 const char* cdata = fn->GetTitle();
2978 printf("fTmp Size %d\n",fTmp.Sizeof());
2980 int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2981 printf("calculated size %d\n",size);
2982 ofstream out(fName.Data(),ios::out | ios::binary);
2983 out.write(cdata,size);
2988 //_____________________________________________________________________________
2989 AliQADataMakerRec * AliReconstruction::GetQADataMaker(Int_t iDet)
2991 // get the quality assurance data maker object and the loader for a detector
2993 if (fQADataMaker[iDet])
2994 return fQADataMaker[iDet];
2996 AliQADataMakerRec * qadm = NULL;
2997 if (iDet == fgkNDetectors) { //Global QA
2998 qadm = new AliGlobalQADataMaker();
2999 fQADataMaker[iDet] = qadm;
3003 // load the QA data maker object
3004 TPluginManager* pluginManager = gROOT->GetPluginManager();
3005 TString detName = fgkDetectorName[iDet];
3006 TString qadmName = "Ali" + detName + "QADataMakerRec";
3007 if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT"))
3010 // first check if a plugin is defined for the quality assurance data maker
3011 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
3012 // if not, add a plugin for it
3013 if (!pluginHandler) {
3014 AliDebug(1, Form("defining plugin for %s", qadmName.Data()));
3015 TString libs = gSystem->GetLibraries();
3016 if (libs.Contains("lib" + detName + "base.so") ||
3017 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
3018 pluginManager->AddHandler("AliQADataMakerRec", detName,
3019 qadmName, detName + "qadm", qadmName + "()");
3021 pluginManager->AddHandler("AliQADataMakerRec", detName,
3022 qadmName, detName, qadmName + "()");
3024 pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
3026 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
3027 qadm = (AliQADataMakerRec *) pluginHandler->ExecPlugin(0);
3030 fQADataMaker[iDet] = qadm;
3035 //_____________________________________________________________________________
3036 Bool_t AliReconstruction::RunQA(const char* detectors, AliESDEvent *& esd)
3038 // run the Quality Assurance data producer
3040 AliCodeTimerAuto("")
3041 TString detStr = detectors;
3042 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3043 if (!IsSelected(fgkDetectorName[iDet], detStr))
3045 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
3048 AliCodeTimerStart(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
3049 AliInfo(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
3051 qadm->Exec(AliQA::kESDS, esd) ;
3054 AliCodeTimerStop(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
3056 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
3057 AliError(Form("the following detectors were not found: %s",
3067 //_____________________________________________________________________________
3068 void AliReconstruction::CheckQA()
3070 // check the QA of SIM for this run and remove the detectors
3071 // with status Fatal
3073 TString newRunLocalReconstruction ;
3074 TString newRunTracking ;
3075 TString newFillESD ;
3077 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
3078 TString detName(AliQA::GetDetName(iDet)) ;
3079 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX(iDet)) ;
3080 if ( qa->IsSet(AliQA::DETECTORINDEX(iDet), AliQA::kSIM, AliQA::kFATAL)) {
3081 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
3083 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
3084 fRunLocalReconstruction.Contains("ALL") ) {
3085 newRunLocalReconstruction += detName ;
3086 newRunLocalReconstruction += " " ;
3088 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
3089 fRunTracking.Contains("ALL") ) {
3090 newRunTracking += detName ;
3091 newRunTracking += " " ;
3093 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
3094 fFillESD.Contains("ALL") ) {
3095 newFillESD += detName ;
3100 fRunLocalReconstruction = newRunLocalReconstruction ;
3101 fRunTracking = newRunTracking ;
3102 fFillESD = newFillESD ;
3105 //_____________________________________________________________________________
3106 Int_t AliReconstruction::GetDetIndex(const char* detector)
3108 // return the detector index corresponding to detector
3110 for (index = 0; index < fgkNDetectors ; index++) {
3111 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
3116 //_____________________________________________________________________________
3117 Bool_t AliReconstruction::FinishPlaneEff() {
3119 // Here execute all the necessary operationis, at the end of the tracking phase,
3120 // in case that evaluation of PlaneEfficiencies was required for some detector.
3121 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
3123 // This Preliminary version works only FOR ITS !!!!!
3124 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3127 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3130 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3131 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
3132 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3133 if(fTracker[iDet]) ret=fTracker[iDet]->GetPlaneEff()->WriteIntoCDB();
3137 //_____________________________________________________________________________
3138 Bool_t AliReconstruction::InitPlaneEff() {
3140 // Here execute all the necessary operations, before of the tracking phase,
3141 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
3142 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
3143 // which should be updated/recalculated.
3145 // This Preliminary version will work only FOR ITS !!!!!
3146 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3149 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3151 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));