1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // class for running the reconstruction //
22 // Clusters and tracks are created for all detectors and all events by //
25 // AliReconstruction rec; //
28 // The Run method returns kTRUE in case of successful execution. //
30 // If the input to the reconstruction are not simulated digits but raw data, //
31 // this can be specified by an argument of the Run method or by the method //
33 // rec.SetInput("..."); //
35 // The input formats and the corresponding argument are: //
36 // - DDL raw data files: directory name, ends with "/" //
37 // - raw data root file: root file name, extension ".root" //
38 // - raw data DATE file: DATE file name, any other non-empty string //
39 // - MC root files : empty string, default //
41 // By default all events are reconstructed. The reconstruction can be //
42 // limited to a range of events by giving the index of the first and the //
43 // last event as an argument to the Run method or by calling //
45 // rec.SetEventRange(..., ...); //
47 // The index -1 (default) can be used for the last event to indicate no //
48 // upper limit of the event range. //
50 // In case of raw-data reconstruction the user can modify the default //
51 // number of events per digits/clusters/tracks file. In case the option //
52 // is not used the number is set 1. In case the user provides 0, than //
53 // the number of events is equal to the number of events inside the //
54 // raw-data file (i.e. one digits/clusters/tracks file): //
56 // rec.SetNumberOfEventsPerFile(...); //
59 // The name of the galice file can be changed from the default //
60 // "galice.root" by passing it as argument to the AliReconstruction //
61 // constructor or by //
63 // rec.SetGAliceFile("..."); //
65 // The local reconstruction can be switched on or off for individual //
68 // rec.SetRunLocalReconstruction("..."); //
70 // The argument is a (case sensitive) string with the names of the //
71 // detectors separated by a space. The special string "ALL" selects all //
72 // available detectors. This is the default. //
74 // The reconstruction of the primary vertex position can be switched off by //
76 // rec.SetRunVertexFinder(kFALSE); //
78 // The tracking and the creation of ESD tracks can be switched on for //
79 // selected detectors by //
81 // rec.SetRunTracking("..."); //
83 // Uniform/nonuniform field tracking switches (default: uniform field) //
85 // rec.SetUniformFieldTracking(); ( rec.SetUniformFieldTracking(kFALSE); ) //
87 // The filling of additional ESD information can be steered by //
89 // rec.SetFillESD("..."); //
91 // Again, for both methods the string specifies the list of detectors. //
92 // The default is "ALL". //
94 // The call of the shortcut method //
96 // rec.SetRunReconstruction("..."); //
98 // is equivalent to calling SetRunLocalReconstruction, SetRunTracking and //
99 // SetFillESD with the same detector selecting string as argument. //
101 // The reconstruction requires digits or raw data as input. For the creation //
102 // of digits and raw data have a look at the class AliSimulation. //
104 // The input data of a detector can be replaced by the corresponding HLT //
105 // data by calling (usual detector string) //
106 // SetUseHLTData("..."); //
108 // For debug purposes the method SetCheckPointLevel can be used. If the //
109 // argument is greater than 0, files with ESD events will be written after //
110 // selected steps of the reconstruction for each event: //
111 // level 1: after tracking and after filling of ESD (final) //
112 // level 2: in addition after each tracking step //
113 // level 3: in addition after the filling of ESD for each detector //
114 // If a final check point file exists for an event, this event will be //
115 // skipped in the reconstruction. The tracking and the filling of ESD for //
116 // a detector will be skipped as well, if the corresponding check point //
117 // file exists. The ESD event will then be loaded from the file instead. //
119 ///////////////////////////////////////////////////////////////////////////////
126 #include <TPluginManager.h>
127 #include <TGeoManager.h>
128 #include <TLorentzVector.h>
131 #include <TObjArray.h>
133 #include "AliReconstruction.h"
134 #include "AliCodeTimer.h"
135 #include "AliReconstructor.h"
137 #include "AliRunLoader.h"
139 #include "AliRawReaderFile.h"
140 #include "AliRawReaderDate.h"
141 #include "AliRawReaderRoot.h"
142 #include "AliRawEventHeaderBase.h"
143 #include "AliESDEvent.h"
144 #include "AliESDMuonTrack.h"
145 #include "AliESDfriend.h"
146 #include "AliESDVertex.h"
147 #include "AliESDcascade.h"
148 #include "AliESDkink.h"
149 #include "AliESDtrack.h"
150 #include "AliESDCaloCluster.h"
151 #include "AliESDCaloCells.h"
152 #include "AliMultiplicity.h"
153 #include "AliTracker.h"
154 #include "AliVertexer.h"
155 #include "AliVertexerTracks.h"
156 #include "AliV0vertexer.h"
157 #include "AliCascadeVertexer.h"
158 #include "AliHeader.h"
159 #include "AliGenEventHeader.h"
161 #include "AliESDpid.h"
162 #include "AliESDtrack.h"
163 #include "AliESDPmdTrack.h"
165 #include "AliESDTagCreator.h"
166 #include "AliAODTagCreator.h"
168 #include "AliGeomManager.h"
169 #include "AliTrackPointArray.h"
170 #include "AliCDBManager.h"
171 #include "AliCDBStorage.h"
172 #include "AliCDBEntry.h"
173 #include "AliAlignObj.h"
175 #include "AliCentralTrigger.h"
176 #include "AliTriggerConfiguration.h"
177 #include "AliTriggerClass.h"
178 #include "AliCTPRawStream.h"
180 #include "AliAODEvent.h"
181 #include "AliAODHeader.h"
182 #include "AliAODTrack.h"
183 #include "AliAODVertex.h"
184 #include "AliAODv0.h"
185 #include "AliAODJet.h"
186 #include "AliAODCaloCells.h"
187 #include "AliAODCaloCluster.h"
188 #include "AliAODPmdCluster.h"
189 #include "AliAODFmdCluster.h"
190 #include "AliAODTracklets.h"
192 #include "AliQADataMakerRec.h"
193 #include "AliGlobalQADataMaker.h"
195 #include "AliQADataMakerSteer.h"
197 #include "AliPlaneEff.h"
199 #include "AliSysInfo.h" // memory snapshots
200 #include "AliRawHLTManager.h"
203 ClassImp(AliReconstruction)
206 //_____________________________________________________________________________
207 const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
209 //_____________________________________________________________________________
210 AliReconstruction::AliReconstruction(const char* gAliceFilename,
211 const char* name, const char* title) :
214 fUniformField(kTRUE),
215 fRunVertexFinder(kTRUE),
216 fRunVertexFinderTracks(kTRUE),
217 fRunHLTTracking(kFALSE),
218 fRunMuonTracking(kFALSE),
220 fRunCascadeFinder(kTRUE),
221 fStopOnError(kFALSE),
222 fWriteAlignmentData(kFALSE),
223 fWriteESDfriend(kFALSE),
225 fFillTriggerESD(kTRUE),
233 fRunLocalReconstruction("ALL"),
236 fUseTrackingErrorsForAlignment(""),
237 fGAliceFileName(gAliceFilename),
242 fNumberOfEventsPerFile(1),
245 fLoadAlignFromCDB(kTRUE),
246 fLoadAlignData("ALL"),
252 fParentRawReader(NULL),
255 fDiamondProfile(NULL),
256 fMeanVertexConstraint(kTRUE),
260 fAlignObjArray(NULL),
263 fInitCDBCalled(kFALSE),
264 fSetRunNumberFromDataCalled(kFALSE),
271 // create reconstruction object with default parameters
273 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
274 fReconstructor[iDet] = NULL;
275 fLoader[iDet] = NULL;
276 fTracker[iDet] = NULL;
277 fQADataMaker[iDet] = NULL;
278 fQACycles[iDet] = 999999;
280 fQADataMaker[fgkNDetectors]=NULL; //Global QA
284 //_____________________________________________________________________________
285 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
288 fUniformField(rec.fUniformField),
289 fRunVertexFinder(rec.fRunVertexFinder),
290 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
291 fRunHLTTracking(rec.fRunHLTTracking),
292 fRunMuonTracking(rec.fRunMuonTracking),
293 fRunV0Finder(rec.fRunV0Finder),
294 fRunCascadeFinder(rec.fRunCascadeFinder),
295 fStopOnError(rec.fStopOnError),
296 fWriteAlignmentData(rec.fWriteAlignmentData),
297 fWriteESDfriend(rec.fWriteESDfriend),
298 fWriteAOD(rec.fWriteAOD),
299 fFillTriggerESD(rec.fFillTriggerESD),
301 fCleanESD(rec.fCleanESD),
302 fV0DCAmax(rec.fV0DCAmax),
303 fV0CsPmin(rec.fV0CsPmin),
307 fRunLocalReconstruction(rec.fRunLocalReconstruction),
308 fRunTracking(rec.fRunTracking),
309 fFillESD(rec.fFillESD),
310 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
311 fGAliceFileName(rec.fGAliceFileName),
313 fEquipIdMap(rec.fEquipIdMap),
314 fFirstEvent(rec.fFirstEvent),
315 fLastEvent(rec.fLastEvent),
316 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
319 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
320 fLoadAlignData(rec.fLoadAlignData),
321 fESDPar(rec.fESDPar),
322 fUseHLTData(rec.fUseHLTData),
326 fParentRawReader(NULL),
329 fDiamondProfile(NULL),
330 fMeanVertexConstraint(rec.fMeanVertexConstraint),
334 fAlignObjArray(rec.fAlignObjArray),
335 fCDBUri(rec.fCDBUri),
337 fInitCDBCalled(rec.fInitCDBCalled),
338 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
340 fRunGlobalQA(rec.fRunGlobalQA),
341 fInLoopQA(rec.fInLoopQA),
342 fRunPlaneEff(rec.fRunPlaneEff)
346 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
347 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
349 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
350 fReconstructor[iDet] = NULL;
351 fLoader[iDet] = NULL;
352 fTracker[iDet] = NULL;
353 fQADataMaker[iDet] = NULL;
354 fQACycles[iDet] = rec.fQACycles[iDet];
356 fQADataMaker[fgkNDetectors]=NULL; //Global QA
357 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
358 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
362 //_____________________________________________________________________________
363 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
365 // assignment operator
367 this->~AliReconstruction();
368 new(this) AliReconstruction(rec);
372 //_____________________________________________________________________________
373 AliReconstruction::~AliReconstruction()
379 fSpecCDBUri.Delete();
381 AliCodeTimer::Instance()->Print();
384 //_____________________________________________________________________________
385 void AliReconstruction::InitCDB()
387 // activate a default CDB storage
388 // First check if we have any CDB storage set, because it is used
389 // to retrieve the calibration and alignment constants
391 if (fInitCDBCalled) return;
392 fInitCDBCalled = kTRUE;
394 AliCDBManager* man = AliCDBManager::Instance();
395 if (man->IsDefaultStorageSet())
397 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
398 AliWarning("Default CDB storage has been already set !");
399 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
400 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
401 fCDBUri = man->GetDefaultStorage()->GetURI();
404 if (fCDBUri.Length() > 0)
406 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
407 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
408 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
410 fCDBUri="local://$ALICE_ROOT";
411 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
412 AliWarning("Default CDB storage not yet set !!!!");
413 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
414 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
417 man->SetDefaultStorage(fCDBUri);
420 // Now activate the detector specific CDB storage locations
421 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
422 TObject* obj = fSpecCDBUri[i];
424 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
425 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
426 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
427 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
432 //_____________________________________________________________________________
433 void AliReconstruction::SetDefaultStorage(const char* uri) {
434 // Store the desired default CDB storage location
435 // Activate it later within the Run() method
441 //_____________________________________________________________________________
442 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
443 // Store a detector-specific CDB storage location
444 // Activate it later within the Run() method
446 AliCDBPath aPath(calibType);
447 if(!aPath.IsValid()){
448 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
449 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
450 if(!strcmp(calibType, fgkDetectorName[iDet])) {
451 aPath.SetPath(Form("%s/*", calibType));
452 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
456 if(!aPath.IsValid()){
457 AliError(Form("Not a valid path or detector: %s", calibType));
462 // // check that calibType refers to a "valid" detector name
463 // Bool_t isDetector = kFALSE;
464 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
465 // TString detName = fgkDetectorName[iDet];
466 // if(aPath.GetLevel0() == detName) {
467 // isDetector = kTRUE;
473 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
477 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
478 if (obj) fSpecCDBUri.Remove(obj);
479 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
483 //_____________________________________________________________________________
484 Bool_t AliReconstruction::SetRunNumberFromData()
486 // The method is called in Run() in order
487 // to set a correct run number.
488 // In case of raw data reconstruction the
489 // run number is taken from the raw data header
491 if (fSetRunNumberFromDataCalled) return kTRUE;
492 fSetRunNumberFromDataCalled = kTRUE;
494 AliCDBManager* man = AliCDBManager::Instance();
496 if(man->GetRun() > 0) {
497 AliWarning("Run number is taken from event header! Ignoring settings in AliCDBManager!");
501 AliError("No run loader is found !");
504 // read run number from gAlice
505 if(fRunLoader->GetAliRun())
506 AliCDBManager::Instance()->SetRun(fRunLoader->GetHeader()->GetRun());
509 if(fRawReader->NextEvent()) {
510 AliCDBManager::Instance()->SetRun(fRawReader->GetRunNumber());
511 fRawReader->RewindEvents();
514 AliError("No raw-data events found !");
519 AliError("Neither gAlice nor RawReader objects are found !");
529 //_____________________________________________________________________________
530 void AliReconstruction::SetCDBLock() {
531 // Set CDB lock: from now on it is forbidden to reset the run number
532 // or the default storage or to activate any further storage!
534 AliCDBManager::Instance()->SetLock(1);
537 //_____________________________________________________________________________
538 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
540 // Read the alignment objects from CDB.
541 // Each detector is supposed to have the
542 // alignment objects in DET/Align/Data CDB path.
543 // All the detector objects are then collected,
544 // sorted by geometry level (starting from ALIC) and
545 // then applied to the TGeo geometry.
546 // Finally an overlaps check is performed.
548 // Load alignment data from CDB and fill fAlignObjArray
549 if(fLoadAlignFromCDB){
551 TString detStr = detectors;
552 TString loadAlObjsListOfDets = "";
554 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
555 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
556 loadAlObjsListOfDets += fgkDetectorName[iDet];
557 loadAlObjsListOfDets += " ";
558 } // end loop over detectors
559 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
560 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
562 // Check if the array with alignment objects was
563 // provided by the user. If yes, apply the objects
564 // to the present TGeo geometry
565 if (fAlignObjArray) {
566 if (gGeoManager && gGeoManager->IsClosed()) {
567 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
568 AliError("The misalignment of one or more volumes failed!"
569 "Compare the list of simulated detectors and the list of detector alignment data!");
574 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
580 delete fAlignObjArray; fAlignObjArray=0;
585 //_____________________________________________________________________________
586 void AliReconstruction::SetGAliceFile(const char* fileName)
588 // set the name of the galice file
590 fGAliceFileName = fileName;
593 //_____________________________________________________________________________
594 void AliReconstruction::SetOption(const char* detector, const char* option)
596 // set options for the reconstruction of a detector
598 TObject* obj = fOptions.FindObject(detector);
599 if (obj) fOptions.Remove(obj);
600 fOptions.Add(new TNamed(detector, option));
604 //_____________________________________________________________________________
605 Bool_t AliReconstruction::Run(const char* input, Bool_t IsOnline)
607 // run the reconstruction
613 if (!input) input = fInput.Data();
614 TString fileName(input);
615 if (fileName.EndsWith("/")) {
616 fRawReader = new AliRawReaderFile(fileName);
617 } else if (fileName.EndsWith(".root")) {
618 fRawReader = new AliRawReaderRoot(fileName);
619 } else if (!fileName.IsNull()) {
620 fRawReader = new AliRawReaderDate(fileName);
625 AliError("Null pointer to the event structure!");
628 fRawReader = new AliRawReaderDate((void *)input);
631 if (!fEquipIdMap.IsNull() && fRawReader)
632 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
634 if (!fUseHLTData.IsNull()) {
635 // create the RawReaderHLT which performs redirection of HLT input data for
636 // the specified detectors
637 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
639 fParentRawReader=fRawReader;
640 fRawReader=pRawReader;
642 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
646 AliSysInfo::AddStamp("Start");
647 // get the run loader
648 if (!InitRunLoader()) return kFALSE;
649 AliSysInfo::AddStamp("LoadLoader");
651 // Initialize the CDB storage
654 AliSysInfo::AddStamp("LoadCDB");
656 // Set run number in CDBManager (if it is not already set by the user)
657 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
659 // Set CDB lock: from now on it is forbidden to reset the run number
660 // or the default storage or to activate any further storage!
663 // Import ideal TGeo geometry and apply misalignment
665 TString geom(gSystem->DirName(fGAliceFileName));
666 geom += "/geometry.root";
667 AliGeomManager::LoadGeometry(geom.Data());
668 if (!gGeoManager) if (fStopOnError) return kFALSE;
671 if (!MisalignGeometry(fLoadAlignData)) if (fStopOnError) return kFALSE;
672 AliSysInfo::AddStamp("LoadGeom");
675 AliQADataMakerSteer qas ;
676 if (fRunQA && fRawReader) qas.Run(fRunLocalReconstruction, fRawReader) ;
677 // checking the QA of previous steps
681 // local reconstruction
682 if (!fRunLocalReconstruction.IsNull()) {
683 if (!RunLocalReconstruction(fRunLocalReconstruction)) {
684 if (fStopOnError) {CleanUp(); return kFALSE;}
690 if (fRunVertexFinder && !CreateVertexer()) {
696 AliSysInfo::AddStamp("Vertexer");
699 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
705 AliSysInfo::AddStamp("LoadTrackers");
707 // get the possibly already existing ESD file and tree
708 AliESDEvent* esd = new AliESDEvent(); AliESDEvent* hltesd = new AliESDEvent();
709 TFile* fileOld = NULL;
710 TTree* treeOld = NULL; TTree *hlttreeOld = NULL;
711 if (!gSystem->AccessPathName("AliESDs.root")){
712 gSystem->CopyFile("AliESDs.root", "AliESDs.old.root", kTRUE);
713 fileOld = TFile::Open("AliESDs.old.root");
714 if (fileOld && fileOld->IsOpen()) {
715 treeOld = (TTree*) fileOld->Get("esdTree");
716 if (treeOld)esd->ReadFromTree(treeOld);
717 hlttreeOld = (TTree*) fileOld->Get("HLTesdTree");
718 if (hlttreeOld) hltesd->ReadFromTree(hlttreeOld);
722 // create the ESD output file and tree
723 TFile* file = TFile::Open("AliESDs.root", "RECREATE");
724 file->SetCompressionLevel(2);
725 if (!file->IsOpen()) {
726 AliError("opening AliESDs.root failed");
727 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
730 TTree* tree = new TTree("esdTree", "Tree with ESD objects");
731 esd = new AliESDEvent();
732 esd->CreateStdContent();
733 esd->WriteToTree(tree);
735 TTree* hlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
736 hltesd = new AliESDEvent();
737 hltesd->CreateStdContent();
738 hltesd->WriteToTree(hlttree);
741 delete esd; delete hltesd;
742 esd = NULL; hltesd = NULL;
744 // create the branch with ESD additions
748 AliESDfriend *esdf = 0;
749 if (fWriteESDfriend) {
750 esdf = new AliESDfriend();
751 TBranch *br=tree->Branch("ESDfriend.","AliESDfriend", &esdf);
752 br->SetFile("AliESDfriends.root");
753 esd->AddObject(esdf);
757 // Get the GRP CDB entry
758 AliCDBEntry* entryGRP = AliCDBManager::Instance()->Get("GRP/GRP/Data");
761 fGRPList = dynamic_cast<TList*> (entryGRP->GetObject());
763 AliError("No GRP entry found in OCDB!");
766 // Get the diamond profile from OCDB
767 AliCDBEntry* entry = AliCDBManager::Instance()
768 ->Get("GRP/Calib/MeanVertex");
771 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
773 AliError("No diamond profile found in OCDB!");
776 AliVertexerTracks tVertexer(AliTracker::GetBz());
777 if(fDiamondProfile && fMeanVertexConstraint) tVertexer.SetVtxStart(fDiamondProfile);
779 if (fRawReader) fRawReader->RewindEvents();
782 gSystem->GetProcInfo(&ProcInfo);
783 AliInfo(Form("Current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual));
786 //Initialize the QA and start of cycle for out-of-cycle QA
788 TString detStr(fFillESD);
789 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
790 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
791 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
793 AliInfo(Form("Initializing the QA data maker for %s",
794 fgkDetectorName[iDet]));
795 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
796 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
798 qadm->StartOfCycle(AliQA::kRECPOINTS);
799 qadm->StartOfCycle(AliQA::kESDS,"same");
804 AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
805 AliInfo(Form("Initializing the global QA data maker"));
807 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
808 AliTracker::SetResidualsArray(arr);
809 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
811 qadm->StartOfCycle(AliQA::kRECPOINTS);
812 qadm->StartOfCycle(AliQA::kESDS,"same");
816 //Initialize the Plane Efficiency framework
817 if (fRunPlaneEff && !InitPlaneEff()) {
818 if(fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
821 //******* The loop over events
822 for (Int_t iEvent = 0; iEvent < fRunLoader->GetNumberOfEvents(); iEvent++) {
823 if (fRawReader) fRawReader->NextEvent();
824 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
825 // copy old ESD to the new one
827 esd->ReadFromTree(treeOld);
828 treeOld->GetEntry(iEvent);
832 esd->ReadFromTree(hlttreeOld);
833 hlttreeOld->GetEntry(iEvent);
839 AliInfo(Form("processing event %d", iEvent));
841 //Start of cycle for the in-loop QA
844 TString detStr(fFillESD);
845 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
846 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
847 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
849 qadm->StartOfCycle(AliQA::kRECPOINTS);
850 qadm->StartOfCycle(AliQA::kESDS, "same") ;
854 AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
855 qadm->StartOfCycle(AliQA::kRECPOINTS);
856 qadm->StartOfCycle(AliQA::kESDS,"same");
860 fRunLoader->GetEvent(iEvent);
863 sprintf(aFileName, "ESD_%d.%d_final.root",
864 fRunLoader->GetHeader()->GetRun(),
865 fRunLoader->GetHeader()->GetEventNrInRun());
866 if (!gSystem->AccessPathName(aFileName)) continue;
868 // local signle event reconstruction
869 if (!fRunLocalReconstruction.IsNull()) {
870 if (!RunLocalEventReconstruction(fRunLocalReconstruction)) {
871 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
875 esd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
876 hltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
877 esd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
878 hltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
880 // Set magnetic field from the tracker
881 esd->SetMagneticField(AliTracker::GetBz());
882 hltesd->SetMagneticField(AliTracker::GetBz());
886 // Fill raw-data error log into the ESD
887 if (fRawReader) FillRawDataErrorLog(iEvent,esd);
890 if (fRunVertexFinder) {
891 if (!ReadESD(esd, "vertex")) {
892 if (!RunVertexFinder(esd)) {
893 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
895 if (fCheckPointLevel > 0) WriteESD(esd, "vertex");
900 if (!fRunTracking.IsNull()) {
901 if (fRunHLTTracking) {
902 hltesd->SetPrimaryVertexSPD(esd->GetVertex());
903 if (!RunHLTTracking(hltesd)) {
904 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
910 if (!fRunTracking.IsNull()) {
911 if (fRunMuonTracking) {
912 if (!RunMuonTracking(esd)) {
913 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
919 if (!fRunTracking.IsNull()) {
920 if (!ReadESD(esd, "tracking")) {
921 if (!RunTracking(esd)) {
922 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
924 if (fCheckPointLevel > 0) WriteESD(esd, "tracking");
929 if (!fFillESD.IsNull()) {
930 if (!FillESD(esd, fFillESD)) {
931 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
935 // fill Event header information from the RawEventHeader
936 if (fRawReader){FillRawEventHeaderESD(esd);}
939 AliESDpid::MakePID(esd);
940 if (fCheckPointLevel > 1) WriteESD(esd, "PID");
942 if (fFillTriggerESD) {
943 if (!ReadESD(esd, "trigger")) {
944 if (!FillTriggerESD(esd)) {
945 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
947 if (fCheckPointLevel > 1) WriteESD(esd, "trigger");
954 // Propagate track to the beam pipe (if not laready done by ITS)
956 const Int_t ntracks = esd->GetNumberOfTracks();
957 const Double_t kBz = esd->GetMagneticField();
958 const Double_t kRadius = 2.8; //something less than the beam pipe radius
961 UShort_t *selectedIdx=new UShort_t[ntracks];
963 for (Int_t itrack=0; itrack<ntracks; itrack++){
964 const Double_t kMaxStep = 5; //max step over the material
967 AliESDtrack *track = esd->GetTrack(itrack);
968 if (!track) continue;
970 AliExternalTrackParam *tpcTrack =
971 (AliExternalTrackParam *)track->GetTPCInnerParam();
975 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kTRUE);
977 Int_t n=trkArray.GetEntriesFast();
978 selectedIdx[n]=track->GetID();
979 trkArray.AddLast(tpcTrack);
982 if (track->GetX() < kRadius) continue;
985 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kTRUE);
987 track->RelateToVertex(esd->GetPrimaryVertexSPD(), kBz, kRadius);
992 // Improve the reconstructed primary vertex position using the tracks
994 TObject *obj = fOptions.FindObject("ITS");
996 TString optITS = obj->GetTitle();
997 if (optITS.Contains("cosmics") || optITS.Contains("COSMICS"))
998 fRunVertexFinderTracks=kFALSE;
1000 if (fRunVertexFinderTracks) {
1001 // TPC + ITS primary vertex
1002 AliESDVertex *pvtx=tVertexer.FindPrimaryVertex(esd);
1004 if (pvtx->GetStatus()) {
1005 esd->SetPrimaryVertex(pvtx);
1006 for (Int_t i=0; i<ntracks; i++) {
1007 AliESDtrack *t = esd->GetTrack(i);
1008 t->RelateToVertex(pvtx, kBz, kRadius);
1013 // TPC-only primary vertex
1014 pvtx=tVertexer.FindPrimaryVertex(&trkArray,selectedIdx);
1016 if (pvtx->GetStatus()) {
1017 esd->SetPrimaryVertexTPC(pvtx);
1018 Int_t nsel=trkArray.GetEntriesFast();
1019 for (Int_t i=0; i<nsel; i++) {
1020 AliExternalTrackParam *t =
1021 (AliExternalTrackParam *)trkArray.UncheckedAt(i);
1022 t->PropagateToDCA(pvtx, kBz, kRadius);
1028 delete[] selectedIdx;
1030 if(fDiamondProfile) esd->SetDiamond(fDiamondProfile);
1035 AliV0vertexer vtxer;
1036 vtxer.Tracks2V0vertices(esd);
1038 if (fRunCascadeFinder) {
1040 AliCascadeVertexer cvtxer;
1041 cvtxer.V0sTracks2CascadeVertices(esd);
1046 if (fCleanESD) CleanESD(esd);
1049 AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
1050 if (qadm) qadm->Exec(AliQA::kESDS, esd);
1053 if (fWriteESDfriend) {
1054 esdf->~AliESDfriend();
1055 new (esdf) AliESDfriend(); // Reset...
1056 esd->GetESDfriend(esdf);
1063 if (fCheckPointLevel > 0) WriteESD(esd, "final");
1066 if (fWriteESDfriend) {
1067 esdf->~AliESDfriend();
1068 new (esdf) AliESDfriend(); // Reset...
1071 gSystem->GetProcInfo(&ProcInfo);
1072 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1075 // End of cycle for the in-loop QA
1078 RunQA(fFillESD.Data(), esd);
1079 TString detStr(fFillESD);
1080 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1081 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1082 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1083 if (!qadm) continue;
1084 qadm->EndOfCycle(AliQA::kRECPOINTS);
1085 qadm->EndOfCycle(AliQA::kESDS);
1090 AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
1092 qadm->EndOfCycle(AliQA::kRECPOINTS);
1093 qadm->EndOfCycle(AliQA::kESDS);
1099 //******** End of the loop over events
1103 tree->GetUserInfo()->Add(esd);
1104 hlttree->GetUserInfo()->Add(hltesd);
1106 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1107 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1109 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1110 cdbMapCopy->SetOwner(1);
1111 cdbMapCopy->SetName("cdbMap");
1112 TIter iter(cdbMap->GetTable());
1115 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1116 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1117 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1118 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1121 TList *cdbListCopy = new TList();
1122 cdbListCopy->SetOwner(1);
1123 cdbListCopy->SetName("cdbList");
1125 TIter iter2(cdbList);
1128 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1129 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1132 tree->GetUserInfo()->Add(cdbMapCopy);
1133 tree->GetUserInfo()->Add(cdbListCopy);
1136 if(fESDPar.Contains("ESD.par")){
1137 AliInfo("Attaching ESD.par to Tree");
1138 TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1139 tree->GetUserInfo()->Add(fn);
1145 if (fWriteESDfriend)
1146 tree->SetBranchStatus("ESDfriend*",0);
1147 // we want to have only one tree version number
1148 tree->Write(tree->GetName(),TObject::kOverwrite);
1151 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1152 if (fRunPlaneEff && !FinishPlaneEff()) {
1153 AliWarning("Finish PlaneEff evaluation failed");
1157 CleanUp(file, fileOld);
1160 TFile *esdFile = TFile::Open("AliESDs.root", "READONLY");
1161 TFile *aodFile = TFile::Open("AliAOD.root", "RECREATE");
1162 ESDFile2AODFile(esdFile, aodFile);
1167 // Create tags for the events in the ESD tree (the ESD tree is always present)
1168 // In case of empty events the tags will contain dummy values
1169 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1170 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPList);
1172 AliAODTagCreator *aodtagCreator = new AliAODTagCreator();
1173 aodtagCreator->CreateAODTags(fFirstEvent,fLastEvent,fGRPList);
1176 //Finish QA and end of cycle for out-of-loop QA
1179 qas.Run(fRunLocalReconstruction.Data(), AliQA::kRECPOINTS);
1181 qas.Run(fRunTracking.Data(), AliQA::kESDS);
1184 AliQADataMakerRec *qadm = GetQADataMaker(fgkNDetectors);
1186 qadm->EndOfCycle(AliQA::kRECPOINTS);
1187 qadm->EndOfCycle(AliQA::kESDS);
1193 // Cleanup of CDB manager: cache and active storages!
1194 AliCDBManager::Instance()->ClearCache();
1201 //_____________________________________________________________________________
1202 Bool_t AliReconstruction::RunLocalReconstruction(const TString& /*detectors*/)
1204 // run the local reconstruction
1205 static Int_t eventNr=0;
1206 AliCodeTimerAuto("")
1208 // AliCDBManager* man = AliCDBManager::Instance();
1209 // Bool_t origCache = man->GetCacheFlag();
1211 // TString detStr = detectors;
1212 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1213 // if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1214 // AliReconstructor* reconstructor = GetReconstructor(iDet);
1215 // if (!reconstructor) continue;
1216 // if (reconstructor->HasLocalReconstruction()) continue;
1218 // AliCodeTimerStart(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1219 // AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1221 // AliCodeTimerStart(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1222 // AliInfo(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1224 // man->SetCacheFlag(kTRUE);
1225 // TString calibPath = Form("%s/Calib/*", fgkDetectorName[iDet]);
1226 // man->GetAll(calibPath); // entries are cached!
1228 // AliCodeTimerStop(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1230 // if (fRawReader) {
1231 // fRawReader->RewindEvents();
1232 // reconstructor->Reconstruct(fRunLoader, fRawReader);
1234 // reconstructor->Reconstruct(fRunLoader);
1237 // AliCodeTimerStop(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1238 // AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr));
1240 // // unload calibration data
1241 // man->UnloadFromCache(calibPath);
1242 // //man->ClearCache();
1245 // man->SetCacheFlag(origCache);
1247 // if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1248 // AliError(Form("the following detectors were not found: %s",
1250 // if (fStopOnError) return kFALSE;
1257 //_____________________________________________________________________________
1258 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1260 // run the local reconstruction
1262 static Int_t eventNr=0;
1263 AliCodeTimerAuto("")
1265 TString detStr = detectors;
1266 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1267 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1268 AliReconstructor* reconstructor = GetReconstructor(iDet);
1269 if (!reconstructor) continue;
1270 AliLoader* loader = fLoader[iDet];
1272 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1275 // conversion of digits
1276 if (fRawReader && reconstructor->HasDigitConversion()) {
1277 AliInfo(Form("converting raw data digits into root objects for %s",
1278 fgkDetectorName[iDet]));
1279 AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1280 fgkDetectorName[iDet]));
1281 loader->LoadDigits("update");
1282 loader->CleanDigits();
1283 loader->MakeDigitsContainer();
1284 TTree* digitsTree = loader->TreeD();
1285 reconstructor->ConvertDigits(fRawReader, digitsTree);
1286 loader->WriteDigits("OVERWRITE");
1287 loader->UnloadDigits();
1289 // local reconstruction
1290 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1291 AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1292 loader->LoadRecPoints("update");
1293 loader->CleanRecPoints();
1294 loader->MakeRecPointsContainer();
1295 TTree* clustersTree = loader->TreeR();
1296 if (fRawReader && !reconstructor->HasDigitConversion()) {
1297 reconstructor->Reconstruct(fRawReader, clustersTree);
1299 loader->LoadDigits("read");
1300 TTree* digitsTree = loader->TreeD();
1302 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1303 if (fStopOnError) return kFALSE;
1305 reconstructor->Reconstruct(digitsTree, clustersTree);
1307 loader->UnloadDigits();
1310 // In-loop QA for local reconstrucion
1311 if (fRunQA && fInLoopQA) {
1312 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1315 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1317 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1319 qadm->Exec(AliQA::kRECPOINTS, clustersTree) ;
1322 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1326 loader->WriteRecPoints("OVERWRITE");
1327 loader->UnloadRecPoints();
1328 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1331 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1332 AliError(Form("the following detectors were not found: %s",
1334 if (fStopOnError) return kFALSE;
1340 //_____________________________________________________________________________
1341 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1343 // run the barrel tracking
1345 AliCodeTimerAuto("")
1347 AliESDVertex* vertex = NULL;
1348 Double_t vtxPos[3] = {0, 0, 0};
1349 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1350 TArrayF mcVertex(3);
1351 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1352 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1353 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1357 if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1358 AliInfo("running the ITS vertex finder");
1359 if (fLoader[0]) fLoader[0]->LoadRecPoints();
1360 vertex = fVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber());
1361 if (fLoader[0]) fLoader[0]->UnloadRecPoints();
1363 AliWarning("Vertex not found");
1364 vertex = new AliESDVertex();
1365 vertex->SetName("default");
1368 vertex->SetName("reconstructed");
1372 AliInfo("getting the primary vertex from MC");
1373 vertex = new AliESDVertex(vtxPos, vtxErr);
1377 vertex->GetXYZ(vtxPos);
1378 vertex->GetSigmaXYZ(vtxErr);
1380 AliWarning("no vertex reconstructed");
1381 vertex = new AliESDVertex(vtxPos, vtxErr);
1383 esd->SetPrimaryVertexSPD(vertex);
1384 // if SPD multiplicity has been determined, it is stored in the ESD
1385 AliMultiplicity *mult = fVertexer->GetMultiplicity();
1386 if(mult)esd->SetMultiplicity(mult);
1388 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1389 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1396 //_____________________________________________________________________________
1397 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
1399 // run the HLT barrel tracking
1401 AliCodeTimerAuto("")
1404 AliError("Missing runLoader!");
1408 AliInfo("running HLT tracking");
1410 // Get a pointer to the HLT reconstructor
1411 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
1412 if (!reconstructor) return kFALSE;
1415 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1416 TString detName = fgkDetectorName[iDet];
1417 AliDebug(1, Form("%s HLT tracking", detName.Data()));
1418 reconstructor->SetOption(detName.Data());
1419 AliTracker *tracker = reconstructor->CreateTracker();
1421 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
1422 if (fStopOnError) return kFALSE;
1426 Double_t vtxErr[3]={0.005,0.005,0.010};
1427 const AliESDVertex *vertex = esd->GetVertex();
1428 vertex->GetXYZ(vtxPos);
1429 tracker->SetVertex(vtxPos,vtxErr);
1431 fLoader[iDet]->LoadRecPoints("read");
1432 TTree* tree = fLoader[iDet]->TreeR();
1434 AliError(Form("Can't get the %s cluster tree", detName.Data()));
1437 tracker->LoadClusters(tree);
1439 if (tracker->Clusters2Tracks(esd) != 0) {
1440 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
1444 tracker->UnloadClusters();
1452 //_____________________________________________________________________________
1453 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
1455 // run the muon spectrometer tracking
1457 AliCodeTimerAuto("")
1460 AliError("Missing runLoader!");
1463 Int_t iDet = 7; // for MUON
1465 AliInfo("is running...");
1467 // Get a pointer to the MUON reconstructor
1468 AliReconstructor *reconstructor = GetReconstructor(iDet);
1469 if (!reconstructor) return kFALSE;
1472 TString detName = fgkDetectorName[iDet];
1473 AliDebug(1, Form("%s tracking", detName.Data()));
1474 AliTracker *tracker = reconstructor->CreateTracker();
1476 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1481 fLoader[iDet]->LoadRecPoints("read");
1483 tracker->LoadClusters(fLoader[iDet]->TreeR());
1485 Int_t rv = tracker->Clusters2Tracks(esd);
1489 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1493 fLoader[iDet]->UnloadRecPoints();
1495 tracker->UnloadClusters();
1503 //_____________________________________________________________________________
1504 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
1506 // run the barrel tracking
1507 static Int_t eventNr=0;
1508 AliCodeTimerAuto("")
1510 AliInfo("running tracking");
1512 //Fill the ESD with the T0 info (will be used by the TOF)
1513 if (fReconstructor[11] && fLoader[11]) {
1514 fLoader[11]->LoadRecPoints("READ");
1515 TTree *treeR = fLoader[11]->TreeR();
1516 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
1519 // pass 1: TPC + ITS inwards
1520 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1521 if (!fTracker[iDet]) continue;
1522 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
1525 fLoader[iDet]->LoadRecPoints("read");
1526 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
1527 TTree* tree = fLoader[iDet]->TreeR();
1529 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1532 fTracker[iDet]->LoadClusters(tree);
1533 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1535 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
1536 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1539 if (fCheckPointLevel > 1) {
1540 WriteESD(esd, Form("%s.tracking", fgkDetectorName[iDet]));
1542 // preliminary PID in TPC needed by the ITS tracker
1544 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1545 AliESDpid::MakePID(esd);
1547 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
1550 // pass 2: ALL backwards
1552 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1553 if (!fTracker[iDet]) continue;
1554 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
1557 if (iDet > 1) { // all except ITS, TPC
1559 fLoader[iDet]->LoadRecPoints("read");
1560 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
1561 tree = fLoader[iDet]->TreeR();
1563 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1566 fTracker[iDet]->LoadClusters(tree);
1567 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1571 if (fTracker[iDet]->PropagateBack(esd) != 0) {
1572 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
1575 if (fCheckPointLevel > 1) {
1576 WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
1580 if (iDet > 2) { // all except ITS, TPC, TRD
1581 fTracker[iDet]->UnloadClusters();
1582 fLoader[iDet]->UnloadRecPoints();
1584 // updated PID in TPC needed by the ITS tracker -MI
1586 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1587 AliESDpid::MakePID(esd);
1589 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1592 // write space-points to the ESD in case alignment data output
1594 if (fWriteAlignmentData)
1595 WriteAlignmentData(esd);
1597 // pass 3: TRD + TPC + ITS refit inwards
1599 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1601 for (Int_t iDet = 2; iDet >= 0; iDet--) {
1602 if (!fTracker[iDet]) continue;
1603 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
1606 if (fTracker[iDet]->RefitInward(esd) != 0) {
1607 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
1610 // run postprocessing
1611 if (fTracker[iDet]->PostProcess(esd) != 0) {
1612 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
1615 if (fCheckPointLevel > 1) {
1616 WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
1618 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1620 fTracker[iDet]->UnloadClusters();
1621 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
1622 fLoader[iDet]->UnloadRecPoints();
1623 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
1626 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1632 //_____________________________________________________________________________
1633 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
1635 // Remove the data which are not needed for the physics analysis.
1638 Int_t nTracks=esd->GetNumberOfTracks();
1639 Int_t nV0s=esd->GetNumberOfV0s();
1641 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
1643 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
1644 Bool_t rc=esd->Clean(cleanPars);
1646 nTracks=esd->GetNumberOfTracks();
1647 nV0s=esd->GetNumberOfV0s();
1649 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
1654 //_____________________________________________________________________________
1655 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
1657 // fill the event summary data
1659 AliCodeTimerAuto("")
1660 static Int_t eventNr=0;
1661 TString detStr = detectors;
1663 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1664 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1665 AliReconstructor* reconstructor = GetReconstructor(iDet);
1666 if (!reconstructor) continue;
1667 if (!ReadESD(esd, fgkDetectorName[iDet])) {
1668 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
1669 TTree* clustersTree = NULL;
1670 if (fLoader[iDet]) {
1671 fLoader[iDet]->LoadRecPoints("read");
1672 clustersTree = fLoader[iDet]->TreeR();
1673 if (!clustersTree) {
1674 AliError(Form("Can't get the %s clusters tree",
1675 fgkDetectorName[iDet]));
1676 if (fStopOnError) return kFALSE;
1679 if (fRawReader && !reconstructor->HasDigitConversion()) {
1680 reconstructor->FillESD(fRawReader, clustersTree, esd);
1682 TTree* digitsTree = NULL;
1683 if (fLoader[iDet]) {
1684 fLoader[iDet]->LoadDigits("read");
1685 digitsTree = fLoader[iDet]->TreeD();
1687 AliError(Form("Can't get the %s digits tree",
1688 fgkDetectorName[iDet]));
1689 if (fStopOnError) return kFALSE;
1692 reconstructor->FillESD(digitsTree, clustersTree, esd);
1693 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
1695 if (fLoader[iDet]) {
1696 fLoader[iDet]->UnloadRecPoints();
1699 if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
1703 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1704 AliError(Form("the following detectors were not found: %s",
1706 if (fStopOnError) return kFALSE;
1708 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
1713 //_____________________________________________________________________________
1714 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
1716 // Reads the trigger decision which is
1717 // stored in Trigger.root file and fills
1718 // the corresponding esd entries
1720 AliCodeTimerAuto("")
1722 AliInfo("Filling trigger information into the ESD");
1724 AliCentralTrigger *aCTP = NULL;
1727 AliCTPRawStream input(fRawReader);
1728 if (!input.Next()) {
1729 AliError("No valid CTP (trigger) DDL raw data is found ! The trigger information is not stored in the ESD !");
1732 esd->SetTriggerMask(input.GetClassMask());
1733 esd->SetTriggerCluster(input.GetClusterMask());
1735 aCTP = new AliCentralTrigger();
1736 TString configstr("");
1737 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
1738 AliError("No trigger configuration found in OCDB! The trigger classes information will no be stored in ESD!");
1743 AliRunLoader *runloader = AliRunLoader::GetRunLoader();
1745 if (!runloader->LoadTrigger()) {
1746 aCTP = runloader->GetTrigger();
1747 esd->SetTriggerMask(aCTP->GetClassMask());
1748 esd->SetTriggerCluster(aCTP->GetClusterMask());
1751 AliWarning("No trigger can be loaded! The trigger information is not stored in the ESD !");
1756 AliError("No run loader is available! The trigger information is not stored in the ESD !");
1761 // Now fill the trigger class names into AliESDRun object
1762 AliTriggerConfiguration *config = aCTP->GetConfiguration();
1764 AliError("No trigger configuration has been found! The trigger classes information will no be stored in ESD!");
1768 const TObjArray& classesArray = config->GetClasses();
1769 Int_t nclasses = classesArray.GetEntriesFast();
1770 for( Int_t j=0; j<nclasses; j++ ) {
1771 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At( j );
1772 Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
1773 esd->SetTriggerClass(trclass->GetName(),trindex);
1783 //_____________________________________________________________________________
1784 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
1787 // Filling information from RawReader Header
1790 AliInfo("Filling information from RawReader Header");
1791 esd->SetBunchCrossNumber(0);
1792 esd->SetOrbitNumber(0);
1793 esd->SetPeriodNumber(0);
1794 esd->SetTimeStamp(0);
1795 esd->SetEventType(0);
1796 const AliRawEventHeaderBase * eventHeader = fRawReader->GetEventHeader();
1799 const UInt_t *id = eventHeader->GetP("Id");
1800 esd->SetBunchCrossNumber((id)[1]&0x00000fff);
1801 esd->SetOrbitNumber((((id)[0]<<20)&0xf00000)|(((id)[1]>>12)&0xfffff));
1802 esd->SetPeriodNumber(((id)[0]>>4)&0x0fffffff);
1804 esd->SetTimeStamp((eventHeader->Get("Timestamp")));
1805 esd->SetEventType((eventHeader->Get("Type")));
1812 //_____________________________________________________________________________
1813 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
1815 // check whether detName is contained in detectors
1816 // if yes, it is removed from detectors
1818 // check if all detectors are selected
1819 if ((detectors.CompareTo("ALL") == 0) ||
1820 detectors.BeginsWith("ALL ") ||
1821 detectors.EndsWith(" ALL") ||
1822 detectors.Contains(" ALL ")) {
1827 // search for the given detector
1828 Bool_t result = kFALSE;
1829 if ((detectors.CompareTo(detName) == 0) ||
1830 detectors.BeginsWith(detName+" ") ||
1831 detectors.EndsWith(" "+detName) ||
1832 detectors.Contains(" "+detName+" ")) {
1833 detectors.ReplaceAll(detName, "");
1837 // clean up the detectors string
1838 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
1839 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
1840 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
1845 //_____________________________________________________________________________
1846 Bool_t AliReconstruction::InitRunLoader()
1848 // get or create the run loader
1850 if (gAlice) delete gAlice;
1853 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
1854 // load all base libraries to get the loader classes
1855 TString libs = gSystem->GetLibraries();
1856 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1857 TString detName = fgkDetectorName[iDet];
1858 if (detName == "HLT") continue;
1859 if (libs.Contains("lib" + detName + "base.so")) continue;
1860 gSystem->Load("lib" + detName + "base.so");
1862 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
1864 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
1868 fRunLoader->CdGAFile();
1869 if (gFile->GetKey(AliRunLoader::GetGAliceName())) {
1870 if (fRunLoader->LoadgAlice() == 0) {
1871 gAlice = fRunLoader->GetAliRun();
1872 AliTracker::SetFieldMap(gAlice->Field(),fUniformField);
1875 if (!gAlice && !fRawReader) {
1876 AliError(Form("no gAlice object found in file %s",
1877 fGAliceFileName.Data()));
1882 //PH This is a temporary fix to give access to the kinematics
1883 //PH that is needed for the labels of ITS clusters
1884 fRunLoader->LoadHeader();
1885 fRunLoader->LoadKinematics();
1887 } else { // galice.root does not exist
1889 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
1893 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
1894 AliConfig::GetDefaultEventFolderName(),
1897 AliError(Form("could not create run loader in file %s",
1898 fGAliceFileName.Data()));
1902 fRunLoader->MakeTree("E");
1904 while (fRawReader->NextEvent()) {
1905 fRunLoader->SetEventNumber(iEvent);
1906 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1908 fRunLoader->MakeTree("H");
1909 fRunLoader->TreeE()->Fill();
1912 fRawReader->RewindEvents();
1913 if (fNumberOfEventsPerFile > 0)
1914 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
1916 fRunLoader->SetNumberOfEventsPerFile(iEvent);
1917 fRunLoader->WriteHeader("OVERWRITE");
1918 fRunLoader->CdGAFile();
1919 fRunLoader->Write(0, TObject::kOverwrite);
1920 // AliTracker::SetFieldMap(???);
1926 //_____________________________________________________________________________
1927 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
1929 // get the reconstructor object and the loader for a detector
1931 if (fReconstructor[iDet]) return fReconstructor[iDet];
1933 // load the reconstructor object
1934 TPluginManager* pluginManager = gROOT->GetPluginManager();
1935 TString detName = fgkDetectorName[iDet];
1936 TString recName = "Ali" + detName + "Reconstructor";
1937 if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT")) return NULL;
1939 AliReconstructor* reconstructor = NULL;
1940 // first check if a plugin is defined for the reconstructor
1941 TPluginHandler* pluginHandler =
1942 pluginManager->FindHandler("AliReconstructor", detName);
1943 // if not, add a plugin for it
1944 if (!pluginHandler) {
1945 AliDebug(1, Form("defining plugin for %s", recName.Data()));
1946 TString libs = gSystem->GetLibraries();
1947 if (libs.Contains("lib" + detName + "base.so") ||
1948 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
1949 pluginManager->AddHandler("AliReconstructor", detName,
1950 recName, detName + "rec", recName + "()");
1952 pluginManager->AddHandler("AliReconstructor", detName,
1953 recName, detName, recName + "()");
1955 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
1957 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
1958 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
1960 if (reconstructor) {
1961 TObject* obj = fOptions.FindObject(detName.Data());
1962 if (obj) reconstructor->SetOption(obj->GetTitle());
1963 reconstructor->Init();
1964 fReconstructor[iDet] = reconstructor;
1967 // get or create the loader
1968 if (detName != "HLT") {
1969 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
1970 if (!fLoader[iDet]) {
1971 AliConfig::Instance()
1972 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
1974 // first check if a plugin is defined for the loader
1976 pluginManager->FindHandler("AliLoader", detName);
1977 // if not, add a plugin for it
1978 if (!pluginHandler) {
1979 TString loaderName = "Ali" + detName + "Loader";
1980 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
1981 pluginManager->AddHandler("AliLoader", detName,
1982 loaderName, detName + "base",
1983 loaderName + "(const char*, TFolder*)");
1984 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
1986 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
1988 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
1989 fRunLoader->GetEventFolder());
1991 if (!fLoader[iDet]) { // use default loader
1992 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
1994 if (!fLoader[iDet]) {
1995 AliWarning(Form("couldn't get loader for %s", detName.Data()));
1996 if (fStopOnError) return NULL;
1998 fRunLoader->AddLoader(fLoader[iDet]);
1999 fRunLoader->CdGAFile();
2000 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2001 fRunLoader->Write(0, TObject::kOverwrite);
2006 return reconstructor;
2009 //_____________________________________________________________________________
2010 Bool_t AliReconstruction::CreateVertexer()
2012 // create the vertexer
2015 AliReconstructor* itsReconstructor = GetReconstructor(0);
2016 if (itsReconstructor) {
2017 fVertexer = itsReconstructor->CreateVertexer();
2020 AliWarning("couldn't create a vertexer for ITS");
2021 if (fStopOnError) return kFALSE;
2027 //_____________________________________________________________________________
2028 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2030 // create the trackers
2032 TString detStr = detectors;
2033 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2034 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2035 AliReconstructor* reconstructor = GetReconstructor(iDet);
2036 if (!reconstructor) continue;
2037 TString detName = fgkDetectorName[iDet];
2038 if (detName == "HLT") {
2039 fRunHLTTracking = kTRUE;
2042 if (detName == "MUON") {
2043 fRunMuonTracking = kTRUE;
2048 fTracker[iDet] = reconstructor->CreateTracker();
2049 if (!fTracker[iDet] && (iDet < 7)) {
2050 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2051 if (fStopOnError) return kFALSE;
2053 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2059 //_____________________________________________________________________________
2060 void AliReconstruction::CleanUp(TFile* file, TFile* fileOld)
2062 // delete trackers and the run loader and close and delete the file
2064 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2065 delete fReconstructor[iDet];
2066 fReconstructor[iDet] = NULL;
2067 fLoader[iDet] = NULL;
2068 delete fTracker[iDet];
2069 fTracker[iDet] = NULL;
2070 // delete fQADataMaker[iDet];
2071 // fQADataMaker[iDet] = NULL;
2076 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2077 delete fDiamondProfile;
2078 fDiamondProfile = NULL;
2088 if (fParentRawReader) delete fParentRawReader;
2089 fParentRawReader=NULL;
2099 gSystem->Unlink("AliESDs.old.root");
2103 //_____________________________________________________________________________
2105 Bool_t AliReconstruction::ReadESD(AliESDEvent*& esd, const char* recStep) const
2107 // read the ESD event from a file
2109 if (!esd) return kFALSE;
2111 sprintf(fileName, "ESD_%d.%d_%s.root",
2112 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2113 if (gSystem->AccessPathName(fileName)) return kFALSE;
2115 AliInfo(Form("reading ESD from file %s", fileName));
2116 AliDebug(1, Form("reading ESD from file %s", fileName));
2117 TFile* file = TFile::Open(fileName);
2118 if (!file || !file->IsOpen()) {
2119 AliError(Form("opening %s failed", fileName));
2126 esd = (AliESDEvent*) file->Get("ESD");
2135 //_____________________________________________________________________________
2136 void AliReconstruction::WriteESD(AliESDEvent* esd, const char* recStep) const
2138 // write the ESD event to a file
2142 sprintf(fileName, "ESD_%d.%d_%s.root",
2143 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2145 AliDebug(1, Form("writing ESD to file %s", fileName));
2146 TFile* file = TFile::Open(fileName, "recreate");
2147 if (!file || !file->IsOpen()) {
2148 AliError(Form("opening %s failed", fileName));
2160 //_____________________________________________________________________________
2161 void AliReconstruction::ESDFile2AODFile(TFile* esdFile, TFile* aodFile)
2163 // write all files from the given esd file to an aod file
2165 // create an AliAOD object
2166 AliAODEvent *aod = new AliAODEvent();
2167 aod->CreateStdContent();
2173 TTree *aodTree = new TTree("aodTree", "AliAOD tree");
2174 aodTree->Branch(aod->GetList());
2177 TTree *t = (TTree*) esdFile->Get("esdTree");
2178 AliESDEvent *esd = new AliESDEvent();
2179 esd->ReadFromTree(t);
2181 Int_t nEvents = t->GetEntries();
2183 // set arrays and pointers
2193 // loop over events and fill them
2194 for (Int_t iEvent = 0; iEvent < nEvents; ++iEvent) {
2195 //cout << "event: " << iEvent << endl;
2196 t->GetEntry(iEvent);
2198 // Multiplicity information needed by the header (to be revised!)
2199 Int_t nTracks = esd->GetNumberOfTracks();
2200 Int_t nPosTracks = 0;
2201 for (Int_t iTrack=0; iTrack<nTracks; ++iTrack)
2202 if (esd->GetTrack(iTrack)->Charge()> 0) nPosTracks++;
2204 // Access the header
2205 AliAODHeader *header = aod->GetHeader();
2208 header->SetRunNumber (esd->GetRunNumber() );
2209 header->SetBunchCrossNumber(esd->GetBunchCrossNumber());
2210 header->SetOrbitNumber (esd->GetOrbitNumber() );
2211 header->SetPeriodNumber (esd->GetPeriodNumber() );
2212 header->SetTriggerMask (esd->GetTriggerMask() );
2213 header->SetTriggerCluster (esd->GetTriggerCluster() );
2214 header->SetEventType (esd->GetEventType() );
2215 header->SetMagneticField (esd->GetMagneticField() );
2216 header->SetZDCN1Energy (esd->GetZDCN1Energy() );
2217 header->SetZDCP1Energy (esd->GetZDCP1Energy() );
2218 header->SetZDCN2Energy (esd->GetZDCN2Energy() );
2219 header->SetZDCP2Energy (esd->GetZDCP2Energy() );
2220 header->SetZDCEMEnergy (esd->GetZDCEMEnergy(0),esd->GetZDCEMEnergy(1));
2221 header->SetRefMultiplicity (nTracks);
2222 header->SetRefMultiplicityPos(nPosTracks);
2223 header->SetRefMultiplicityNeg(nTracks - nPosTracks);
2224 header->SetMuonMagFieldScale(-999.); // FIXME
2225 header->SetCentrality(-999.); // FIXME
2227 Int_t nV0s = esd->GetNumberOfV0s();
2228 Int_t nCascades = esd->GetNumberOfCascades();
2229 Int_t nKinks = esd->GetNumberOfKinks();
2230 Int_t nVertices = nV0s + 2*nCascades /*could lead to two vertices, one V0 and the Xi */+ nKinks + 1 /* = prim. vtx*/;
2232 Int_t nCaloClus = esd->GetNumberOfCaloClusters();
2234 Int_t nPmdClus = esd->GetNumberOfPmdTracks();
2236 aod->ResetStd(nTracks, nVertices, nV0s+nCascades, nJets, nCaloClus, nFmdClus, nPmdClus);
2238 // Array to take into account the tracks already added to the AOD
2239 Bool_t * usedTrack = NULL;
2241 usedTrack = new Bool_t[nTracks];
2242 for (Int_t iTrack=0; iTrack<nTracks; ++iTrack) usedTrack[iTrack]=kFALSE;
2244 // Array to take into account the V0s already added to the AOD
2245 Bool_t * usedV0 = NULL;
2247 usedV0 = new Bool_t[nV0s];
2248 for (Int_t iV0=0; iV0<nV0s; ++iV0) usedV0[iV0]=kFALSE;
2250 // Array to take into account the kinks already added to the AOD
2251 Bool_t * usedKink = NULL;
2253 usedKink = new Bool_t[nKinks];
2254 for (Int_t iKink=0; iKink<nKinks; ++iKink) usedKink[iKink]=kFALSE;
2257 // Access to the AOD container of vertices
2258 TClonesArray &vertices = *(aod->GetVertices());
2261 // Access to the AOD container of tracks
2262 TClonesArray &tracks = *(aod->GetTracks());
2265 // Access to the AOD container of V0s
2266 TClonesArray &V0s = *(aod->GetV0s());
2269 // Add primary vertex. The primary tracks will be defined
2270 // after the loops on the composite objects (V0, cascades, kinks)
2271 const AliESDVertex *vtx = esd->GetPrimaryVertex();
2273 vtx->GetXYZ(pos); // position
2274 vtx->GetCovMatrix(covVtx); //covariance matrix
2276 AliAODVertex * primary = new(vertices[jVertices++])
2277 AliAODVertex(pos, covVtx, vtx->GetChi2toNDF(), NULL, -1, AliAODVertex::kPrimary);
2280 AliAODTrack *aodTrack = 0x0;
2282 // Create vertices starting from the most complex objects
2285 for (Int_t nCascade = 0; nCascade < nCascades; ++nCascade) {
2286 AliESDcascade *cascade = esd->GetCascade(nCascade);
2288 cascade->GetXYZ(pos[0], pos[1], pos[2]);
2289 cascade->GetPosCovXi(covVtx);
2291 // Add the cascade vertex
2292 AliAODVertex * vcascade = new(vertices[jVertices++]) AliAODVertex(pos,
2294 cascade->GetChi2Xi(), // = chi2/NDF since NDF = 2*2-3
2297 AliAODVertex::kCascade);
2299 primary->AddDaughter(vcascade); // the cascade 'particle' (represented by a vertex) is added as a daughter to the primary vertex
2301 // Add the V0 from the cascade. The ESD class have to be optimized...
2302 // Now we have to search for the corresponding V0 in the list of V0s
2303 // using the indeces of the positive and negative tracks
2305 Int_t posFromV0 = cascade->GetPindex();
2306 Int_t negFromV0 = cascade->GetNindex();
2309 AliESDv0 * v0 = 0x0;
2312 for (Int_t iV0=0; iV0<nV0s; ++iV0) {
2314 v0 = esd->GetV0(iV0);
2315 Int_t posV0 = v0->GetPindex();
2316 Int_t negV0 = v0->GetNindex();
2318 if (posV0==posFromV0 && negV0==negFromV0) {
2324 AliAODVertex * vV0FromCascade = 0x0;
2326 if (indV0>-1 && !usedV0[indV0]) {
2328 // the V0 exists in the array of V0s and is not used
2330 usedV0[indV0] = kTRUE;
2332 v0->GetXYZ(pos[0], pos[1], pos[2]);
2333 v0->GetPosCov(covVtx);
2335 vV0FromCascade = new(vertices[jVertices++]) AliAODVertex(pos,
2337 v0->GetChi2V0(), // = chi2/NDF since NDF = 2*2-3
2343 // the V0 doesn't exist in the array of V0s or was used
2344 cerr << "Error: event " << iEvent << " cascade " << nCascade
2345 << " The V0 " << indV0
2346 << " doesn't exist in the array of V0s or was used!" << endl;
2348 cascade->GetXYZ(pos[0], pos[1], pos[2]);
2349 cascade->GetPosCov(covVtx);
2351 vV0FromCascade = new(vertices[jVertices++]) AliAODVertex(pos,
2353 v0->GetChi2V0(), // = chi2/NDF since NDF = 2*2-3
2357 vcascade->AddDaughter(vV0FromCascade);
2361 // Add the positive tracks from the V0
2363 if (! usedTrack[posFromV0]) {
2365 usedTrack[posFromV0] = kTRUE;
2367 AliESDtrack *esdTrack = esd->GetTrack(posFromV0);
2368 esdTrack->GetPxPyPz(p_pos);
2369 esdTrack->GetXYZ(pos);
2370 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2371 esdTrack->GetESDpid(pid);
2373 vV0FromCascade->AddDaughter(aodTrack =
2374 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2375 esdTrack->GetLabel(),
2381 (Short_t)esdTrack->Charge(),
2382 esdTrack->GetITSClusterMap(),
2385 kTRUE, // check if this is right
2386 kFALSE, // check if this is right
2387 AliAODTrack::kSecondary)
2389 aodTrack->ConvertAliPIDtoAODPID();
2392 cerr << "Error: event " << iEvent << " cascade " << nCascade
2393 << " track " << posFromV0 << " has already been used!" << endl;
2396 // Add the negative tracks from the V0
2398 if (!usedTrack[negFromV0]) {
2400 usedTrack[negFromV0] = kTRUE;
2402 AliESDtrack *esdTrack = esd->GetTrack(negFromV0);
2403 esdTrack->GetPxPyPz(p_neg);
2404 esdTrack->GetXYZ(pos);
2405 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2406 esdTrack->GetESDpid(pid);
2408 vV0FromCascade->AddDaughter(aodTrack =
2409 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2410 esdTrack->GetLabel(),
2416 (Short_t)esdTrack->Charge(),
2417 esdTrack->GetITSClusterMap(),
2420 kTRUE, // check if this is right
2421 kFALSE, // check if this is right
2422 AliAODTrack::kSecondary)
2424 aodTrack->ConvertAliPIDtoAODPID();
2427 cerr << "Error: event " << iEvent << " cascade " << nCascade
2428 << " track " << negFromV0 << " has already been used!" << endl;
2431 // add it to the V0 array as well
2432 Double_t d0[2] = { -999., -99.};
2433 // counting is probably wrong
2434 new(V0s[jV0s++]) AliAODv0(vV0FromCascade, -999., -99., p_pos, p_neg, d0); // to be refined
2436 // Add the bachelor track from the cascade
2438 Int_t bachelor = cascade->GetBindex();
2440 if(!usedTrack[bachelor]) {
2442 usedTrack[bachelor] = kTRUE;
2444 AliESDtrack *esdTrack = esd->GetTrack(bachelor);
2445 esdTrack->GetPxPyPz(p);
2446 esdTrack->GetXYZ(pos);
2447 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2448 esdTrack->GetESDpid(pid);
2450 vcascade->AddDaughter(aodTrack =
2451 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2452 esdTrack->GetLabel(),
2458 (Short_t)esdTrack->Charge(),
2459 esdTrack->GetITSClusterMap(),
2462 kTRUE, // check if this is right
2463 kFALSE, // check if this is right
2464 AliAODTrack::kSecondary)
2466 aodTrack->ConvertAliPIDtoAODPID();
2469 cerr << "Error: event " << iEvent << " cascade " << nCascade
2470 << " track " << bachelor << " has already been used!" << endl;
2473 // Add the primary track of the cascade (if any)
2475 } // end of the loop on cascades
2479 for (Int_t nV0 = 0; nV0 < nV0s; ++nV0) {
2481 if (usedV0[nV0]) continue; // skip if aready added to the AOD
2483 AliESDv0 *v0 = esd->GetV0(nV0);
2485 v0->GetXYZ(pos[0], pos[1], pos[2]);
2486 v0->GetPosCov(covVtx);
2488 AliAODVertex * vV0 =
2489 new(vertices[jVertices++]) AliAODVertex(pos,
2491 v0->GetChi2V0(), // = chi2/NDF since NDF = 2*2-3
2495 primary->AddDaughter(vV0);
2497 Int_t posFromV0 = v0->GetPindex();
2498 Int_t negFromV0 = v0->GetNindex();
2500 // Add the positive tracks from the V0
2502 if (!usedTrack[posFromV0]) {
2504 usedTrack[posFromV0] = kTRUE;
2506 AliESDtrack *esdTrack = esd->GetTrack(posFromV0);
2507 esdTrack->GetPxPyPz(p_pos);
2508 esdTrack->GetXYZ(pos);
2509 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2510 esdTrack->GetESDpid(pid);
2512 vV0->AddDaughter(aodTrack =
2513 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2514 esdTrack->GetLabel(),
2520 (Short_t)esdTrack->Charge(),
2521 esdTrack->GetITSClusterMap(),
2524 kTRUE, // check if this is right
2525 kFALSE, // check if this is right
2526 AliAODTrack::kSecondary)
2528 aodTrack->ConvertAliPIDtoAODPID();
2531 cerr << "Error: event " << iEvent << " V0 " << nV0
2532 << " track " << posFromV0 << " has already been used!" << endl;
2535 // Add the negative tracks from the V0
2537 if (!usedTrack[negFromV0]) {
2539 usedTrack[negFromV0] = kTRUE;
2541 AliESDtrack *esdTrack = esd->GetTrack(negFromV0);
2542 esdTrack->GetPxPyPz(p_neg);
2543 esdTrack->GetXYZ(pos);
2544 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2545 esdTrack->GetESDpid(pid);
2547 vV0->AddDaughter(aodTrack =
2548 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2549 esdTrack->GetLabel(),
2555 (Short_t)esdTrack->Charge(),
2556 esdTrack->GetITSClusterMap(),
2559 kTRUE, // check if this is right
2560 kFALSE, // check if this is right
2561 AliAODTrack::kSecondary)
2563 aodTrack->ConvertAliPIDtoAODPID();
2566 cerr << "Error: event " << iEvent << " V0 " << nV0
2567 << " track " << negFromV0 << " has already been used!" << endl;
2570 // add it to the V0 array as well
2571 Double_t d0[2] = { 999., 99.};
2572 new(V0s[jV0s++]) AliAODv0(vV0, 999., 99., p_pos, p_neg, d0); // to be refined
2575 // end of the loop on V0s
2577 // Kinks: it is a big mess the access to the information in the kinks
2578 // The loop is on the tracks in order to find the mother and daugther of each kink
2581 for (Int_t iTrack=0; iTrack<nTracks; ++iTrack) {
2583 AliESDtrack * esdTrack = esd->GetTrack(iTrack);
2585 Int_t ikink = esdTrack->GetKinkIndex(0);
2588 // Negative kink index: mother, positive: daughter
2590 // Search for the second track of the kink
2592 for (Int_t jTrack = iTrack+1; jTrack<nTracks; ++jTrack) {
2594 AliESDtrack * esdTrack1 = esd->GetTrack(jTrack);
2596 Int_t jkink = esdTrack1->GetKinkIndex(0);
2598 if ( TMath::Abs(ikink)==TMath::Abs(jkink) ) {
2600 // The two tracks are from the same kink
2602 if (usedKink[TMath::Abs(ikink)-1]) continue; // skip used kinks
2605 Int_t idaughter = -1;
2607 if (ikink<0 && jkink>0) {
2612 else if (ikink>0 && jkink<0) {
2618 cerr << "Error: Wrong combination of kink indexes: "
2619 << ikink << " " << jkink << endl;
2623 // Add the mother track
2625 AliAODTrack * mother = NULL;
2627 if (!usedTrack[imother]) {
2629 usedTrack[imother] = kTRUE;
2631 AliESDtrack *esdTrack = esd->GetTrack(imother);
2632 esdTrack->GetPxPyPz(p);
2633 esdTrack->GetXYZ(pos);
2634 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2635 esdTrack->GetESDpid(pid);
2638 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2639 esdTrack->GetLabel(),
2645 (Short_t)esdTrack->Charge(),
2646 esdTrack->GetITSClusterMap(),
2649 kTRUE, // check if this is right
2650 kTRUE, // check if this is right
2651 AliAODTrack::kPrimary);
2652 primary->AddDaughter(mother);
2653 mother->ConvertAliPIDtoAODPID();
2656 cerr << "Error: event " << iEvent << " kink " << TMath::Abs(ikink)-1
2657 << " track " << imother << " has already been used!" << endl;
2660 // Add the kink vertex
2661 AliESDkink * kink = esd->GetKink(TMath::Abs(ikink)-1);
2663 AliAODVertex * vkink =
2664 new(vertices[jVertices++]) AliAODVertex(kink->GetPosition(),
2668 esdTrack->GetID(), // This is the track ID of the mother's track!
2669 AliAODVertex::kKink);
2670 // Add the daughter track
2672 AliAODTrack * daughter = NULL;
2674 if (!usedTrack[idaughter]) {
2676 usedTrack[idaughter] = kTRUE;
2678 AliESDtrack *esdTrack = esd->GetTrack(idaughter);
2679 esdTrack->GetPxPyPz(p);
2680 esdTrack->GetXYZ(pos);
2681 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2682 esdTrack->GetESDpid(pid);
2685 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2686 esdTrack->GetLabel(),
2692 (Short_t)esdTrack->Charge(),
2693 esdTrack->GetITSClusterMap(),
2696 kTRUE, // check if this is right
2697 kTRUE, // check if this is right
2698 AliAODTrack::kPrimary);
2699 vkink->AddDaughter(daughter);
2700 daughter->ConvertAliPIDtoAODPID();
2703 cerr << "Error: event " << iEvent << " kink " << TMath::Abs(ikink)-1
2704 << " track " << idaughter << " has already been used!" << endl;
2710 vertices.Expand(jVertices);
2712 // Tracks (primary and orphan)
2713 for (Int_t nTrack = 0; nTrack < nTracks; ++nTrack) {
2715 if (usedTrack[nTrack]) continue;
2717 AliESDtrack *esdTrack = esd->GetTrack(nTrack);
2718 esdTrack->GetPxPyPz(p);
2719 esdTrack->GetXYZ(pos);
2720 esdTrack->GetCovarianceXYZPxPyPz(covTr);
2721 esdTrack->GetESDpid(pid);
2723 Float_t impactXY, impactZ;
2725 esdTrack->GetImpactParameters(impactXY,impactZ);
2728 // track inside the beam pipe
2730 primary->AddDaughter(aodTrack =
2731 new(tracks[jTracks++]) AliAODTrack(esdTrack->GetID(),
2732 esdTrack->GetLabel(),
2738 (Short_t)esdTrack->Charge(),
2739 esdTrack->GetITSClusterMap(),
2742 kTRUE, // check if this is right
2743 kTRUE, // check if this is right
2744 AliAODTrack::kPrimary)
2746 aodTrack->ConvertAliPIDtoAODPID();
2749 // outside the beam pipe: orphan track
2750 // Don't write them anymore!
2753 } // end of loop on tracks
2756 Int_t nMuTracks = esd->GetNumberOfMuonTracks();
2757 for (Int_t nMuTrack = 0; nMuTrack < nMuTracks; ++nMuTrack) {
2759 AliESDMuonTrack *esdMuTrack = esd->GetMuonTrack(nMuTrack);
2760 p[0] = esdMuTrack->Px();
2761 p[1] = esdMuTrack->Py();
2762 p[2] = esdMuTrack->Pz();
2763 pos[0] = primary->GetX();
2764 pos[1] = primary->GetY();
2765 pos[2] = primary->GetZ();
2767 // has to be changed once the muon pid is provided by the ESD
2768 for (Int_t i = 0; i < 10; pid[i++] = 0.); pid[AliAODTrack::kMuon]=1.;
2770 primary->AddDaughter(aodTrack =
2771 new(tracks[jTracks++]) AliAODTrack(0, // no ID provided
2772 0, // no label provided
2777 NULL, // no covariance matrix provided
2778 esdMuTrack->Charge(),
2779 0, // ITSClusterMap is set below
2782 kFALSE, // muon tracks are not used to fit the primary vtx
2783 kFALSE, // not used for vertex fit
2784 AliAODTrack::kPrimary)
2787 aodTrack->SetHitsPatternInTrigCh(esdMuTrack->GetHitsPatternInTrigCh());
2788 Int_t track2Trigger = esdMuTrack->GetMatchTrigger();
2789 aodTrack->SetMatchTrigger(track2Trigger);
2791 aodTrack->SetChi2MatchTrigger(esdMuTrack->GetChi2MatchTrigger());
2793 aodTrack->SetChi2MatchTrigger(0.);
2795 tracks.Expand(jTracks); // remove 'empty slots' due to unwritten tracks
2797 // Access to the AOD container of PMD clusters
2798 TClonesArray &pmdClusters = *(aod->GetPmdClusters());
2799 Int_t jPmdClusters=0;
2801 for (Int_t iPmd = 0; iPmd < nPmdClus; ++iPmd) {
2802 // file pmd clusters, to be revised!
2803 AliESDPmdTrack *pmdTrack = esd->GetPmdTrack(iPmd);
2806 Double_t pos[3] = { pmdTrack->GetClusterX(), pmdTrack->GetClusterY(), pmdTrack->GetClusterZ() };
2807 Double_t pid[9] = { 0., 0., 0., 0., 0., 0., 0., 0., 0. }; // to be revised!
2809 // assoc cluster not set
2810 new(pmdClusters[jPmdClusters++]) AliAODPmdCluster(iPmd, nLabel, label, pmdTrack->GetClusterADC(), pos, pid);
2813 // Access to the AOD container of clusters
2814 TClonesArray &caloClusters = *(aod->GetCaloClusters());
2817 for (Int_t iClust=0; iClust<nCaloClus; ++iClust) {
2819 AliESDCaloCluster * cluster = esd->GetCaloCluster(iClust);
2821 Int_t id = cluster->GetID();
2824 Float_t energy = cluster->E();
2825 cluster->GetPosition(posF);
2826 Char_t ttype=AliAODCluster::kUndef;
2828 if (cluster->GetClusterType() == AliESDCaloCluster::kPHOSCluster) {
2829 ttype=AliAODCluster::kPHOSNeutral;
2831 else if (cluster->GetClusterType() == AliESDCaloCluster::kEMCALClusterv1) {
2832 ttype = AliAODCluster::kEMCALClusterv1;
2836 AliAODCaloCluster *caloCluster = new(caloClusters[jClusters++]) AliAODCaloCluster(id,
2844 caloCluster->SetCaloCluster(); // to be refined!
2847 caloClusters.Expand(jClusters); // resize TObjArray to 'remove' slots for pseudo clusters
2848 // end of loop on calo clusters
2850 // fill EMCAL cell info
2851 if (esd->GetEMCALCells()) { // protection against missing ESD information
2852 AliESDCaloCells &esdEMcells = *(esd->GetEMCALCells());
2853 Int_t nEMcell = esdEMcells.GetNumberOfCells() ;
2855 AliAODCaloCells &aodEMcells = *(aod->GetEMCALCells());
2856 aodEMcells.CreateContainer(nEMcell);
2857 aodEMcells.SetType(AliAODCaloCells::kEMCAL);
2858 for (Int_t iCell = 0; iCell < nEMcell; iCell++) {
2859 aodEMcells.SetCell(iCell,esdEMcells.GetCellNumber(iCell),esdEMcells.GetAmplitude(iCell));
2864 // fill PHOS cell info
2865 if (esd->GetPHOSCells()) { // protection against missing ESD information
2866 AliESDCaloCells &esdPHcells = *(esd->GetPHOSCells());
2867 Int_t nPHcell = esdPHcells.GetNumberOfCells() ;
2869 AliAODCaloCells &aodPHcells = *(aod->GetPHOSCells());
2870 aodPHcells.CreateContainer(nPHcell);
2871 aodPHcells.SetType(AliAODCaloCells::kPHOS);
2872 for (Int_t iCell = 0; iCell < nPHcell; iCell++) {
2873 aodPHcells.SetCell(iCell,esdPHcells.GetCellNumber(iCell),esdPHcells.GetAmplitude(iCell));
2879 AliAODTracklets &SPDTracklets = *(aod->GetTracklets());
2880 const AliMultiplicity *mult = esd->GetMultiplicity();
2882 if (mult->GetNumberOfTracklets()>0) {
2883 SPDTracklets.CreateContainer(mult->GetNumberOfTracklets());
2885 for (Int_t n=0; n<mult->GetNumberOfTracklets(); n++) {
2886 SPDTracklets.SetTracklet(n, mult->GetTheta(n), mult->GetPhi(n), mult->GetDeltaPhi(n), mult->GetLabel(n));
2890 Printf("ERROR: AliMultiplicity could not be retrieved from ESD");
2893 delete [] usedTrack;
2897 // fill the tree for this event
2899 } // end of event loop
2901 aodTree->GetUserInfo()->Add(aod);
2903 // write the tree to the specified file
2904 aodFile = aodTree->GetCurrentFile();
2911 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2913 // Write space-points which are then used in the alignment procedures
2914 // For the moment only ITS, TRD and TPC
2916 // Load TOF clusters
2918 fLoader[3]->LoadRecPoints("read");
2919 TTree* tree = fLoader[3]->TreeR();
2921 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[3]));
2924 fTracker[3]->LoadClusters(tree);
2926 Int_t ntracks = esd->GetNumberOfTracks();
2927 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2929 AliESDtrack *track = esd->GetTrack(itrack);
2932 for (Int_t iDet = 3; iDet >= 0; iDet--)
2933 nsp += track->GetNcls(iDet);
2935 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2936 track->SetTrackPointArray(sp);
2938 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2939 AliTracker *tracker = fTracker[iDet];
2940 if (!tracker) continue;
2941 Int_t nspdet = track->GetNcls(iDet);
2942 if (nspdet <= 0) continue;
2943 track->GetClusters(iDet,idx);
2947 while (isp2 < nspdet) {
2949 TString dets = fgkDetectorName[iDet];
2950 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2951 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2952 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2953 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2954 isvalid = tracker->GetTrackPointTrackingError(idx[isp2],p,track);
2956 isvalid = tracker->GetTrackPoint(idx[isp2],p);
2959 const Int_t kNTPCmax = 159;
2960 if (iDet==1 && isp2>kNTPCmax) break; // to be fixed
2961 if (!isvalid) continue;
2962 sp->AddPoint(isptrack,&p); isptrack++; isp++;
2968 fTracker[3]->UnloadClusters();
2969 fLoader[3]->UnloadRecPoints();
2973 //_____________________________________________________________________________
2974 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2976 // The method reads the raw-data error log
2977 // accumulated within the rawReader.
2978 // It extracts the raw-data errors related to
2979 // the current event and stores them into
2980 // a TClonesArray inside the esd object.
2982 if (!fRawReader) return;
2984 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2986 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2988 if (iEvent != log->GetEventNumber()) continue;
2990 esd->AddRawDataErrorLog(log);
2995 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString fName){
2996 // Dump a file content into a char in TNamed
2998 in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2999 Int_t kBytes = (Int_t)in.tellg();
3000 printf("Size: %d \n",kBytes);
3003 char* memblock = new char [kBytes];
3004 in.seekg (0, ios::beg);
3005 in.read (memblock, kBytes);
3007 TString fData(memblock,kBytes);
3008 fn = new TNamed(fName,fData);
3009 printf("fData Size: %d \n",fData.Sizeof());
3010 printf("fName Size: %d \n",fName.Sizeof());
3011 printf("fn Size: %d \n",fn->Sizeof());
3015 AliInfo(Form("Could not Open %s\n",fPath.Data()));
3021 void AliReconstruction::TNamedToFile(TTree* fTree, TString fName){
3022 // This is not really needed in AliReconstruction at the moment
3023 // but can serve as a template
3025 TList *fList = fTree->GetUserInfo();
3026 TNamed *fn = (TNamed*)fList->FindObject(fName.Data());
3027 printf("fn Size: %d \n",fn->Sizeof());
3029 TString fTmp(fn->GetName()); // to be 100% sure in principle fName also works
3030 const char* cdata = fn->GetTitle();
3031 printf("fTmp Size %d\n",fTmp.Sizeof());
3033 int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
3034 printf("calculated size %d\n",size);
3035 ofstream out(fName.Data(),ios::out | ios::binary);
3036 out.write(cdata,size);
3041 //_____________________________________________________________________________
3042 AliQADataMakerRec * AliReconstruction::GetQADataMaker(Int_t iDet)
3044 // get the quality assurance data maker object and the loader for a detector
3046 if (fQADataMaker[iDet])
3047 return fQADataMaker[iDet];
3049 AliQADataMakerRec * qadm = NULL;
3050 if (iDet == fgkNDetectors) { //Global QA
3051 qadm = new AliGlobalQADataMaker();
3052 fQADataMaker[iDet] = qadm;
3056 // load the QA data maker object
3057 TPluginManager* pluginManager = gROOT->GetPluginManager();
3058 TString detName = fgkDetectorName[iDet];
3059 TString qadmName = "Ali" + detName + "QADataMakerRec";
3060 if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT"))
3063 // first check if a plugin is defined for the quality assurance data maker
3064 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
3065 // if not, add a plugin for it
3066 if (!pluginHandler) {
3067 AliDebug(1, Form("defining plugin for %s", qadmName.Data()));
3068 TString libs = gSystem->GetLibraries();
3069 if (libs.Contains("lib" + detName + "base.so") ||
3070 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
3071 pluginManager->AddHandler("AliQADataMakerRec", detName,
3072 qadmName, detName + "qadm", qadmName + "()");
3074 pluginManager->AddHandler("AliQADataMakerRec", detName,
3075 qadmName, detName, qadmName + "()");
3077 pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
3079 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
3080 qadm = (AliQADataMakerRec *) pluginHandler->ExecPlugin(0);
3083 fQADataMaker[iDet] = qadm;
3088 //_____________________________________________________________________________
3089 Bool_t AliReconstruction::RunQA(const char* detectors, AliESDEvent *& esd)
3091 // run the Quality Assurance data producer
3093 AliCodeTimerAuto("")
3094 TString detStr = detectors;
3095 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3096 if (!IsSelected(fgkDetectorName[iDet], detStr))
3098 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
3101 AliCodeTimerStart(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
3102 AliInfo(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
3104 qadm->Exec(AliQA::kESDS, esd) ;
3107 AliCodeTimerStop(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
3109 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
3110 AliError(Form("the following detectors were not found: %s",
3120 //_____________________________________________________________________________
3121 void AliReconstruction::CheckQA()
3123 // check the QA of SIM for this run and remove the detectors
3124 // with status Fatal
3126 TString newRunLocalReconstruction ;
3127 TString newRunTracking ;
3128 TString newFillESD ;
3130 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
3131 TString detName(AliQA::GetDetName(iDet)) ;
3132 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX(iDet)) ;
3133 if ( qa->IsSet(AliQA::DETECTORINDEX(iDet), AliQA::kSIM, AliQA::kFATAL)) {
3134 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
3136 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
3137 fRunLocalReconstruction.Contains("ALL") ) {
3138 newRunLocalReconstruction += detName ;
3139 newRunLocalReconstruction += " " ;
3141 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
3142 fRunTracking.Contains("ALL") ) {
3143 newRunTracking += detName ;
3144 newRunTracking += " " ;
3146 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
3147 fFillESD.Contains("ALL") ) {
3148 newFillESD += detName ;
3153 fRunLocalReconstruction = newRunLocalReconstruction ;
3154 fRunTracking = newRunTracking ;
3155 fFillESD = newFillESD ;
3158 //_____________________________________________________________________________
3159 Int_t AliReconstruction::GetDetIndex(const char* detector)
3161 // return the detector index corresponding to detector
3163 for (index = 0; index < fgkNDetectors ; index++) {
3164 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
3169 //_____________________________________________________________________________
3170 Bool_t AliReconstruction::FinishPlaneEff() {
3172 // Here execute all the necessary operationis, at the end of the tracking phase,
3173 // in case that evaluation of PlaneEfficiencies was required for some detector.
3174 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
3176 // This Preliminary version works only FOR ITS !!!!!
3177 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3180 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3183 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3184 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
3185 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3186 if(fTracker[iDet]) {
3187 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
3188 ret=planeeff->WriteIntoCDB();
3189 if(planeeff->GetCreateHistos()) {
3190 TString name="PlaneEffHisto";
3191 name+=fgkDetectorName[iDet];
3193 ret*=planeeff->WriteHistosToFile(name,"RECREATE");
3199 //_____________________________________________________________________________
3200 Bool_t AliReconstruction::InitPlaneEff() {
3202 // Here execute all the necessary operations, before of the tracking phase,
3203 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
3204 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
3205 // which should be updated/recalculated.
3207 // This Preliminary version will work only FOR ITS !!!!!
3208 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3211 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3213 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));