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 // For debug purposes the method SetCheckPointLevel can be used. If the //
105 // argument is greater than 0, files with ESD events will be written after //
106 // selected steps of the reconstruction for each event: //
107 // level 1: after tracking and after filling of ESD (final) //
108 // level 2: in addition after each tracking step //
109 // level 3: in addition after the filling of ESD for each detector //
110 // If a final check point file exists for an event, this event will be //
111 // skipped in the reconstruction. The tracking and the filling of ESD for //
112 // a detector will be skipped as well, if the corresponding check point //
113 // file exists. The ESD event will then be loaded from the file instead. //
115 ///////////////////////////////////////////////////////////////////////////////
122 #include <TPluginManager.h>
123 #include <TGeoManager.h>
124 #include <TLorentzVector.h>
128 #include "AliReconstruction.h"
129 #include "AliCodeTimer.h"
130 #include "AliReconstructor.h"
132 #include "AliRunLoader.h"
134 #include "AliRawReaderFile.h"
135 #include "AliRawReaderDate.h"
136 #include "AliRawReaderRoot.h"
137 #include "AliRawEventHeaderBase.h"
138 #include "AliESDEvent.h"
139 #include "AliESDMuonTrack.h"
140 #include "AliESDfriend.h"
141 #include "AliESDVertex.h"
142 #include "AliESDcascade.h"
143 #include "AliESDkink.h"
144 #include "AliESDtrack.h"
145 #include "AliESDCaloCluster.h"
146 #include "AliESDCaloCells.h"
147 #include "AliMultiplicity.h"
148 #include "AliTracker.h"
149 #include "AliVertexer.h"
150 #include "AliVertexerTracks.h"
151 #include "AliV0vertexer.h"
152 #include "AliCascadeVertexer.h"
153 #include "AliHeader.h"
154 #include "AliGenEventHeader.h"
156 #include "AliESDpid.h"
157 #include "AliESDtrack.h"
158 #include "AliESDPmdTrack.h"
160 #include "AliESDTagCreator.h"
161 #include "AliAODTagCreator.h"
163 #include "AliGeomManager.h"
164 #include "AliTrackPointArray.h"
165 #include "AliCDBManager.h"
166 #include "AliCDBStorage.h"
167 #include "AliCDBEntry.h"
168 #include "AliAlignObj.h"
170 #include "AliCentralTrigger.h"
171 #include "AliTriggerConfiguration.h"
172 #include "AliTriggerClass.h"
173 #include "AliCTPRawStream.h"
175 #include "AliAODEvent.h"
176 #include "AliAODHeader.h"
177 #include "AliAODTrack.h"
178 #include "AliAODVertex.h"
179 #include "AliAODv0.h"
180 #include "AliAODJet.h"
181 #include "AliAODCaloCells.h"
182 #include "AliAODCaloCluster.h"
183 #include "AliAODPmdCluster.h"
184 #include "AliAODFmdCluster.h"
185 #include "AliAODTracklets.h"
187 #include "AliQADataMakerRec.h"
188 #include "AliGlobalQADataMaker.h"
190 #include "AliQADataMakerSteer.h"
192 #include "AliPlaneEff.h"
194 #include "AliSysInfo.h" // memory snapshots
197 ClassImp(AliReconstruction)
200 //_____________________________________________________________________________
201 const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
203 //_____________________________________________________________________________
204 AliReconstruction::AliReconstruction(const char* gAliceFilename,
205 const char* name, const char* title) :
208 fUniformField(kTRUE),
209 fRunVertexFinder(kTRUE),
210 fRunHLTTracking(kFALSE),
211 fRunMuonTracking(kFALSE),
213 fRunCascadeFinder(kTRUE),
214 fStopOnError(kFALSE),
215 fWriteAlignmentData(kFALSE),
216 fWriteESDfriend(kFALSE),
218 fFillTriggerESD(kTRUE),
226 fRunLocalReconstruction("ALL"),
229 fUseTrackingErrorsForAlignment(""),
230 fGAliceFileName(gAliceFilename),
235 fNumberOfEventsPerFile(1),
238 fLoadAlignFromCDB(kTRUE),
239 fLoadAlignData("ALL"),
246 fDiamondProfile(NULL),
247 fMeanVertexConstraint(kTRUE),
251 fAlignObjArray(NULL),
254 fInitCDBCalled(kFALSE),
255 fSetRunNumberFromDataCalled(kFALSE),
257 fRunGlobalQA(kFALSE),
262 // create reconstruction object with default parameters
264 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
265 fReconstructor[iDet] = NULL;
266 fLoader[iDet] = NULL;
267 fTracker[iDet] = NULL;
268 fQADataMaker[iDet] = NULL;
269 fQACycles[iDet] = 999999;
271 fQADataMaker[fgkNDetectors]=NULL; //Global QA
275 //_____________________________________________________________________________
276 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
279 fUniformField(rec.fUniformField),
280 fRunVertexFinder(rec.fRunVertexFinder),
281 fRunHLTTracking(rec.fRunHLTTracking),
282 fRunMuonTracking(rec.fRunMuonTracking),
283 fRunV0Finder(rec.fRunV0Finder),
284 fRunCascadeFinder(rec.fRunCascadeFinder),
285 fStopOnError(rec.fStopOnError),
286 fWriteAlignmentData(rec.fWriteAlignmentData),
287 fWriteESDfriend(rec.fWriteESDfriend),
288 fWriteAOD(rec.fWriteAOD),
289 fFillTriggerESD(rec.fFillTriggerESD),
291 fCleanESD(rec.fCleanESD),
292 fV0DCAmax(rec.fV0DCAmax),
293 fV0CsPmin(rec.fV0CsPmin),
297 fRunLocalReconstruction(rec.fRunLocalReconstruction),
298 fRunTracking(rec.fRunTracking),
299 fFillESD(rec.fFillESD),
300 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
301 fGAliceFileName(rec.fGAliceFileName),
303 fEquipIdMap(rec.fEquipIdMap),
304 fFirstEvent(rec.fFirstEvent),
305 fLastEvent(rec.fLastEvent),
306 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
309 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
310 fLoadAlignData(rec.fLoadAlignData),
311 fESDPar(rec.fESDPar),
317 fDiamondProfile(NULL),
318 fMeanVertexConstraint(rec.fMeanVertexConstraint),
322 fAlignObjArray(rec.fAlignObjArray),
323 fCDBUri(rec.fCDBUri),
325 fInitCDBCalled(rec.fInitCDBCalled),
326 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
328 fRunGlobalQA(rec.fRunGlobalQA),
329 fInLoopQA(rec.fInLoopQA),
330 fRunPlaneEff(rec.fRunPlaneEff)
334 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
335 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
337 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
338 fReconstructor[iDet] = NULL;
339 fLoader[iDet] = NULL;
340 fTracker[iDet] = NULL;
341 fQADataMaker[iDet] = NULL;
342 fQACycles[iDet] = rec.fQACycles[iDet];
344 fQADataMaker[fgkNDetectors]=NULL; //Global QA
345 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
346 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
350 //_____________________________________________________________________________
351 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
353 // assignment operator
355 this->~AliReconstruction();
356 new(this) AliReconstruction(rec);
360 //_____________________________________________________________________________
361 AliReconstruction::~AliReconstruction()
367 fSpecCDBUri.Delete();
369 AliCodeTimer::Instance()->Print();
372 //_____________________________________________________________________________
373 void AliReconstruction::InitCDB()
375 // activate a default CDB storage
376 // First check if we have any CDB storage set, because it is used
377 // to retrieve the calibration and alignment constants
379 if (fInitCDBCalled) return;
380 fInitCDBCalled = kTRUE;
382 AliCDBManager* man = AliCDBManager::Instance();
383 if (man->IsDefaultStorageSet())
385 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
386 AliWarning("Default CDB storage has been already set !");
387 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
388 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
389 fCDBUri = man->GetDefaultStorage()->GetURI();
392 if (fCDBUri.Length() > 0)
394 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
395 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
396 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
398 fCDBUri="local://$ALICE_ROOT";
399 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
400 AliWarning("Default CDB storage not yet set !!!!");
401 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
402 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
405 man->SetDefaultStorage(fCDBUri);
408 // Now activate the detector specific CDB storage locations
409 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
410 TObject* obj = fSpecCDBUri[i];
412 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
413 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
414 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
415 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
420 //_____________________________________________________________________________
421 void AliReconstruction::SetDefaultStorage(const char* uri) {
422 // Store the desired default CDB storage location
423 // Activate it later within the Run() method
429 //_____________________________________________________________________________
430 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
431 // Store a detector-specific CDB storage location
432 // Activate it later within the Run() method
434 AliCDBPath aPath(calibType);
435 if(!aPath.IsValid()){
436 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
437 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
438 if(!strcmp(calibType, fgkDetectorName[iDet])) {
439 aPath.SetPath(Form("%s/*", calibType));
440 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
444 if(!aPath.IsValid()){
445 AliError(Form("Not a valid path or detector: %s", calibType));
450 // // check that calibType refers to a "valid" detector name
451 // Bool_t isDetector = kFALSE;
452 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
453 // TString detName = fgkDetectorName[iDet];
454 // if(aPath.GetLevel0() == detName) {
455 // isDetector = kTRUE;
461 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
465 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
466 if (obj) fSpecCDBUri.Remove(obj);
467 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
471 //_____________________________________________________________________________
472 Bool_t AliReconstruction::SetRunNumberFromData()
474 // The method is called in Run() in order
475 // to set a correct run number.
476 // In case of raw data reconstruction the
477 // run number is taken from the raw data header
479 if (fSetRunNumberFromDataCalled) return kTRUE;
480 fSetRunNumberFromDataCalled = kTRUE;
482 AliCDBManager* man = AliCDBManager::Instance();
484 if(man->GetRun() > 0) {
485 AliWarning("Run number is taken from event header! Ignoring settings in AliCDBManager!");
489 AliError("No run loader is found !");
492 // read run number from gAlice
493 if(fRunLoader->GetAliRun())
494 AliCDBManager::Instance()->SetRun(fRunLoader->GetHeader()->GetRun());
497 if(fRawReader->NextEvent()) {
498 AliCDBManager::Instance()->SetRun(fRawReader->GetRunNumber());
499 fRawReader->RewindEvents();
502 AliError("No raw-data events found !");
507 AliError("Neither gAlice nor RawReader objects are found !");
517 //_____________________________________________________________________________
518 void AliReconstruction::SetCDBLock() {
519 // Set CDB lock: from now on it is forbidden to reset the run number
520 // or the default storage or to activate any further storage!
522 AliCDBManager::Instance()->SetLock(1);
525 //_____________________________________________________________________________
526 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
528 // Read the alignment objects from CDB.
529 // Each detector is supposed to have the
530 // alignment objects in DET/Align/Data CDB path.
531 // All the detector objects are then collected,
532 // sorted by geometry level (starting from ALIC) and
533 // then applied to the TGeo geometry.
534 // Finally an overlaps check is performed.
536 // Load alignment data from CDB and fill fAlignObjArray
537 if(fLoadAlignFromCDB){
539 TString detStr = detectors;
540 TString loadAlObjsListOfDets = "";
542 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
543 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
544 loadAlObjsListOfDets += fgkDetectorName[iDet];
545 loadAlObjsListOfDets += " ";
546 } // end loop over detectors
547 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
548 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
550 // Check if the array with alignment objects was
551 // provided by the user. If yes, apply the objects
552 // to the present TGeo geometry
553 if (fAlignObjArray) {
554 if (gGeoManager && gGeoManager->IsClosed()) {
555 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
556 AliError("The misalignment of one or more volumes failed!"
557 "Compare the list of simulated detectors and the list of detector alignment data!");
562 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
568 delete fAlignObjArray; fAlignObjArray=0;
573 //_____________________________________________________________________________
574 void AliReconstruction::SetGAliceFile(const char* fileName)
576 // set the name of the galice file
578 fGAliceFileName = fileName;
581 //_____________________________________________________________________________
582 void AliReconstruction::SetOption(const char* detector, const char* option)
584 // set options for the reconstruction of a detector
586 TObject* obj = fOptions.FindObject(detector);
587 if (obj) fOptions.Remove(obj);
588 fOptions.Add(new TNamed(detector, option));
592 //_____________________________________________________________________________
593 Bool_t AliReconstruction::Run(const char* input, Bool_t IsOnline)
595 // run the reconstruction
601 if (!input) input = fInput.Data();
602 TString fileName(input);
603 if (fileName.EndsWith("/")) {
604 fRawReader = new AliRawReaderFile(fileName);
605 } else if (fileName.EndsWith(".root")) {
606 fRawReader = new AliRawReaderRoot(fileName);
607 } else if (!fileName.IsNull()) {
608 fRawReader = new AliRawReaderDate(fileName);
613 AliError("Null pointer to the event structure!");
616 fRawReader = new AliRawReaderDate((void *)input);
619 if (!fEquipIdMap.IsNull() && fRawReader)
620 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
622 AliSysInfo::AddStamp("Start");
623 // get the run loader
624 if (!InitRunLoader()) return kFALSE;
625 AliSysInfo::AddStamp("LoadLoader");
627 // Initialize the CDB storage
630 AliSysInfo::AddStamp("LoadCDB");
632 // Set run number in CDBManager (if it is not already set by the user)
633 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
635 // Set CDB lock: from now on it is forbidden to reset the run number
636 // or the default storage or to activate any further storage!
639 // Import ideal TGeo geometry and apply misalignment
641 TString geom(gSystem->DirName(fGAliceFileName));
642 geom += "/geometry.root";
643 AliGeomManager::LoadGeometry(geom.Data());
644 if (!gGeoManager) if (fStopOnError) return kFALSE;
647 if (!MisalignGeometry(fLoadAlignData)) if (fStopOnError) return kFALSE;
648 AliSysInfo::AddStamp("LoadGeom");
651 AliQADataMakerSteer qas ;
652 if (fRunQA && fRawReader) qas.Run(fRunLocalReconstruction, fRawReader) ;
653 // checking the QA of previous steps
657 // local reconstruction
658 if (!fRunLocalReconstruction.IsNull()) {
659 if (!RunLocalReconstruction(fRunLocalReconstruction)) {
660 if (fStopOnError) {CleanUp(); return kFALSE;}
666 if (fRunVertexFinder && !CreateVertexer()) {
672 AliSysInfo::AddStamp("Vertexer");
675 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
681 AliSysInfo::AddStamp("LoadTrackers");
683 // get the possibly already existing ESD file and tree
684 AliESDEvent* esd = new AliESDEvent(); AliESDEvent* hltesd = new AliESDEvent();
685 TFile* fileOld = NULL;
686 TTree* treeOld = NULL; TTree *hlttreeOld = NULL;
687 if (!gSystem->AccessPathName("AliESDs.root")){
688 gSystem->CopyFile("AliESDs.root", "AliESDs.old.root", kTRUE);
689 fileOld = TFile::Open("AliESDs.old.root");
690 if (fileOld && fileOld->IsOpen()) {
691 treeOld = (TTree*) fileOld->Get("esdTree");
692 if (treeOld)esd->ReadFromTree(treeOld);
693 hlttreeOld = (TTree*) fileOld->Get("HLTesdTree");
694 if (hlttreeOld) hltesd->ReadFromTree(hlttreeOld);
698 // create the ESD output file and tree
699 TFile* file = TFile::Open("AliESDs.root", "RECREATE");
700 file->SetCompressionLevel(2);
701 if (!file->IsOpen()) {
702 AliError("opening AliESDs.root failed");
703 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
706 TTree* tree = new TTree("esdTree", "Tree with ESD objects");
707 esd = new AliESDEvent();
708 esd->CreateStdContent();
709 esd->WriteToTree(tree);
711 TTree* hlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
712 hltesd = new AliESDEvent();
713 hltesd->CreateStdContent();
714 hltesd->WriteToTree(hlttree);
717 delete esd; delete hltesd;
718 esd = NULL; hltesd = NULL;
720 // create the branch with ESD additions
724 AliESDfriend *esdf = 0;
725 if (fWriteESDfriend) {
726 esdf = new AliESDfriend();
727 TBranch *br=tree->Branch("ESDfriend.","AliESDfriend", &esdf);
728 br->SetFile("AliESDfriends.root");
729 esd->AddObject(esdf);
733 // Get the GRP CDB entry
734 AliCDBEntry* entryGRP = AliCDBManager::Instance()->Get("GRP/GRP/Data");
737 fGRPList = dynamic_cast<TList*> (entryGRP->GetObject());
739 AliError("No GRP entry found in OCDB!");
742 // Get the diamond profile from OCDB
743 AliCDBEntry* entry = AliCDBManager::Instance()
744 ->Get("GRP/Calib/MeanVertex");
747 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
749 AliError("No diamond profile found in OCDB!");
752 AliVertexerTracks tVertexer(AliTracker::GetBz());
753 if(fDiamondProfile && fMeanVertexConstraint) tVertexer.SetVtxStart(fDiamondProfile);
755 if (fRawReader) fRawReader->RewindEvents();
758 gSystem->GetProcInfo(&ProcInfo);
759 AliInfo(Form("Current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual));
762 //Initialize the QA and start of cycle for out-of-cycle QA
764 TString detStr(fFillESD);
765 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
766 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
767 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
769 AliInfo(Form("Initializing the QA data maker for %s",
770 fgkDetectorName[iDet]));
771 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
772 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
774 qadm->StartOfCycle(AliQA::kRECPOINTS);
775 qadm->StartOfCycle(AliQA::kESDS,"same");
779 AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
780 AliInfo(Form("Initializing the global QA data maker"));
782 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
783 AliTracker::SetResidualsArray(arr);
784 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
786 qadm->StartOfCycle(AliQA::kRECPOINTS);
787 qadm->StartOfCycle(AliQA::kESDS);
792 //Initialize the Plane Efficiency framework
793 if (fRunPlaneEff && !InitPlaneEff()) {
794 if(fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
797 //******* The loop over events
798 for (Int_t iEvent = 0; iEvent < fRunLoader->GetNumberOfEvents(); iEvent++) {
799 if (fRawReader) fRawReader->NextEvent();
800 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
801 // copy old ESD to the new one
803 esd->ReadFromTree(treeOld);
804 treeOld->GetEntry(iEvent);
808 esd->ReadFromTree(hlttreeOld);
809 hlttreeOld->GetEntry(iEvent);
815 AliInfo(Form("processing event %d", iEvent));
817 //Start of cycle for the in-loop QA
818 if (fRunQA && fInLoopQA) {
819 TString detStr(fFillESD);
820 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
821 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
822 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
824 qadm->StartOfCycle(AliQA::kRECPOINTS);
825 qadm->StartOfCycle(AliQA::kESDS, "same") ;
828 AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
829 qadm->StartOfCycle(AliQA::kRECPOINTS);
830 qadm->StartOfCycle(AliQA::kESDS);
834 fRunLoader->GetEvent(iEvent);
837 sprintf(aFileName, "ESD_%d.%d_final.root",
838 fRunLoader->GetHeader()->GetRun(),
839 fRunLoader->GetHeader()->GetEventNrInRun());
840 if (!gSystem->AccessPathName(aFileName)) continue;
842 // local signle event reconstruction
843 if (!fRunLocalReconstruction.IsNull()) {
844 if (!RunLocalEventReconstruction(fRunLocalReconstruction)) {
845 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
849 esd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
850 hltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
851 esd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
852 hltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
854 // Set magnetic field from the tracker
855 esd->SetMagneticField(AliTracker::GetBz());
856 hltesd->SetMagneticField(AliTracker::GetBz());
860 // Fill raw-data error log into the ESD
861 if (fRawReader) FillRawDataErrorLog(iEvent,esd);
864 if (fRunVertexFinder) {
865 if (!ReadESD(esd, "vertex")) {
866 if (!RunVertexFinder(esd)) {
867 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
869 if (fCheckPointLevel > 0) WriteESD(esd, "vertex");
874 if (!fRunTracking.IsNull()) {
875 if (fRunHLTTracking) {
876 hltesd->SetVertex(esd->GetVertex());
877 if (!RunHLTTracking(hltesd)) {
878 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
884 if (!fRunTracking.IsNull()) {
885 if (fRunMuonTracking) {
886 if (!RunMuonTracking(esd)) {
887 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
893 if (!fRunTracking.IsNull()) {
894 if (!ReadESD(esd, "tracking")) {
895 if (!RunTracking(esd)) {
896 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
898 if (fCheckPointLevel > 0) WriteESD(esd, "tracking");
903 if (!fFillESD.IsNull()) {
904 if (!FillESD(esd, fFillESD)) {
905 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
909 // fill Event header information from the RawEventHeader
910 if (fRawReader){FillRawEventHeaderESD(esd);}
913 AliESDpid::MakePID(esd);
914 if (fCheckPointLevel > 1) WriteESD(esd, "PID");
916 if (fFillTriggerESD) {
917 if (!ReadESD(esd, "trigger")) {
918 if (!FillTriggerESD(esd)) {
919 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
921 if (fCheckPointLevel > 1) WriteESD(esd, "trigger");
927 //Try to improve the reconstructed primary vertex position using the tracks
928 AliESDVertex *pvtx=0;
929 Bool_t dovertex=kTRUE;
930 TObject* obj = fOptions.FindObject("ITS");
932 TString optITS = obj->GetTitle();
933 if (optITS.Contains("cosmics") || optITS.Contains("COSMICS"))
936 if(dovertex) pvtx=tVertexer.FindPrimaryVertex(esd);
937 if(fDiamondProfile) esd->SetDiamond(fDiamondProfile);
940 if (pvtx->GetStatus()) {
941 // Store the improved primary vertex
942 esd->SetPrimaryVertex(pvtx);
943 // Propagate the tracks to the DCA to the improved primary vertex
944 Double_t somethingbig = 777.;
945 Double_t bz = esd->GetMagneticField();
946 Int_t nt=esd->GetNumberOfTracks();
948 AliESDtrack *t = esd->GetTrack(nt);
949 t->RelateToVertex(pvtx, bz, somethingbig);
956 vtxer.Tracks2V0vertices(esd);
958 if (fRunCascadeFinder) {
960 AliCascadeVertexer cvtxer;
961 cvtxer.V0sTracks2CascadeVertices(esd);
966 if (fCleanESD) CleanESD(esd);
967 if (fWriteESDfriend) {
968 esdf->~AliESDfriend();
969 new (esdf) AliESDfriend(); // Reset...
970 esd->GetESDfriend(esdf);
977 if (fCheckPointLevel > 0) WriteESD(esd, "final");
980 if (fWriteESDfriend) {
981 esdf->~AliESDfriend();
982 new (esdf) AliESDfriend(); // Reset...
985 gSystem->GetProcInfo(&ProcInfo);
986 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
989 // End of cycle for the in-loop QA
990 if (fRunQA && fInLoopQA) {
991 RunQA(fFillESD.Data(), esd);
992 TString detStr(fFillESD);
993 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
994 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
995 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
997 qadm->EndOfCycle(AliQA::kRECPOINTS);
998 qadm->EndOfCycle(AliQA::kESDS);
1002 AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
1004 qadm->EndOfCycle(AliQA::kRECPOINTS);
1005 qadm->EndOfCycle(AliQA::kESDS);
1011 //******** End of the loop over events
1015 tree->GetUserInfo()->Add(esd);
1016 hlttree->GetUserInfo()->Add(hltesd);
1018 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1019 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1021 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1022 cdbMapCopy->SetOwner(1);
1023 cdbMapCopy->SetName("cdbMap");
1024 TIter iter(cdbMap->GetTable());
1027 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1028 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1029 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1030 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1033 TList *cdbListCopy = new TList();
1034 cdbListCopy->SetOwner(1);
1035 cdbListCopy->SetName("cdbList");
1037 TIter iter2(cdbList);
1040 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1041 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1044 tree->GetUserInfo()->Add(cdbMapCopy);
1045 tree->GetUserInfo()->Add(cdbListCopy);
1048 if(fESDPar.Contains("ESD.par")){
1049 AliInfo("Attaching ESD.par to Tree");
1050 TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1051 tree->GetUserInfo()->Add(fn);
1057 if (fWriteESDfriend)
1058 tree->SetBranchStatus("ESDfriend*",0);
1059 // we want to have only one tree version number
1060 tree->Write(tree->GetName(),TObject::kOverwrite);
1063 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1064 if (fRunPlaneEff && !FinishPlaneEff()) {
1065 AliWarning("Finish PlaneEff evaluation failed");
1069 CleanUp(file, fileOld);
1072 TFile *esdFile = TFile::Open("AliESDs.root", "READONLY");
1073 TFile *aodFile = TFile::Open("AliAOD.root", "RECREATE");
1074 ESDFile2AODFile(esdFile, aodFile);
1079 // Create tags for the events in the ESD tree (the ESD tree is always present)
1080 // In case of empty events the tags will contain dummy values
1081 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1082 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPList);
1084 AliAODTagCreator *aodtagCreator = new AliAODTagCreator();
1085 aodtagCreator->CreateAODTags(fFirstEvent,fLastEvent,fGRPList);
1088 //Finish QA and end of cycle for out-of-loop QA
1089 if (fRunQA && !fInLoopQA) {
1090 qas.Run(fRunLocalReconstruction.Data(), AliQA::kRECPOINTS);
1092 qas.Run(fRunTracking.Data(), AliQA::kESDS);
1095 AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
1097 qadm->EndOfCycle(AliQA::kRECPOINTS);
1098 qadm->EndOfCycle(AliQA::kESDS);
1104 // Cleanup of CDB manager: cache and active storages!
1105 AliCDBManager::Instance()->ClearCache();
1112 //_____________________________________________________________________________
1113 Bool_t AliReconstruction::RunLocalReconstruction(const TString& /*detectors*/)
1115 // run the local reconstruction
1116 static Int_t eventNr=0;
1117 AliCodeTimerAuto("")
1119 // AliCDBManager* man = AliCDBManager::Instance();
1120 // Bool_t origCache = man->GetCacheFlag();
1122 // TString detStr = detectors;
1123 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1124 // if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1125 // AliReconstructor* reconstructor = GetReconstructor(iDet);
1126 // if (!reconstructor) continue;
1127 // if (reconstructor->HasLocalReconstruction()) continue;
1129 // AliCodeTimerStart(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1130 // AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1132 // AliCodeTimerStart(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1133 // AliInfo(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1135 // man->SetCacheFlag(kTRUE);
1136 // TString calibPath = Form("%s/Calib/*", fgkDetectorName[iDet]);
1137 // man->GetAll(calibPath); // entries are cached!
1139 // AliCodeTimerStop(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1141 // if (fRawReader) {
1142 // fRawReader->RewindEvents();
1143 // reconstructor->Reconstruct(fRunLoader, fRawReader);
1145 // reconstructor->Reconstruct(fRunLoader);
1148 // AliCodeTimerStop(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1149 // AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr));
1151 // // unload calibration data
1152 // man->UnloadFromCache(calibPath);
1153 // //man->ClearCache();
1156 // man->SetCacheFlag(origCache);
1158 // if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1159 // AliError(Form("the following detectors were not found: %s",
1161 // if (fStopOnError) return kFALSE;
1168 //_____________________________________________________________________________
1169 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1171 // run the local reconstruction
1173 static Int_t eventNr=0;
1174 AliCodeTimerAuto("")
1176 TString detStr = detectors;
1177 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1178 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1179 AliReconstructor* reconstructor = GetReconstructor(iDet);
1180 if (!reconstructor) continue;
1181 AliLoader* loader = fLoader[iDet];
1183 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1186 // conversion of digits
1187 if (fRawReader && reconstructor->HasDigitConversion()) {
1188 AliInfo(Form("converting raw data digits into root objects for %s",
1189 fgkDetectorName[iDet]));
1190 AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1191 fgkDetectorName[iDet]));
1192 loader->LoadDigits("update");
1193 loader->CleanDigits();
1194 loader->MakeDigitsContainer();
1195 TTree* digitsTree = loader->TreeD();
1196 reconstructor->ConvertDigits(fRawReader, digitsTree);
1197 loader->WriteDigits("OVERWRITE");
1198 loader->UnloadDigits();
1200 // local reconstruction
1201 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1202 AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1203 loader->LoadRecPoints("update");
1204 loader->CleanRecPoints();
1205 loader->MakeRecPointsContainer();
1206 TTree* clustersTree = loader->TreeR();
1207 if (fRawReader && !reconstructor->HasDigitConversion()) {
1208 reconstructor->Reconstruct(fRawReader, clustersTree);
1210 loader->LoadDigits("read");
1211 TTree* digitsTree = loader->TreeD();
1213 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1214 if (fStopOnError) return kFALSE;
1216 reconstructor->Reconstruct(digitsTree, clustersTree);
1218 loader->UnloadDigits();
1221 // In-loop QA for local reconstrucion
1222 if (fRunQA && fInLoopQA) {
1223 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1226 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1228 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1230 qadm->Exec(AliQA::kRECPOINTS, clustersTree) ;
1233 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1237 loader->WriteRecPoints("OVERWRITE");
1238 loader->UnloadRecPoints();
1239 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1242 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1243 AliError(Form("the following detectors were not found: %s",
1245 if (fStopOnError) return kFALSE;
1251 //_____________________________________________________________________________
1252 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1254 // run the barrel tracking
1256 AliCodeTimerAuto("")
1258 AliESDVertex* vertex = NULL;
1259 Double_t vtxPos[3] = {0, 0, 0};
1260 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1261 TArrayF mcVertex(3);
1262 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1263 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1264 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1268 if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1269 AliInfo("running the ITS vertex finder");
1270 if (fLoader[0]) fLoader[0]->LoadRecPoints();
1271 vertex = fVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber());
1272 if (fLoader[0]) fLoader[0]->UnloadRecPoints();
1274 AliWarning("Vertex not found");
1275 vertex = new AliESDVertex();
1276 vertex->SetName("default");
1279 vertex->SetName("reconstructed");
1283 AliInfo("getting the primary vertex from MC");
1284 vertex = new AliESDVertex(vtxPos, vtxErr);
1288 vertex->GetXYZ(vtxPos);
1289 vertex->GetSigmaXYZ(vtxErr);
1291 AliWarning("no vertex reconstructed");
1292 vertex = new AliESDVertex(vtxPos, vtxErr);
1294 esd->SetVertex(vertex);
1295 // if SPD multiplicity has been determined, it is stored in the ESD
1296 AliMultiplicity *mult = fVertexer->GetMultiplicity();
1297 if(mult)esd->SetMultiplicity(mult);
1299 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1300 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1307 //_____________________________________________________________________________
1308 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
1310 // run the HLT barrel tracking
1312 AliCodeTimerAuto("")
1315 AliError("Missing runLoader!");
1319 AliInfo("running HLT tracking");
1321 // Get a pointer to the HLT reconstructor
1322 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
1323 if (!reconstructor) return kFALSE;
1326 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1327 TString detName = fgkDetectorName[iDet];
1328 AliDebug(1, Form("%s HLT tracking", detName.Data()));
1329 reconstructor->SetOption(detName.Data());
1330 AliTracker *tracker = reconstructor->CreateTracker();
1332 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
1333 if (fStopOnError) return kFALSE;
1337 Double_t vtxErr[3]={0.005,0.005,0.010};
1338 const AliESDVertex *vertex = esd->GetVertex();
1339 vertex->GetXYZ(vtxPos);
1340 tracker->SetVertex(vtxPos,vtxErr);
1342 fLoader[iDet]->LoadRecPoints("read");
1343 TTree* tree = fLoader[iDet]->TreeR();
1345 AliError(Form("Can't get the %s cluster tree", detName.Data()));
1348 tracker->LoadClusters(tree);
1350 if (tracker->Clusters2Tracks(esd) != 0) {
1351 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
1355 tracker->UnloadClusters();
1363 //_____________________________________________________________________________
1364 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
1366 // run the muon spectrometer tracking
1368 AliCodeTimerAuto("")
1371 AliError("Missing runLoader!");
1374 Int_t iDet = 7; // for MUON
1376 AliInfo("is running...");
1378 // Get a pointer to the MUON reconstructor
1379 AliReconstructor *reconstructor = GetReconstructor(iDet);
1380 if (!reconstructor) return kFALSE;
1383 TString detName = fgkDetectorName[iDet];
1384 AliDebug(1, Form("%s tracking", detName.Data()));
1385 AliTracker *tracker = reconstructor->CreateTracker();
1387 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1392 fLoader[iDet]->LoadRecPoints("read");
1394 tracker->LoadClusters(fLoader[iDet]->TreeR());
1396 Int_t rv = tracker->Clusters2Tracks(esd);
1400 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1404 fLoader[iDet]->UnloadRecPoints();
1406 tracker->UnloadClusters();
1414 //_____________________________________________________________________________
1415 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
1417 // run the barrel tracking
1418 static Int_t eventNr=0;
1419 AliCodeTimerAuto("")
1421 AliInfo("running tracking");
1423 //Fill the ESD with the T0 info (will be used by the TOF)
1424 if (fReconstructor[11] && fLoader[11]) {
1425 fLoader[11]->LoadRecPoints("READ");
1426 TTree *treeR = fLoader[11]->TreeR();
1427 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
1430 // pass 1: TPC + ITS inwards
1431 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1432 if (!fTracker[iDet]) continue;
1433 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
1436 fLoader[iDet]->LoadRecPoints("read");
1437 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
1438 TTree* tree = fLoader[iDet]->TreeR();
1440 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1443 fTracker[iDet]->LoadClusters(tree);
1444 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1446 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
1447 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1450 if (fCheckPointLevel > 1) {
1451 WriteESD(esd, Form("%s.tracking", fgkDetectorName[iDet]));
1453 // preliminary PID in TPC needed by the ITS tracker
1455 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1456 AliESDpid::MakePID(esd);
1458 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
1461 // pass 2: ALL backwards
1463 if (fRunQA && fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1465 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1466 if (!fTracker[iDet]) continue;
1467 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
1470 if (iDet > 1) { // all except ITS, TPC
1472 fLoader[iDet]->LoadRecPoints("read");
1473 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
1474 tree = fLoader[iDet]->TreeR();
1476 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1479 fTracker[iDet]->LoadClusters(tree);
1480 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1484 if (fTracker[iDet]->PropagateBack(esd) != 0) {
1485 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
1488 if (fCheckPointLevel > 1) {
1489 WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
1493 if (iDet > 2) { // all except ITS, TPC, TRD
1494 fTracker[iDet]->UnloadClusters();
1495 fLoader[iDet]->UnloadRecPoints();
1497 // updated PID in TPC needed by the ITS tracker -MI
1499 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1500 AliESDpid::MakePID(esd);
1502 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1505 if (fRunQA && fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1507 // write space-points to the ESD in case alignment data output
1509 if (fWriteAlignmentData)
1510 WriteAlignmentData(esd);
1512 // pass 3: TRD + TPC + ITS refit inwards
1515 for (Int_t iDet = 2; iDet >= 0; iDet--) {
1516 if (!fTracker[iDet]) continue;
1517 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
1520 if (fTracker[iDet]->RefitInward(esd) != 0) {
1521 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
1524 // run postprocessing
1525 if (fTracker[iDet]->PostProcess(esd) != 0) {
1526 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
1529 if (fCheckPointLevel > 1) {
1530 WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
1532 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1534 fTracker[iDet]->UnloadClusters();
1535 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
1536 fLoader[iDet]->UnloadRecPoints();
1537 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
1541 // Propagate track to the vertex - if not done by ITS
1543 Int_t ntracks = esd->GetNumberOfTracks();
1544 for (Int_t itrack=0; itrack<ntracks; itrack++){
1545 const Double_t kRadius = 3; // beam pipe radius
1546 const Double_t kMaxStep = 5; // max step
1547 const Double_t kMaxD = 123456; // max distance to prim vertex
1548 Double_t fieldZ = AliTracker::GetBz(); //
1549 AliESDtrack * track = esd->GetTrack(itrack);
1550 if (!track) continue;
1551 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1552 AliTracker::PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kTRUE);
1553 track->RelateToVertex(esd->GetVertex(),fieldZ, kMaxD);
1559 //_____________________________________________________________________________
1560 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
1562 // Remove the data which are not needed for the physics analysis.
1565 Int_t nTracks=esd->GetNumberOfTracks();
1566 Int_t nV0s=esd->GetNumberOfV0s();
1568 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
1570 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
1571 Bool_t rc=esd->Clean(cleanPars);
1573 nTracks=esd->GetNumberOfTracks();
1574 nV0s=esd->GetNumberOfV0s();
1576 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
1581 //_____________________________________________________________________________
1582 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
1584 // fill the event summary data
1586 AliCodeTimerAuto("")
1587 static Int_t eventNr=0;
1588 TString detStr = detectors;
1590 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1591 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1592 AliReconstructor* reconstructor = GetReconstructor(iDet);
1593 if (!reconstructor) continue;
1594 if (!ReadESD(esd, fgkDetectorName[iDet])) {
1595 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
1596 TTree* clustersTree = NULL;
1597 if (fLoader[iDet]) {
1598 fLoader[iDet]->LoadRecPoints("read");
1599 clustersTree = fLoader[iDet]->TreeR();
1600 if (!clustersTree) {
1601 AliError(Form("Can't get the %s clusters tree",
1602 fgkDetectorName[iDet]));
1603 if (fStopOnError) return kFALSE;
1606 if (fRawReader && !reconstructor->HasDigitConversion()) {
1607 reconstructor->FillESD(fRawReader, clustersTree, esd);
1609 TTree* digitsTree = NULL;
1610 if (fLoader[iDet]) {
1611 fLoader[iDet]->LoadDigits("read");
1612 digitsTree = fLoader[iDet]->TreeD();
1614 AliError(Form("Can't get the %s digits tree",
1615 fgkDetectorName[iDet]));
1616 if (fStopOnError) return kFALSE;
1619 reconstructor->FillESD(digitsTree, clustersTree, esd);
1620 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
1622 if (fLoader[iDet]) {
1623 fLoader[iDet]->UnloadRecPoints();
1626 if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
1630 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1631 AliError(Form("the following detectors were not found: %s",
1633 if (fStopOnError) return kFALSE;
1635 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
1640 //_____________________________________________________________________________
1641 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
1643 // Reads the trigger decision which is
1644 // stored in Trigger.root file and fills
1645 // the corresponding esd entries
1647 AliCodeTimerAuto("")
1649 AliInfo("Filling trigger information into the ESD");
1651 AliCentralTrigger *aCTP = NULL;
1654 AliCTPRawStream input(fRawReader);
1655 if (!input.Next()) {
1656 AliError("No valid CTP (trigger) DDL raw data is found ! The trigger information is not stored in the ESD !");
1659 esd->SetTriggerMask(input.GetClassMask());
1660 esd->SetTriggerCluster(input.GetClusterMask());
1662 aCTP = new AliCentralTrigger();
1663 TString configstr("");
1664 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
1665 AliError("No trigger configuration found in OCDB! The trigger classes information will no be stored in ESD!");
1670 AliRunLoader *runloader = AliRunLoader::GetRunLoader();
1672 if (!runloader->LoadTrigger()) {
1673 AliCentralTrigger *aCTP = runloader->GetTrigger();
1674 esd->SetTriggerMask(aCTP->GetClassMask());
1675 esd->SetTriggerCluster(aCTP->GetClusterMask());
1678 AliWarning("No trigger can be loaded! The trigger information is not stored in the ESD !");
1683 AliError("No run loader is available! The trigger information is not stored in the ESD !");
1688 // Now fill the trigger class names into AliESDRun object
1689 AliTriggerConfiguration *config = aCTP->GetConfiguration();
1691 AliError("No trigger configuration has been found! The trigger classes information will no be stored in ESD!");
1695 const TObjArray& classesArray = config->GetClasses();
1696 Int_t nclasses = classesArray.GetEntriesFast();
1697 for( Int_t j=0; j<nclasses; j++ ) {
1698 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At( j );
1699 Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
1700 esd->SetTriggerClass(trclass->GetName(),trindex);
1710 //_____________________________________________________________________________
1711 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
1714 // Filling information from RawReader Header
1717 AliInfo("Filling information from RawReader Header");
1718 esd->SetBunchCrossNumber(0);
1719 esd->SetOrbitNumber(0);
1720 esd->SetPeriodNumber(0);
1721 esd->SetTimeStamp(0);
1722 esd->SetEventType(0);
1723 const AliRawEventHeaderBase * eventHeader = fRawReader->GetEventHeader();
1726 const UInt_t *id = eventHeader->GetP("Id");
1727 esd->SetBunchCrossNumber((id)[1]&0x00000fff);
1728 esd->SetOrbitNumber((((id)[0]<<20)&0xf00000)|(((id)[1]>>12)&0xfffff));
1729 esd->SetPeriodNumber(((id)[0]>>4)&0x0fffffff);
1731 esd->SetTimeStamp((eventHeader->Get("Timestamp")));
1732 esd->SetEventType((eventHeader->Get("Type")));
1739 //_____________________________________________________________________________
1740 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
1742 // check whether detName is contained in detectors
1743 // if yes, it is removed from detectors
1745 // check if all detectors are selected
1746 if ((detectors.CompareTo("ALL") == 0) ||
1747 detectors.BeginsWith("ALL ") ||
1748 detectors.EndsWith(" ALL") ||
1749 detectors.Contains(" ALL ")) {
1754 // search for the given detector
1755 Bool_t result = kFALSE;
1756 if ((detectors.CompareTo(detName) == 0) ||
1757 detectors.BeginsWith(detName+" ") ||
1758 detectors.EndsWith(" "+detName) ||
1759 detectors.Contains(" "+detName+" ")) {
1760 detectors.ReplaceAll(detName, "");
1764 // clean up the detectors string
1765 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1766 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1767 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1772 //_____________________________________________________________________________
1773 Bool_t AliReconstruction::InitRunLoader()
1775 // get or create the run loader
1777 if (gAlice) delete gAlice;
1780 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
1781 // load all base libraries to get the loader classes
1782 TString libs = gSystem->GetLibraries();
1783 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1784 TString detName = fgkDetectorName[iDet];
1785 if (detName == "HLT") continue;
1786 if (libs.Contains("lib" + detName + "base.so")) continue;
1787 gSystem->Load("lib" + detName + "base.so");
1789 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
1791 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1795 fRunLoader->CdGAFile();
1796 if (gFile->GetKey(AliRunLoader::GetGAliceName())) {
1797 if (fRunLoader->LoadgAlice() == 0) {
1798 gAlice = fRunLoader->GetAliRun();
1799 AliTracker::SetFieldMap(gAlice->Field(),fUniformField);
1802 if (!gAlice && !fRawReader) {
1803 AliError(Form("no gAlice object found in file %s",
1804 fGAliceFileName.Data()));
1809 //PH This is a temporary fix to give access to the kinematics
1810 //PH that is needed for the labels of ITS clusters
1811 fRunLoader->LoadHeader();
1812 fRunLoader->LoadKinematics();
1814 } else { // galice.root does not exist
1816 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
1820 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
1821 AliConfig::GetDefaultEventFolderName(),
1824 AliError(Form("could not create run loader in file %s",
1825 fGAliceFileName.Data()));
1829 fRunLoader->MakeTree("E");
1831 while (fRawReader->NextEvent()) {
1832 fRunLoader->SetEventNumber(iEvent);
1833 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1835 fRunLoader->MakeTree("H");
1836 fRunLoader->TreeE()->Fill();
1839 fRawReader->RewindEvents();
1840 if (fNumberOfEventsPerFile > 0)
1841 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
1843 fRunLoader->SetNumberOfEventsPerFile(iEvent);
1844 fRunLoader->WriteHeader("OVERWRITE");
1845 fRunLoader->CdGAFile();
1846 fRunLoader->Write(0, TObject::kOverwrite);
1847 // AliTracker::SetFieldMap(???);
1853 //_____________________________________________________________________________
1854 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
1856 // get the reconstructor object and the loader for a detector
1858 if (fReconstructor[iDet]) return fReconstructor[iDet];
1860 // load the reconstructor object
1861 TPluginManager* pluginManager = gROOT->GetPluginManager();
1862 TString detName = fgkDetectorName[iDet];
1863 TString recName = "Ali" + detName + "Reconstructor";
1864 if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT")) return NULL;
1866 AliReconstructor* reconstructor = NULL;
1867 // first check if a plugin is defined for the reconstructor
1868 TPluginHandler* pluginHandler =
1869 pluginManager->FindHandler("AliReconstructor", detName);
1870 // if not, add a plugin for it
1871 if (!pluginHandler) {
1872 AliDebug(1, Form("defining plugin for %s", recName.Data()));
1873 TString libs = gSystem->GetLibraries();
1874 if (libs.Contains("lib" + detName + "base.so") ||
1875 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
1876 pluginManager->AddHandler("AliReconstructor", detName,
1877 recName, detName + "rec", recName + "()");
1879 pluginManager->AddHandler("AliReconstructor", detName,
1880 recName, detName, recName + "()");
1882 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
1884 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
1885 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
1887 if (reconstructor) {
1888 TObject* obj = fOptions.FindObject(detName.Data());
1889 if (obj) reconstructor->SetOption(obj->GetTitle());
1890 reconstructor->Init();
1891 fReconstructor[iDet] = reconstructor;
1894 // get or create the loader
1895 if (detName != "HLT") {
1896 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
1897 if (!fLoader[iDet]) {
1898 AliConfig::Instance()
1899 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
1901 // first check if a plugin is defined for the loader
1903 pluginManager->FindHandler("AliLoader", detName);
1904 // if not, add a plugin for it
1905 if (!pluginHandler) {
1906 TString loaderName = "Ali" + detName + "Loader";
1907 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
1908 pluginManager->AddHandler("AliLoader", detName,
1909 loaderName, detName + "base",
1910 loaderName + "(const char*, TFolder*)");
1911 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
1913 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
1915 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
1916 fRunLoader->GetEventFolder());
1918 if (!fLoader[iDet]) { // use default loader
1919 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
1921 if (!fLoader[iDet]) {
1922 AliWarning(Form("couldn't get loader for %s", detName.Data()));
1923 if (fStopOnError) return NULL;
1925 fRunLoader->AddLoader(fLoader[iDet]);
1926 fRunLoader->CdGAFile();
1927 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
1928 fRunLoader->Write(0, TObject::kOverwrite);
1933 return reconstructor;
1936 //_____________________________________________________________________________
1937 Bool_t AliReconstruction::CreateVertexer()
1939 // create the vertexer
1942 AliReconstructor* itsReconstructor = GetReconstructor(0);
1943 if (itsReconstructor) {
1944 fVertexer = itsReconstructor->CreateVertexer();
1947 AliWarning("couldn't create a vertexer for ITS");
1948 if (fStopOnError) return kFALSE;
1954 //_____________________________________________________________________________
1955 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
1957 // create the trackers
1959 TString detStr = detectors;
1960 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1961 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1962 AliReconstructor* reconstructor = GetReconstructor(iDet);
1963 if (!reconstructor) continue;
1964 TString detName = fgkDetectorName[iDet];
1965 if (detName == "HLT") {
1966 fRunHLTTracking = kTRUE;
1969 if (detName == "MUON") {
1970 fRunMuonTracking = kTRUE;
1975 fTracker[iDet] = reconstructor->CreateTracker();
1976 if (!fTracker[iDet] && (iDet < 7)) {
1977 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1978 if (fStopOnError) return kFALSE;
1980 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
1986 //_____________________________________________________________________________
1987 void AliReconstruction::CleanUp(TFile* file, TFile* fileOld)
1989 // delete trackers and the run loader and close and delete the file
1991 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1992 delete fReconstructor[iDet];
1993 fReconstructor[iDet] = NULL;
1994 fLoader[iDet] = NULL;
1995 delete fTracker[iDet];
1996 fTracker[iDet] = NULL;
1997 // delete fQADataMaker[iDet];
1998 // fQADataMaker[iDet] = NULL;
2003 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2004 delete fDiamondProfile;
2005 fDiamondProfile = NULL;
2024 gSystem->Unlink("AliESDs.old.root");
2028 //_____________________________________________________________________________
2030 Bool_t AliReconstruction::ReadESD(AliESDEvent*& esd, const char* recStep) const
2032 // read the ESD event from a file
2034 if (!esd) return kFALSE;
2036 sprintf(fileName, "ESD_%d.%d_%s.root",
2037 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2038 if (gSystem->AccessPathName(fileName)) return kFALSE;
2040 AliInfo(Form("reading ESD from file %s", fileName));
2041 AliDebug(1, Form("reading ESD from file %s", fileName));
2042 TFile* file = TFile::Open(fileName);
2043 if (!file || !file->IsOpen()) {
2044 AliError(Form("opening %s failed", fileName));
2051 esd = (AliESDEvent*) file->Get("ESD");
2060 //_____________________________________________________________________________
2061 void AliReconstruction::WriteESD(AliESDEvent* esd, const char* recStep) const
2063 // write the ESD event to a file
2067 sprintf(fileName, "ESD_%d.%d_%s.root",
2068 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2070 AliDebug(1, Form("writing ESD to file %s", fileName));
2071 TFile* file = TFile::Open(fileName, "recreate");
2072 if (!file || !file->IsOpen()) {
2073 AliError(Form("opening %s failed", fileName));
2085 //_____________________________________________________________________________
2086 void AliReconstruction::ESDFile2AODFile(TFile* esdFile, TFile* aodFile)
2088 // write all files from the given esd file to an aod file
2090 // create an AliAOD object
2091 AliAODEvent *aod = new AliAODEvent();
2092 aod->CreateStdContent();
2098 TTree *aodTree = new TTree("aodTree", "AliAOD tree");
2099 aodTree->Branch(aod->GetList());
2102 TTree *t = (TTree*) esdFile->Get("esdTree");
2103 AliESDEvent *esd = new AliESDEvent();
2104 esd->ReadFromTree(t);
2106 Int_t nEvents = t->GetEntries();
2108 // set arrays and pointers
2118 // loop over events and fill them
2119 for (Int_t iEvent = 0; iEvent < nEvents; ++iEvent) {
2120 //cout << "event: " << iEvent << endl;
2121 t->GetEntry(iEvent);
2123 // Multiplicity information needed by the header (to be revised!)
2124 Int_t nTracks = esd->GetNumberOfTracks();
2125 Int_t nPosTracks = 0;
2126 for (Int_t iTrack=0; iTrack<nTracks; ++iTrack)
2127 if (esd->GetTrack(iTrack)->Charge()> 0) nPosTracks++;
2129 // Access the header
2130 AliAODHeader *header = aod->GetHeader();
2133 header->SetRunNumber (esd->GetRunNumber() );
2134 header->SetBunchCrossNumber(esd->GetBunchCrossNumber());
2135 header->SetOrbitNumber (esd->GetOrbitNumber() );
2136 header->SetPeriodNumber (esd->GetPeriodNumber() );
2137 header->SetTriggerMask (esd->GetTriggerMask() );
2138 header->SetTriggerCluster (esd->GetTriggerCluster() );
2139 header->SetEventType (esd->GetEventType() );
2140 header->SetMagneticField (esd->GetMagneticField() );
2141 header->SetZDCN1Energy (esd->GetZDCN1Energy() );
2142 header->SetZDCP1Energy (esd->GetZDCP1Energy() );
2143 header->SetZDCN2Energy (esd->GetZDCN2Energy() );
2144 header->SetZDCP2Energy (esd->GetZDCP2Energy() );
2145 header->SetZDCEMEnergy (esd->GetZDCEMEnergy(0),esd->GetZDCEMEnergy(1));
2146 header->SetRefMultiplicity (nTracks);
2147 header->SetRefMultiplicityPos(nPosTracks);
2148 header->SetRefMultiplicityNeg(nTracks - nPosTracks);
2149 header->SetMuonMagFieldScale(-999.); // FIXME
2150 header->SetCentrality(-999.); // FIXME
2152 Int_t nV0s = esd->GetNumberOfV0s();
2153 Int_t nCascades = esd->GetNumberOfCascades();
2154 Int_t nKinks = esd->GetNumberOfKinks();
2155 Int_t nVertices = nV0s + 2*nCascades /*could lead to two vertices, one V0 and the Xi */+ nKinks + 1 /* = prim. vtx*/;
2157 Int_t nCaloClus = esd->GetNumberOfCaloClusters();
2159 Int_t nPmdClus = esd->GetNumberOfPmdTracks();
2161 aod->ResetStd(nTracks, nVertices, nV0s+nCascades, nJets, nCaloClus, nFmdClus, nPmdClus);
2163 // Array to take into account the tracks already added to the AOD
2164 Bool_t * usedTrack = NULL;
2166 usedTrack = new Bool_t[nTracks];
2167 for (Int_t iTrack=0; iTrack<nTracks; ++iTrack) usedTrack[iTrack]=kFALSE;
2169 // Array to take into account the V0s already added to the AOD
2170 Bool_t * usedV0 = NULL;
2172 usedV0 = new Bool_t[nV0s];
2173 for (Int_t iV0=0; iV0<nV0s; ++iV0) usedV0[iV0]=kFALSE;
2175 // Array to take into account the kinks already added to the AOD
2176 Bool_t * usedKink = NULL;
2178 usedKink = new Bool_t[nKinks];
2179 for (Int_t iKink=0; iKink<nKinks; ++iKink) usedKink[iKink]=kFALSE;
2182 // Access to the AOD container of vertices
2183 TClonesArray &vertices = *(aod->GetVertices());
2186 // Access to the AOD container of tracks
2187 TClonesArray &tracks = *(aod->GetTracks());
2190 // Access to the AOD container of V0s
2191 TClonesArray &V0s = *(aod->GetV0s());
2194 // Add primary vertex. The primary tracks will be defined
2195 // after the loops on the composite objects (V0, cascades, kinks)
2196 const AliESDVertex *vtx = esd->GetPrimaryVertex();
2198 vtx->GetXYZ(pos); // position
2199 vtx->GetCovMatrix(covVtx); //covariance matrix
2201 AliAODVertex * primary = new(vertices[jVertices++])
2202 AliAODVertex(pos, covVtx, vtx->GetChi2toNDF(), NULL, -1, AliAODVertex::kPrimary);
2205 AliAODTrack *aodTrack = 0x0;
2207 // Create vertices starting from the most complex objects
2210 for (Int_t nCascade = 0; nCascade < nCascades; ++nCascade) {
2211 AliESDcascade *cascade = esd->GetCascade(nCascade);
2213 cascade->GetXYZ(pos[0], pos[1], pos[2]);
2214 cascade->GetPosCovXi(covVtx);
2216 // Add the cascade vertex
2217 AliAODVertex * vcascade = new(vertices[jVertices++]) AliAODVertex(pos,
2219 cascade->GetChi2Xi(), // = chi2/NDF since NDF = 2*2-3
2222 AliAODVertex::kCascade);
2224 primary->AddDaughter(vcascade); // the cascade 'particle' (represented by a vertex) is added as a daughter to the primary vertex
2226 // Add the V0 from the cascade. The ESD class have to be optimized...
2227 // Now we have to search for the corresponding V0 in the list of V0s
2228 // using the indeces of the positive and negative tracks
2230 Int_t posFromV0 = cascade->GetPindex();
2231 Int_t negFromV0 = cascade->GetNindex();
2234 AliESDv0 * v0 = 0x0;
2237 for (Int_t iV0=0; iV0<nV0s; ++iV0) {
2239 v0 = esd->GetV0(iV0);
2240 Int_t posV0 = v0->GetPindex();
2241 Int_t negV0 = v0->GetNindex();
2243 if (posV0==posFromV0 && negV0==negFromV0) {
2249 AliAODVertex * vV0FromCascade = 0x0;
2251 if (indV0>-1 && !usedV0[indV0]) {
2253 // the V0 exists in the array of V0s and is not used
2255 usedV0[indV0] = kTRUE;
2257 v0->GetXYZ(pos[0], pos[1], pos[2]);
2258 v0->GetPosCov(covVtx);
2260 vV0FromCascade = new(vertices[jVertices++]) AliAODVertex(pos,
2262 v0->GetChi2V0(), // = chi2/NDF since NDF = 2*2-3
2268 // the V0 doesn't exist in the array of V0s or was used
2269 cerr << "Error: event " << iEvent << " cascade " << nCascade
2270 << " The V0 " << indV0
2271 << " doesn't exist in the array of V0s or was used!" << endl;
2273 cascade->GetXYZ(pos[0], pos[1], pos[2]);
2274 cascade->GetPosCov(covVtx);
2276 vV0FromCascade = new(vertices[jVertices++]) AliAODVertex(pos,
2278 v0->GetChi2V0(), // = chi2/NDF since NDF = 2*2-3
2282 vcascade->AddDaughter(vV0FromCascade);
2286 // Add the positive tracks from the V0
2288 if (! usedTrack[posFromV0]) {
2290 usedTrack[posFromV0] = kTRUE;
2292 AliESDtrack *esdTrack = esd->GetTrack(posFromV0);
2293 esdTrack->GetPxPyPz(p_pos);
2294 esdTrack->GetXYZ(pos);
2295 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2296 esdTrack->GetESDpid(pid);
2298 vV0FromCascade->AddDaughter(aodTrack =
2299 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2300 esdTrack->GetLabel(),
2306 (Short_t)esdTrack->Charge(),
2307 esdTrack->GetITSClusterMap(),
2310 kTRUE, // check if this is right
2311 kFALSE, // check if this is right
2312 AliAODTrack::kSecondary)
2314 aodTrack->ConvertAliPIDtoAODPID();
2317 cerr << "Error: event " << iEvent << " cascade " << nCascade
2318 << " track " << posFromV0 << " has already been used!" << endl;
2321 // Add the negative tracks from the V0
2323 if (!usedTrack[negFromV0]) {
2325 usedTrack[negFromV0] = kTRUE;
2327 AliESDtrack *esdTrack = esd->GetTrack(negFromV0);
2328 esdTrack->GetPxPyPz(p_neg);
2329 esdTrack->GetXYZ(pos);
2330 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2331 esdTrack->GetESDpid(pid);
2333 vV0FromCascade->AddDaughter(aodTrack =
2334 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2335 esdTrack->GetLabel(),
2341 (Short_t)esdTrack->Charge(),
2342 esdTrack->GetITSClusterMap(),
2345 kTRUE, // check if this is right
2346 kFALSE, // check if this is right
2347 AliAODTrack::kSecondary)
2349 aodTrack->ConvertAliPIDtoAODPID();
2352 cerr << "Error: event " << iEvent << " cascade " << nCascade
2353 << " track " << negFromV0 << " has already been used!" << endl;
2356 // add it to the V0 array as well
2357 Double_t d0[2] = { -999., -99.};
2358 // counting is probably wrong
2359 new(V0s[jV0s++]) AliAODv0(vV0FromCascade, -999., -99., p_pos, p_neg, d0); // to be refined
2361 // Add the bachelor track from the cascade
2363 Int_t bachelor = cascade->GetBindex();
2365 if(!usedTrack[bachelor]) {
2367 usedTrack[bachelor] = kTRUE;
2369 AliESDtrack *esdTrack = esd->GetTrack(bachelor);
2370 esdTrack->GetPxPyPz(p);
2371 esdTrack->GetXYZ(pos);
2372 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2373 esdTrack->GetESDpid(pid);
2375 vcascade->AddDaughter(aodTrack =
2376 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2377 esdTrack->GetLabel(),
2383 (Short_t)esdTrack->Charge(),
2384 esdTrack->GetITSClusterMap(),
2387 kTRUE, // check if this is right
2388 kFALSE, // check if this is right
2389 AliAODTrack::kSecondary)
2391 aodTrack->ConvertAliPIDtoAODPID();
2394 cerr << "Error: event " << iEvent << " cascade " << nCascade
2395 << " track " << bachelor << " has already been used!" << endl;
2398 // Add the primary track of the cascade (if any)
2400 } // end of the loop on cascades
2404 for (Int_t nV0 = 0; nV0 < nV0s; ++nV0) {
2406 if (usedV0[nV0]) continue; // skip if aready added to the AOD
2408 AliESDv0 *v0 = esd->GetV0(nV0);
2410 v0->GetXYZ(pos[0], pos[1], pos[2]);
2411 v0->GetPosCov(covVtx);
2413 AliAODVertex * vV0 =
2414 new(vertices[jVertices++]) AliAODVertex(pos,
2416 v0->GetChi2V0(), // = chi2/NDF since NDF = 2*2-3
2420 primary->AddDaughter(vV0);
2422 Int_t posFromV0 = v0->GetPindex();
2423 Int_t negFromV0 = v0->GetNindex();
2425 // Add the positive tracks from the V0
2427 if (!usedTrack[posFromV0]) {
2429 usedTrack[posFromV0] = kTRUE;
2431 AliESDtrack *esdTrack = esd->GetTrack(posFromV0);
2432 esdTrack->GetPxPyPz(p_pos);
2433 esdTrack->GetXYZ(pos);
2434 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2435 esdTrack->GetESDpid(pid);
2437 vV0->AddDaughter(aodTrack =
2438 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2439 esdTrack->GetLabel(),
2445 (Short_t)esdTrack->Charge(),
2446 esdTrack->GetITSClusterMap(),
2449 kTRUE, // check if this is right
2450 kFALSE, // check if this is right
2451 AliAODTrack::kSecondary)
2453 aodTrack->ConvertAliPIDtoAODPID();
2456 cerr << "Error: event " << iEvent << " V0 " << nV0
2457 << " track " << posFromV0 << " has already been used!" << endl;
2460 // Add the negative tracks from the V0
2462 if (!usedTrack[negFromV0]) {
2464 usedTrack[negFromV0] = kTRUE;
2466 AliESDtrack *esdTrack = esd->GetTrack(negFromV0);
2467 esdTrack->GetPxPyPz(p_neg);
2468 esdTrack->GetXYZ(pos);
2469 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2470 esdTrack->GetESDpid(pid);
2472 vV0->AddDaughter(aodTrack =
2473 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2474 esdTrack->GetLabel(),
2480 (Short_t)esdTrack->Charge(),
2481 esdTrack->GetITSClusterMap(),
2484 kTRUE, // check if this is right
2485 kFALSE, // check if this is right
2486 AliAODTrack::kSecondary)
2488 aodTrack->ConvertAliPIDtoAODPID();
2491 cerr << "Error: event " << iEvent << " V0 " << nV0
2492 << " track " << negFromV0 << " has already been used!" << endl;
2495 // add it to the V0 array as well
2496 Double_t d0[2] = { 999., 99.};
2497 new(V0s[jV0s++]) AliAODv0(vV0, 999., 99., p_pos, p_neg, d0); // to be refined
2500 // end of the loop on V0s
2502 // Kinks: it is a big mess the access to the information in the kinks
2503 // The loop is on the tracks in order to find the mother and daugther of each kink
2506 for (Int_t iTrack=0; iTrack<nTracks; ++iTrack) {
2508 AliESDtrack * esdTrack = esd->GetTrack(iTrack);
2510 Int_t ikink = esdTrack->GetKinkIndex(0);
2513 // Negative kink index: mother, positive: daughter
2515 // Search for the second track of the kink
2517 for (Int_t jTrack = iTrack+1; jTrack<nTracks; ++jTrack) {
2519 AliESDtrack * esdTrack1 = esd->GetTrack(jTrack);
2521 Int_t jkink = esdTrack1->GetKinkIndex(0);
2523 if ( TMath::Abs(ikink)==TMath::Abs(jkink) ) {
2525 // The two tracks are from the same kink
2527 if (usedKink[TMath::Abs(ikink)-1]) continue; // skip used kinks
2530 Int_t idaughter = -1;
2532 if (ikink<0 && jkink>0) {
2537 else if (ikink>0 && jkink<0) {
2543 cerr << "Error: Wrong combination of kink indexes: "
2544 << ikink << " " << jkink << endl;
2548 // Add the mother track
2550 AliAODTrack * mother = NULL;
2552 if (!usedTrack[imother]) {
2554 usedTrack[imother] = kTRUE;
2556 AliESDtrack *esdTrack = esd->GetTrack(imother);
2557 esdTrack->GetPxPyPz(p);
2558 esdTrack->GetXYZ(pos);
2559 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2560 esdTrack->GetESDpid(pid);
2563 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2564 esdTrack->GetLabel(),
2570 (Short_t)esdTrack->Charge(),
2571 esdTrack->GetITSClusterMap(),
2574 kTRUE, // check if this is right
2575 kTRUE, // check if this is right
2576 AliAODTrack::kPrimary);
2577 primary->AddDaughter(mother);
2578 mother->ConvertAliPIDtoAODPID();
2581 cerr << "Error: event " << iEvent << " kink " << TMath::Abs(ikink)-1
2582 << " track " << imother << " has already been used!" << endl;
2585 // Add the kink vertex
2586 AliESDkink * kink = esd->GetKink(TMath::Abs(ikink)-1);
2588 AliAODVertex * vkink =
2589 new(vertices[jVertices++]) AliAODVertex(kink->GetPosition(),
2593 esdTrack->GetID(), // This is the track ID of the mother's track!
2594 AliAODVertex::kKink);
2595 // Add the daughter track
2597 AliAODTrack * daughter = NULL;
2599 if (!usedTrack[idaughter]) {
2601 usedTrack[idaughter] = kTRUE;
2603 AliESDtrack *esdTrack = esd->GetTrack(idaughter);
2604 esdTrack->GetPxPyPz(p);
2605 esdTrack->GetXYZ(pos);
2606 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2607 esdTrack->GetESDpid(pid);
2610 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2611 esdTrack->GetLabel(),
2617 (Short_t)esdTrack->Charge(),
2618 esdTrack->GetITSClusterMap(),
2621 kTRUE, // check if this is right
2622 kTRUE, // check if this is right
2623 AliAODTrack::kPrimary);
2624 vkink->AddDaughter(daughter);
2625 daughter->ConvertAliPIDtoAODPID();
2628 cerr << "Error: event " << iEvent << " kink " << TMath::Abs(ikink)-1
2629 << " track " << idaughter << " has already been used!" << endl;
2635 vertices.Expand(jVertices);
2637 // Tracks (primary and orphan)
2638 for (Int_t nTrack = 0; nTrack < nTracks; ++nTrack) {
2640 if (usedTrack[nTrack]) continue;
2642 AliESDtrack *esdTrack = esd->GetTrack(nTrack);
2643 esdTrack->GetPxPyPz(p);
2644 esdTrack->GetXYZ(pos);
2645 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2646 esdTrack->GetESDpid(pid);
2648 Float_t impactXY, impactZ;
2650 esdTrack->GetImpactParameters(impactXY,impactZ);
2653 // track inside the beam pipe
2655 primary->AddDaughter(aodTrack =
2656 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2657 esdTrack->GetLabel(),
2663 (Short_t)esdTrack->Charge(),
2664 esdTrack->GetITSClusterMap(),
2667 kTRUE, // check if this is right
2668 kTRUE, // check if this is right
2669 AliAODTrack::kPrimary)
2671 aodTrack->ConvertAliPIDtoAODPID();
2674 // outside the beam pipe: orphan track
2675 // Don't write them anymore!
2678 } // end of loop on tracks
2681 Int_t nMuTracks = esd->GetNumberOfMuonTracks();
2682 for (Int_t nMuTrack = 0; nMuTrack < nMuTracks; ++nMuTrack) {
2684 AliESDMuonTrack *esdMuTrack = esd->GetMuonTrack(nMuTrack);
2685 p[0] = esdMuTrack->Px();
2686 p[1] = esdMuTrack->Py();
2687 p[2] = esdMuTrack->Pz();
2688 pos[0] = primary->GetX();
2689 pos[1] = primary->GetY();
2690 pos[2] = primary->GetZ();
2692 // has to be changed once the muon pid is provided by the ESD
2693 for (Int_t i = 0; i < 10; pid[i++] = 0.); pid[AliAODTrack::kMuon]=1.;
2695 primary->AddDaughter(aodTrack =
2696 new(tracks[jTracks++]) AliAODTrack(0, // no ID provided
2697 0, // no label provided
2702 NULL, // no covariance matrix provided
2703 esdMuTrack->Charge(),
2704 0, // ITSClusterMap is set below
2707 kFALSE, // muon tracks are not used to fit the primary vtx
2708 kFALSE, // not used for vertex fit
2709 AliAODTrack::kPrimary)
2712 aodTrack->SetHitsPatternInTrigCh(esdMuTrack->GetHitsPatternInTrigCh());
2713 Int_t track2Trigger = esdMuTrack->GetMatchTrigger();
2714 aodTrack->SetMatchTrigger(track2Trigger);
2716 aodTrack->SetChi2MatchTrigger(esdMuTrack->GetChi2MatchTrigger());
2718 aodTrack->SetChi2MatchTrigger(0.);
2720 tracks.Expand(jTracks); // remove 'empty slots' due to unwritten tracks
2722 // Access to the AOD container of PMD clusters
2723 TClonesArray &pmdClusters = *(aod->GetPmdClusters());
2724 Int_t jPmdClusters=0;
2726 for (Int_t iPmd = 0; iPmd < nPmdClus; ++iPmd) {
2727 // file pmd clusters, to be revised!
2728 AliESDPmdTrack *pmdTrack = esd->GetPmdTrack(iPmd);
2731 Double_t pos[3] = { pmdTrack->GetClusterX(), pmdTrack->GetClusterY(), pmdTrack->GetClusterZ() };
2732 Double_t pid[9] = { 0., 0., 0., 0., 0., 0., 0., 0., 0. }; // to be revised!
2734 // assoc cluster not set
2735 new(pmdClusters[jPmdClusters++]) AliAODPmdCluster(iPmd, nLabel, label, pmdTrack->GetClusterADC(), pos, pid);
2738 // Access to the AOD container of clusters
2739 TClonesArray &caloClusters = *(aod->GetCaloClusters());
2742 for (Int_t iClust=0; iClust<nCaloClus; ++iClust) {
2744 AliESDCaloCluster * cluster = esd->GetCaloCluster(iClust);
2746 Int_t id = cluster->GetID();
2749 Float_t energy = cluster->E();
2750 cluster->GetPosition(posF);
2751 Char_t ttype=AliAODCluster::kUndef;
2753 if (cluster->GetClusterType() == AliESDCaloCluster::kPHOSCluster) {
2754 ttype=AliAODCluster::kPHOSNeutral;
2756 else if (cluster->GetClusterType() == AliESDCaloCluster::kEMCALClusterv1) {
2757 ttype = AliAODCluster::kEMCALClusterv1;
2761 AliAODCaloCluster *caloCluster = new(caloClusters[jClusters++]) AliAODCaloCluster(id,
2769 caloCluster->SetCaloCluster(); // to be refined!
2772 caloClusters.Expand(jClusters); // resize TObjArray to 'remove' slots for pseudo clusters
2773 // end of loop on calo clusters
2775 // fill EMCAL cell info
2776 if (esd->GetEMCALCells()) { // protection against missing ESD information
2777 AliESDCaloCells &esdEMcells = *(esd->GetEMCALCells());
2778 Int_t nEMcell = esdEMcells.GetNumberOfCells() ;
2780 AliAODCaloCells &aodEMcells = *(aod->GetEMCALCells());
2781 aodEMcells.CreateContainer(nEMcell);
2782 aodEMcells.SetType(AliAODCaloCells::kEMCAL);
2783 for (Int_t iCell = 0; iCell < nEMcell; iCell++) {
2784 aodEMcells.SetCell(iCell,esdEMcells.GetCellNumber(iCell),esdEMcells.GetAmplitude(iCell));
2789 // fill PHOS cell info
2790 if (esd->GetPHOSCells()) { // protection against missing ESD information
2791 AliESDCaloCells &esdPHcells = *(esd->GetPHOSCells());
2792 Int_t nPHcell = esdPHcells.GetNumberOfCells() ;
2794 AliAODCaloCells &aodPHcells = *(aod->GetPHOSCells());
2795 aodPHcells.CreateContainer(nPHcell);
2796 aodPHcells.SetType(AliAODCaloCells::kPHOS);
2797 for (Int_t iCell = 0; iCell < nPHcell; iCell++) {
2798 aodPHcells.SetCell(iCell,esdPHcells.GetCellNumber(iCell),esdPHcells.GetAmplitude(iCell));
2804 AliAODTracklets &SPDTracklets = *(aod->GetTracklets());
2805 const AliMultiplicity *mult = esd->GetMultiplicity();
2807 if (mult->GetNumberOfTracklets()>0) {
2808 SPDTracklets.CreateContainer(mult->GetNumberOfTracklets());
2810 for (Int_t n=0; n<mult->GetNumberOfTracklets(); n++) {
2811 SPDTracklets.SetTracklet(n, mult->GetTheta(n), mult->GetPhi(n), mult->GetDeltaPhi(n), mult->GetLabel(n));
2815 Printf("ERROR: AliMultiplicity could not be retrieved from ESD");
2818 delete [] usedTrack;
2822 // fill the tree for this event
2824 } // end of event loop
2826 aodTree->GetUserInfo()->Add(aod);
2828 // write the tree to the specified file
2829 aodFile = aodTree->GetCurrentFile();
2836 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2838 // Write space-points which are then used in the alignment procedures
2839 // For the moment only ITS, TRD and TPC
2841 // Load TOF clusters
2843 fLoader[3]->LoadRecPoints("read");
2844 TTree* tree = fLoader[3]->TreeR();
2846 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[3]));
2849 fTracker[3]->LoadClusters(tree);
2851 Int_t ntracks = esd->GetNumberOfTracks();
2852 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2854 AliESDtrack *track = esd->GetTrack(itrack);
2857 for (Int_t iDet = 3; iDet >= 0; iDet--)
2858 nsp += track->GetNcls(iDet);
2860 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2861 track->SetTrackPointArray(sp);
2863 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2864 AliTracker *tracker = fTracker[iDet];
2865 if (!tracker) continue;
2866 Int_t nspdet = track->GetNcls(iDet);
2867 if (nspdet <= 0) continue;
2868 track->GetClusters(iDet,idx);
2872 while (isp2 < nspdet) {
2874 TString dets = fgkDetectorName[iDet];
2875 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2876 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2877 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2878 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2879 isvalid = tracker->GetTrackPointTrackingError(idx[isp2],p,track);
2881 isvalid = tracker->GetTrackPoint(idx[isp2],p);
2884 const Int_t kNTPCmax = 159;
2885 if (iDet==1 && isp2>kNTPCmax) break; // to be fixed
2886 if (!isvalid) continue;
2887 sp->AddPoint(isptrack,&p); isptrack++; isp++;
2893 fTracker[3]->UnloadClusters();
2894 fLoader[3]->UnloadRecPoints();
2898 //_____________________________________________________________________________
2899 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2901 // The method reads the raw-data error log
2902 // accumulated within the rawReader.
2903 // It extracts the raw-data errors related to
2904 // the current event and stores them into
2905 // a TClonesArray inside the esd object.
2907 if (!fRawReader) return;
2909 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2911 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2913 if (iEvent != log->GetEventNumber()) continue;
2915 esd->AddRawDataErrorLog(log);
2920 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString fName){
2921 // Dump a file content into a char in TNamed
2923 in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2924 Int_t kBytes = (Int_t)in.tellg();
2925 printf("Size: %d \n",kBytes);
2928 char* memblock = new char [kBytes];
2929 in.seekg (0, ios::beg);
2930 in.read (memblock, kBytes);
2932 TString fData(memblock,kBytes);
2933 fn = new TNamed(fName,fData);
2934 printf("fData Size: %d \n",fData.Sizeof());
2935 printf("fName Size: %d \n",fName.Sizeof());
2936 printf("fn Size: %d \n",fn->Sizeof());
2940 AliInfo(Form("Could not Open %s\n",fPath.Data()));
2946 void AliReconstruction::TNamedToFile(TTree* fTree, TString fName){
2947 // This is not really needed in AliReconstruction at the moment
2948 // but can serve as a template
2950 TList *fList = fTree->GetUserInfo();
2951 TNamed *fn = (TNamed*)fList->FindObject(fName.Data());
2952 printf("fn Size: %d \n",fn->Sizeof());
2954 TString fTmp(fn->GetName()); // to be 100% sure in principle fName also works
2955 const char* cdata = fn->GetTitle();
2956 printf("fTmp Size %d\n",fTmp.Sizeof());
2958 int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2959 printf("calculated size %d\n",size);
2960 ofstream out(fName.Data(),ios::out | ios::binary);
2961 out.write(cdata,size);
2966 //_____________________________________________________________________________
2967 AliQADataMakerRec * AliReconstruction::GetQADataMaker(Int_t iDet)
2969 // get the quality assurance data maker object and the loader for a detector
2971 if (fQADataMaker[iDet])
2972 return fQADataMaker[iDet];
2974 AliQADataMakerRec * qadm = NULL;
2975 if (iDet == fgkNDetectors) { //Global QA
2976 qadm = new AliGlobalQADataMaker();
2977 fQADataMaker[iDet] = qadm;
2981 // load the QA data maker object
2982 TPluginManager* pluginManager = gROOT->GetPluginManager();
2983 TString detName = fgkDetectorName[iDet];
2984 TString qadmName = "Ali" + detName + "QADataMakerRec";
2985 if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT"))
2988 // first check if a plugin is defined for the quality assurance data maker
2989 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2990 // if not, add a plugin for it
2991 if (!pluginHandler) {
2992 AliDebug(1, Form("defining plugin for %s", qadmName.Data()));
2993 TString libs = gSystem->GetLibraries();
2994 if (libs.Contains("lib" + detName + "base.so") ||
2995 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2996 pluginManager->AddHandler("AliQADataMakerRec", detName,
2997 qadmName, detName + "qadm", qadmName + "()");
2999 pluginManager->AddHandler("AliQADataMakerRec", detName,
3000 qadmName, detName, qadmName + "()");
3002 pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
3004 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
3005 qadm = (AliQADataMakerRec *) pluginHandler->ExecPlugin(0);
3008 fQADataMaker[iDet] = qadm;
3013 //_____________________________________________________________________________
3014 Bool_t AliReconstruction::RunQA(const char* detectors, AliESDEvent *& esd)
3016 // run the Quality Assurance data producer
3018 AliCodeTimerAuto("")
3019 TString detStr = detectors;
3020 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3021 if (!IsSelected(fgkDetectorName[iDet], detStr))
3023 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
3026 AliCodeTimerStart(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
3027 AliInfo(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
3029 qadm->Exec(AliQA::kESDS, esd) ;
3032 AliCodeTimerStop(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
3034 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
3035 AliError(Form("the following detectors were not found: %s",
3045 //_____________________________________________________________________________
3046 void AliReconstruction::CheckQA()
3048 // check the QA of SIM for this run and remove the detectors
3049 // with status Fatal
3051 TString newRunLocalReconstruction ;
3052 TString newRunTracking ;
3053 TString newFillESD ;
3055 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
3056 TString detName(AliQA::GetDetName(iDet)) ;
3057 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX(iDet)) ;
3058 if ( qa->IsSet(AliQA::DETECTORINDEX(iDet), AliQA::kSIM, AliQA::kFATAL)) {
3059 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
3061 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
3062 fRunLocalReconstruction.Contains("ALL") ) {
3063 newRunLocalReconstruction += detName ;
3064 newRunLocalReconstruction += " " ;
3066 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
3067 fRunTracking.Contains("ALL") ) {
3068 newRunTracking += detName ;
3069 newRunTracking += " " ;
3071 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
3072 fFillESD.Contains("ALL") ) {
3073 newFillESD += detName ;
3078 fRunLocalReconstruction = newRunLocalReconstruction ;
3079 fRunTracking = newRunTracking ;
3080 fFillESD = newFillESD ;
3083 //_____________________________________________________________________________
3084 Int_t AliReconstruction::GetDetIndex(const char* detector)
3086 // return the detector index corresponding to detector
3088 for (index = 0; index < fgkNDetectors ; index++) {
3089 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
3094 //_____________________________________________________________________________
3095 Bool_t AliReconstruction::FinishPlaneEff() {
3097 // Here execute all the necessary operationis, at the end of the tracking phase,
3098 // in case that evaluation of PlaneEfficiencies was required for some detector.
3099 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
3101 // This Preliminary version works only FOR ITS !!!!!
3102 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3105 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3108 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3109 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
3110 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3111 if(fTracker[iDet]) ret=fTracker[iDet]->GetPlaneEff()->WriteIntoCDB();
3115 //_____________________________________________________________________________
3116 Bool_t AliReconstruction::InitPlaneEff() {
3118 // Here execute all the necessary operations, before of the tracking phase,
3119 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
3120 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
3121 // which should be updated/recalculated.
3123 // This Preliminary version will work 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
3129 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));