]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - STEER/STEER/AliReconstruction.cxx
Fixes for coverity.
[u/mrichter/AliRoot.git] / STEER / STEER / AliReconstruction.cxx
index b3eaf7b10d2255812541d454d3e8c859d1247085..c21aa1ebbd11d9f4dc1dd9aa3854ffa6f490f626 100644 (file)
@@ -13,7 +13,7 @@
  * provided "as is" without express or implied warranty.                  *
  **************************************************************************/
 
-/* $Id$ */
+/* $Id: AliReconstruction.cxx 63911 2013-08-19 16:46:41Z hristov $ */
 
 ///////////////////////////////////////////////////////////////////////////////
 //                                                                           //
 #include "AliLHCData.h"
 #include "ARVersion.h"
 #include <RVersion.h>
+#include <stdlib.h>
 #include <unistd.h>
 #include <sys/resource.h>
 ClassImp(AliReconstruction)
@@ -211,13 +212,9 @@ using std::endl;
 
 //_____________________________________________________________________________
 const char* AliReconstruction::fgkStopEvFName = "_stopEvent_";
-const char* AliReconstruction::fgkDetectorName[AliReconstruction::kNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE"
-// #ifdef MFT_UPGRADE
-//                                                                                   , "MFT"
-// #endif 
-                                                                                  , "MFT"    // AU
-                                                                                 , "HLT"
-};
+const char* AliReconstruction::fgkDetectorName[AliReconstruction::kNDetectors] = {"ITS", "TPC", "TRD",
+"TOF", "PHOS", 
+"HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE","AD","FIT","MFT", "HLT"};
 
 //_____________________________________________________________________________
 AliReconstruction::AliReconstruction(const char* gAliceFilename) :
@@ -225,6 +222,7 @@ AliReconstruction::AliReconstruction(const char* gAliceFilename) :
   fRunVertexFinder(kTRUE),
   fRunVertexFinderTracks(kTRUE),
   fRunMuonTracking(kFALSE),
+  fRunMFTTrackingMU(kFALSE),
   fRunV0Finder(kTRUE),
   fRunCascadeFinder(kTRUE),
   fRunMultFinder(kTRUE),
@@ -288,6 +286,7 @@ AliReconstruction::AliReconstruction(const char* gAliceFilename) :
   fCDBUri(),
   fQARefUri(),
   fSpecCDBUri(), 
+  fCheckRecoCDBvsSimuCDB(),
   fInitCDBCalled(kFALSE),
   fCDBSnapshotMode(kFALSE),
   fSetRunNumberFromDataCalled(kFALSE),
@@ -300,6 +299,8 @@ AliReconstruction::AliReconstruction(const char* gAliceFilename) :
   fWriteQAExpertData(kTRUE), 
   fRunPlaneEff(kFALSE),
 
+  fESDpid(NULL),
+
   fesd(NULL),
   fhltesd(NULL),
   fesdf(NULL),
@@ -317,6 +318,12 @@ AliReconstruction::AliReconstruction(const char* gAliceFilename) :
   fSspecie(0),
   fNhighPt(0),
   fShighPt(0),
+  //
+  fTreeBuffSize(30000000),
+  fMemCountESD(0),
+  fMemCountESDF(0),
+  fMemCountESDHLT(0),
+  //
   fUpgradeModule(""),
   fAnalysisMacro(),
   fAnalysis(0),
@@ -340,7 +347,9 @@ AliReconstruction::AliReconstruction(const char* gAliceFilename) :
     fQAWriteExpert[iDet] = kFALSE ; 
   }
   fBeamInt[0][0]=fBeamInt[0][1]=fBeamInt[1][0]=fBeamInt[1][1] = -1;
-
+  //
+  AddCheckRecoCDBvsSimuCDB("TPC/Calib/RecoParam"); // check for similarity in the sim and rec
+  //
   AliPID pid;
 }
 
@@ -350,6 +359,7 @@ AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
   fRunVertexFinder(rec.fRunVertexFinder),
   fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
   fRunMuonTracking(rec.fRunMuonTracking),
+  fRunMFTTrackingMU(rec.fRunMFTTrackingMU),
   fRunV0Finder(rec.fRunV0Finder),
   fRunCascadeFinder(rec.fRunCascadeFinder),
   fRunMultFinder(rec.fRunMultFinder),
@@ -413,6 +423,7 @@ AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
   fCDBUri(rec.fCDBUri),
   fQARefUri(rec.fQARefUri),
   fSpecCDBUri(), 
+  fCheckRecoCDBvsSimuCDB(),
   fInitCDBCalled(rec.fInitCDBCalled),
   fCDBSnapshotMode(rec.fCDBSnapshotMode),
   fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
@@ -425,6 +436,8 @@ AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
   fWriteQAExpertData(rec.fWriteQAExpertData), 
   fRunPlaneEff(rec.fRunPlaneEff),
 
+  fESDpid(NULL),
+
   fesd(NULL),
   fhltesd(NULL),
   fesdf(NULL),
@@ -442,6 +455,12 @@ AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
   fSspecie(0),
   fNhighPt(0),
   fShighPt(0),
+  //
+  fTreeBuffSize(rec.fTreeBuffSize),
+  fMemCountESD(0),
+  fMemCountESDF(0),
+  fMemCountESDHLT(0),
+  //
   fUpgradeModule(""),
   fAnalysisMacro(rec.fAnalysisMacro),
   fAnalysis(0),
@@ -472,6 +491,10 @@ AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
     if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
   }
 
+  for (Int_t i = 0; i < rec.fCheckRecoCDBvsSimuCDB.GetEntriesFast(); i++) {
+    if (rec.fCheckRecoCDBvsSimuCDB[i]) fCheckRecoCDBvsSimuCDB.AddLast(rec.fCheckRecoCDBvsSimuCDB[i]->Clone());
+  }
+
   for (int i=2;i--;) for (int j=2;j--;) fBeamInt[i][j] = rec.fBeamInt[i][j];
 
 }
@@ -491,6 +514,7 @@ AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
   fRunVertexFinder       = rec.fRunVertexFinder;
   fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
   fRunMuonTracking       = rec.fRunMuonTracking;
+  fRunMFTTrackingMU      = rec.fRunMFTTrackingMU;
   fRunV0Finder           = rec.fRunV0Finder;
   fRunCascadeFinder      = rec.fRunCascadeFinder;
   fRunMultFinder         = rec.fRunMultFinder;
@@ -586,6 +610,12 @@ AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
   fCDBUri        = "";
   fQARefUri      = rec.fQARefUri;
   fSpecCDBUri.Delete();
+  fCheckRecoCDBvsSimuCDB.Delete();
+  //
+  for (Int_t i = 0; i < rec.fCheckRecoCDBvsSimuCDB.GetEntriesFast(); i++) {
+    if (rec.fCheckRecoCDBvsSimuCDB[i]) fCheckRecoCDBvsSimuCDB.AddLast(rec.fCheckRecoCDBvsSimuCDB[i]->Clone());
+  }
+  //
   fInitCDBCalled               = rec.fInitCDBCalled;
   fCDBSnapshotMode             = rec.fCDBSnapshotMode;
   fSetRunNumberFromDataCalled  = rec.fSetRunNumberFromDataCalled;
@@ -598,6 +628,7 @@ AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
   fWriteQAExpertData           = rec.fWriteQAExpertData;
   fRunPlaneEff                 = rec.fRunPlaneEff;
   for (int i=2;i--;) for (int j=2;j--;) fBeamInt[i][j] = rec.fBeamInt[i][j];
+  fESDpid  = NULL;
   fesd     = NULL;
   fhltesd  = NULL;
   fesdf    = NULL;
@@ -615,6 +646,12 @@ AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
   fSspecie = 0;
   fNhighPt = 0;
   fShighPt = 0;
+  //
+  fTreeBuffSize = rec.fTreeBuffSize;
+  fMemCountESD = 0;
+  fMemCountESDF = 0;
+  fMemCountESDHLT = 0;
+  //
   fUpgradeModule="";
   fAnalysisMacro = rec.fAnalysisMacro;
   fAnalysis = 0;
@@ -644,7 +681,7 @@ AliReconstruction::~AliReconstruction()
     delete fAlignObjArray;
   }
   fSpecCDBUri.Delete();
-
+  fCheckRecoCDBvsSimuCDB.Delete();
   AliCodeTimer::Instance()->Print();
 }
 
@@ -807,6 +844,7 @@ void AliReconstruction::SetQARefDefaultStorage(const char* uri) {
   AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
   
 }
+
 //_____________________________________________________________________________
 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
 // Store a detector-specific CDB storage location
@@ -849,6 +887,48 @@ void AliReconstruction::SetSpecificStorage(const char* calibType, const char* ur
 
 }
 
+//_____________________________________________________________________________
+void AliReconstruction::AddCheckRecoCDBvsSimuCDB(const char* cdbpath,const char* comment) 
+{
+  // require the cdb item to be the same in the rec as in the sim
+  // Activate it later within the Run() method
+  TString newent = cdbpath;
+  if (newent.IsNull()) return;
+  TIter nextit(&fCheckRecoCDBvsSimuCDB);
+  TNamed* cdbent=0;
+  while ((cdbent=(TNamed*)nextit())) {
+    TString str = cdbent->GetName();
+    if (str==newent) {
+      AliInfo(Form("%s is already in the list to check",cdbpath));
+      return;
+    }
+  }
+  fCheckRecoCDBvsSimuCDB.AddLast(new TNamed(cdbpath,comment));
+  //
+}
+
+//_____________________________________________________________________________
+void AliReconstruction::RemCheckRecoCDBvsSimuCDB(const char* cdbpath) 
+{
+  // require the cdb item to be the same in the rec as in the sim
+  // Activate it later within the Run() method
+  TString newent = cdbpath;
+  if (newent.IsNull()) return;
+  TIter nextit(&fCheckRecoCDBvsSimuCDB);
+  TNamed* cdbent=0;
+  while ((cdbent=(TNamed*)nextit())) {
+    TString str = cdbent->GetName();
+    if (str==newent) {
+      AliInfo(Form("Removing %s from the list to check",cdbpath));
+      delete fCheckRecoCDBvsSimuCDB.Remove(cdbent);
+      fCheckRecoCDBvsSimuCDB.Compress();
+      return;
+    }
+  }
+  AliInfo(Form("%s is not in the list to check",cdbpath));
+  //
+}
+
 //_____________________________________________________________________________
 Bool_t AliReconstruction::SetRunNumberFromData()
 {
@@ -1277,6 +1357,7 @@ Bool_t AliReconstruction::LoadCDB()
   // in the trigger or that are needed in order to put correct
   // information in ESD
   AliCDBManager::Instance()->GetAll("TRIGGER/*/*");
+  AliCDBManager::Instance()->GetAll("HLT/*/*");
 
   return kTRUE;
 }
@@ -1437,6 +1518,7 @@ Bool_t AliReconstruction::Run(const char* input)
         Abort("ProcessEvent",TSelector::kAbortFile);
         return kFALSE;
       }
+      CleanProcessedEvent();
       iEvent++;
     }
     if (!iEvent) AliWarning("No events passed trigger selection");
@@ -1529,7 +1611,7 @@ void AliReconstruction::Begin(TTree *)
   }
 
   // Import ideal TGeo geometry and apply misalignment
-  if (!gGeoManager) {
+  if (!AliGeomManager::GetGeometry()) {
     TString geom(gSystem->DirName(fGAliceFileName));
     geom += "/geometry.root";
     AliGeomManager::LoadGeometry(geom.Data());
@@ -1563,7 +1645,10 @@ void AliReconstruction::Begin(TTree *)
     Abort("MisalignGeometry", TSelector::kAbortProcess);
     return;
   }
-  AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
+
+  const TMap* cdbCache = AliCDBManager::Instance()->GetEntryCache();
+  if(cdbCache->Contains("GRP/Geometry/Data"))
+         AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
   if(!toCDBSnapshot) AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
   AliSysInfo::AddStamp("MisalignGeom");
 
@@ -1572,7 +1657,9 @@ void AliReconstruction::Begin(TTree *)
     return;
   }
   AliSysInfo::AddStamp("InitGRP");
-  if(!toCDBSnapshot) AliCDBManager::Instance()->UnloadFromCache("GRP/Calib/CosmicTriggers");
+  if(!toCDBSnapshot)
+      if(cdbCache->Contains("GRP/Calib/CosmicTriggers"))
+         AliCDBManager::Instance()->UnloadFromCache("GRP/Calib/CosmicTriggers");
 
   if(!fCDBSnapshotMode || toCDBSnapshot){
       if (!LoadCDB()) {
@@ -1609,8 +1696,7 @@ void AliReconstruction::Begin(TTree *)
   if(toCDBSnapshot)
   {
       AliCDBManager::Instance()->DumpToSnapshotFile(snapshotFileOut.Data(),kFALSE);
-      AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
-      AliCDBManager::Instance()->UnloadFromCache("GRP/Calib/CosmicTriggers");
+      exit(0);
   }
 
   if (fInput && gProof) {
@@ -1720,6 +1806,8 @@ void AliReconstruction::SlaveBegin(TTree*)
   }
   AliSysInfo::AddStamp("LoadLoader");
  
+  CheckRecoCDBvsSimuCDB();
+
   ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
 
   // get trackers
@@ -1800,6 +1888,9 @@ void AliReconstruction::SlaveBegin(TTree*)
   gSystem->GetProcInfo(&procInfo);
   AliInfo(Form("Current memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
   
+  // PID
+  fESDpid = new AliESDpid();
+
   //QA
   //Initialize the QA and start of cycle 
   if (fRunQA || fRunGlobalQA) 
@@ -1883,14 +1974,12 @@ Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
   static Long_t oldMres=0;
   static Long_t oldMvir=0;
   static Float_t oldCPU=0;
-  static Long_t aveDMres=0;
-  static Long_t aveDMvir=0;
-  static Float_t aveDCPU=0;
+  // static Long_t aveDMres=0;
+  // static Long_t aveDMvir=0;
+  // static Float_t aveDCPU=0;
 
   AliCodeTimerAuto("",0);
 
-  AliESDpid pid;
-
   AliSysInfo::AddStamp(Form("StartEv_%d",iEvent), 0,0,iEvent);
 
   if (iEvent >= fRunLoader->GetNumberOfEvents()) {
@@ -1935,7 +2024,7 @@ Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
       if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
         const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
         reconstructor->SetRecoParam(par);
-       reconstructor->GetPidSettings(&pid);
+       reconstructor->GetPidSettings(fESDpid);
        reconstructor->SetEventInfo(&fEventInfo);
         if (fRunQA) {
           AliQAManager::QAManager()->SetEventInfo(&fEventInfo) ;
@@ -1964,6 +2053,11 @@ Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
     // fill Event header information from the RawEventHeader
     if (fRawReader){FillRawEventHeaderESD(fesd);}
     if (fRawReader){FillRawEventHeaderESD(fhltesd);}
+    if (fRawReader){
+      // Store DAQ detector pattern and attributes
+      fesd->SetDAQDetectorPattern(fRawReader->GetDetectorPattern()[0]);
+      fesd->SetDAQAttributes(fRawReader->GetAttributes()[2]);
+    }
 
     fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
     fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
@@ -2089,9 +2183,23 @@ Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
       AliSysInfo::AddStamp(Form("TrackingMUON_%d",iEvent), 0,0,iEvent);      
     }
 
+    //---------------- AU From here...
+
+    // MFT tracking of MUON tracks
+    if (!fRunTracking.IsNull()) {
+      if (fRunMFTTrackingMU && fRunMuonTracking) {
+       if (!RunMFTTrackingMU(fesd)) {
+         if (fStopOnError) {CleanUp(); return kFALSE;}
+       }
+      }
+      AliSysInfo::AddStamp(Form("TrackingMFT_MUON_%d",iEvent), 0,0,iEvent);      
+    }
+
+    //---------------- ...to here
+
     // barrel tracking
     if (!fRunTracking.IsNull()) {
-      if (!RunTracking(fesd,pid)) {
+      if (!RunTracking(fesd,*fESDpid)) {
        if (fStopOnError) {CleanUp(); return kFALSE;}
       }
     }
@@ -2135,7 +2243,7 @@ Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
       ok = kFALSE;
       if (tpcTrack)
        ok = AliTracker::
-         PropagateTrackToBxByBz(tpcTrack,kRadius,track->GetMass(),kMaxStep,kFALSE);
+         PropagateTrackToBxByBz(tpcTrack,kRadius,track->GetMassForTracking(),kMaxStep,kFALSE);
 
       if (ok) {
        Int_t n=trkArray.GetEntriesFast();
@@ -2147,7 +2255,7 @@ Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
       if (track->IsOn(AliESDtrack::kITSrefit)) continue;
 
       AliTracker::
-         PropagateTrackToBxByBz(track,kRadius,track->GetMass(),kMaxStep,kFALSE);
+         PropagateTrackToBxByBz(track,kRadius,track->GetMassForTracking(),kMaxStep,kFALSE);
       Double_t x[3]; track->GetXYZ(x);
       Double_t b[3]; AliTracker::GetBxByBz(x,b);
       track->RelateToVertexBxByBz(fesd->GetPrimaryVertexSPD(), b, kVeryBig);
@@ -2272,10 +2380,10 @@ Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
 
     // AdC+FN
     if (fReconstructor[3])
-      GetReconstructor(3)->FillEventTimeWithTOF(fesd,&pid);
+      GetReconstructor(3)->FillEventTimeWithTOF(fesd,fESDpid);
 
     // combined PID
-    pid.MakePID(fesd);
+    //    fESDpid->MakePID(fesd);
 
     if (fFillTriggerESD) {
       if (!FillTriggerESD(fesd)) {
@@ -2346,7 +2454,16 @@ Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
   
   }
   //
-  ftree->Fill();
+  Long64_t nbf;
+  nbf = ftree->Fill();
+  if (fTreeBuffSize>0 && ftree->GetAutoFlush()<0 && (fMemCountESD += nbf)>fTreeBuffSize ) { // default limit is still not reached
+    nbf = ftree->GetZipBytes();
+    if (nbf>0) nbf = -nbf;
+    else       nbf = ftree->GetEntries();
+    ftree->SetAutoFlush(nbf);
+    AliInfo(Form("Calling ftree->SetAutoFlush(%lld) | W:%lld T:%lld Z:%lld",
+                nbf,fMemCountESD,ftree->GetTotBytes(),ftree->GetZipBytes()));        
+  }
   AliSysInfo::AddStamp(Form("ESDFill_%d",iEvent), 0,0,iEvent);     
   //
   if (fWriteESDfriend) {
@@ -2360,11 +2477,24 @@ Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
     ftree->AutoSave("SaveSelf");
     if (fWriteESDfriend) ftreeF->AutoSave("SaveSelf");
   }
-    // write HLT ESD
-    fhlttree->Fill();
+  // write HLT ESD
+  
+  nbf = fhlttree->Fill();
+  if (fTreeBuffSize>0 && fhlttree->GetAutoFlush()<0 && (fMemCountESDHLT += nbf)>fTreeBuffSize ) { // default limit is still not reached
+    nbf = fhlttree->GetZipBytes();
+    if (nbf>0) nbf = -nbf;
+    else       nbf = fhlttree->GetEntries();
+    fhlttree->SetAutoFlush(nbf);
+    AliInfo(Form("Calling fhlttree->SetAutoFlush(%lld) | W:%lld T:%lld Z:%lld",
+                nbf,fMemCountESDHLT,fhlttree->GetTotBytes(),fhlttree->GetZipBytes()));        
+  }
+    
+    
+  return kTRUE;
+}
 
-    // call AliEVE
-    if (fRunAliEVE) RunAliEVE();
+void AliReconstruction::CleanProcessedEvent()
+{
     //
     fesd->Reset();
     fhltesd->Reset();
@@ -2377,19 +2507,8 @@ Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
       if (fReconstructor[iDet]) fReconstructor[iDet]->FinishEvent();
     }
  
-    gSystem->GetProcInfo(&procInfo);
-    Long_t dMres=(procInfo.fMemResident-oldMres)/1024;
-    Long_t dMvir=(procInfo.fMemVirtual-oldMvir)/1024;
-    Float_t dCPU=procInfo.fCpuUser+procInfo.fCpuSys-oldCPU;
-    aveDMres+=(dMres-aveDMres)/(iEvent-fFirstEvent+1);
-    aveDMvir+=(dMvir-aveDMvir)/(iEvent-fFirstEvent+1);
-    aveDCPU+=(dCPU-aveDCPU)/(iEvent-fFirstEvent+1);
-    AliInfo(Form("======================= End Event %d: Res %ld(%3ld <%3ld>) Vir %ld(%3ld <%3ld>) CPU %5.2f <%5.2f> ===================",
-                iEvent, procInfo.fMemResident/1024, dMres, aveDMres, procInfo.fMemVirtual/1024, dMvir, aveDMvir, dCPU, aveDCPU));
-    oldMres=procInfo.fMemResident;
-    oldMvir=procInfo.fMemVirtual;
-    oldCPU=procInfo.fCpuUser+procInfo.fCpuSys;
-  
+    AliInfo("======================= End Event ===================");
+    
     fEventInfo.Reset();
     for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
       if (fReconstructor[iDet]) {
@@ -2405,7 +2524,7 @@ Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
   DeleteRecPoints(fDeleteRecPoints);
   DeleteDigits(fDeleteDigits);
   //
-  return kTRUE;
+
 }
 
 //_____________________________________________________________________________
@@ -2461,13 +2580,25 @@ void AliReconstruction::SlaveTerminate()
 
    // Add the AliRoot version that created this file
    TString sVersion("aliroot ");
-   sVersion += ALIROOT_SVN_BRANCH;
+   sVersion += ALIROOT_BRANCH;
    sVersion += ":";
-   sVersion += ALIROOT_SVN_REVISION;
+   sVersion += ALIROOT_REVISION;
    sVersion += "; root ";
+#ifdef ROOT_SVN_BRANCH
    sVersion += ROOT_SVN_BRANCH;
+#elif defined(ROOT_GIT_BRANCH)
+   sVersion += ROOT_GIT_BRANCH;
+#else
+   sVersion += "?";
+#endif
    sVersion += ":";
+#ifdef ROOT_SVN_REVSION
    sVersion += ROOT_SVN_REVISION;
+#elif defined(ROOT_GIT_COMMIT)
+   sVersion += ROOT_GIT_COMMIT;
+#else 
+   sVersion += "?";
+#endif
    sVersion += "; metadata ";
    sVersion += getenv("PRODUCTION_METADATA");
                    
@@ -2866,6 +2997,53 @@ Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
 }
 
 
+//_____________________________________________________________________________
+Bool_t AliReconstruction::RunMFTTrackingMU(AliESDEvent*& esd) {
+
+  // AU
+
+  // run the global muon tracking: matching the MUON tracks with the MFT clusters
+
+  AliCodeTimerAuto("",0)
+
+  if (!fRunLoader) {
+    AliError("Missing runLoader!");
+    return kFALSE;
+  }
+  Int_t iDet = GetDetIndex("MFT"); // for MFT
+
+  // Get a pointer to the MFT reconstructor
+  AliReconstructor *reconstructor = GetReconstructor(iDet);
+  if (!reconstructor) return kFALSE;
+  
+  TString detName = fgkDetectorName[iDet];
+  AliDebug(1, Form("%s tracking for muon tracks", detName.Data()));
+  AliTracker *tracker = reconstructor->CreateTracker();
+  if (!tracker) {
+    AliWarning(Form("couldn't create a Muon tracker for %s", detName.Data()));
+    return kFALSE;
+  }
+     
+  // read RecPoints
+  fLoader[iDet]->LoadRecPoints("read");  
+
+  tracker->LoadClusters(fLoader[iDet]->TreeR());
+  
+  Int_t rv = tracker->Clusters2Tracks(esd);
+  
+  fLoader[iDet]->UnloadRecPoints();
+
+  tracker->UnloadClusters();
+  
+  if (rv) {
+    AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
+    return kFALSE;
+  }
+  
+  return kTRUE;
+
+}
+
 //_____________________________________________________________________________
 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd,AliESDpid &PID)
 {
@@ -2918,7 +3096,7 @@ Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd,AliESDpid &PID)
     // preliminary PID in TPC needed by the ITS tracker
     if (iDet == 1) {
       GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
-      PID.MakePID(esd,kTRUE);
+      PID.MakePIDForTracking(esd);
       AliSysInfo::AddStamp(Form("MakePID0%s_%d",fgkDetectorName[iDet],eventNr), iDet,4,eventNr);
     } 
   }
@@ -2973,7 +3151,7 @@ Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd,AliESDpid &PID)
     if (iDet == 1) {
       //GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
       //AliESDpid::MakePID(esd);
-      PID.MakePID(esd,kTRUE);
+      PID.MakePIDForTracking(esd);
       AliSysInfo::AddStamp(Form("MakePID1%s_%d",fgkDetectorName[iDet],eventNr), iDet,4,eventNr);
     }
 
@@ -3465,6 +3643,10 @@ Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
       fRunMuonTracking = kTRUE;
       continue;
     }
+    if (detName == "MFT") {           // AU    
+      fRunMFTTrackingMU = kTRUE;      // AU
+      continue;                              // AU
+    }                                 // AU
 
     fTracker[iDet] = reconstructor->CreateTracker();
     if (!fTracker[iDet] && (iDet < 7)) {
@@ -3505,6 +3687,9 @@ void AliReconstruction::CleanUp()
   delete fParentRawReader;
   fParentRawReader=NULL;
 
+  delete fESDpid;
+  fESDpid = NULL;
+
   if (ffile) {
     ffile->Close();
     delete ffile;
@@ -3963,6 +4148,9 @@ Bool_t AliReconstruction::GetEventInfo()
   }
 
   // Load trigger aliases and declare the trigger classes included in aliases
+  //PH Why do we do it in each event and not only once in the beginning of the chunk??
+  //PH Temporary fix for #99725: AliReconstruction::GetEventInfo bug
+  fDeclTriggerClasses.Clear();
   AliCDBEntry * entry = AliCDBManager::Instance()->Get("GRP/CTP/Aliases");
   if (entry) {
     THashList * lst = dynamic_cast<THashList*>(entry->GetObject());
@@ -4340,7 +4528,16 @@ void AliReconstruction::WriteESDfriend() {
     fesdf->SetSkipBit(kTRUE);
   }
   //
-  ftreeF->Fill();
+  Long64_t nbf = ftreeF->Fill();
+  if (fTreeBuffSize>0 && ftreeF->GetAutoFlush()<0 && (fMemCountESDF += nbf)>fTreeBuffSize ) { // default limit is still not reached
+    nbf = ftreeF->GetZipBytes();
+    if (nbf>0) nbf = -nbf;
+    else       nbf = ftreeF->GetEntries();
+    ftreeF->SetAutoFlush(nbf);
+    AliInfo(Form("Calling ftreeF->SetAutoFlush(%lld) | W:%lld T:%lld Z:%lld",
+                nbf,fMemCountESDF,ftreeF->GetTotBytes(),ftreeF->GetZipBytes()));        
+  }
+  
 }
 
 //_________________________________________________________________
@@ -4441,3 +4638,125 @@ Bool_t AliReconstruction::HasEnoughResources(int ev)
   }
   return res;
 }
+
+Bool_t AliReconstruction::HasNextEventAfter(Int_t eventId)
+{
+        return ( (eventId < fRunLoader->GetNumberOfEvents()) ||
+          (fRawReader && fRawReader->NextEvent()) );
+}
+
+//_________________________________________________________________
+void AliReconstruction::CheckRecoCDBvsSimuCDB()
+{
+  // if some CDB entries must be the same in the simulation
+  // and reconstruction, check here
+  int nent = fCheckRecoCDBvsSimuCDB.GetEntriesFast();
+  AliInfo(Form("Check %d entries for matching between sim and rec",nent));
+  //
+  // get simulation CDB
+  fRunLoader->CdGAFile();
+  TMap*  cdbMapSim  = (TMap*)gDirectory->Get("cdbMap");
+  TList* cdbListSim = (TList*)gDirectory->Get("cdbList");
+  if (!(cdbMapSim && cdbListSim)) {
+    AliInfo(Form("No CDBMap/List found in %s, nothing to check",fGAliceFileName.Data()));
+    return;
+  }
+  // read the requested objects to make sure they will appear in the reco list
+  for (Int_t i=0;i<nent;i++) {
+    TNamed* cdbent = (TNamed*) fCheckRecoCDBvsSimuCDB[i];
+    if (!cdbent) continue;
+    AliCDBManager::Instance()->Get(cdbent->GetName());
+  }
+  // get default path for simulation
+  TPair* pair;
+  TObjString* stro;
+  pair = (TPair*)cdbMapSim->FindObject("default");
+  if (!pair) {AliFatal("Did not find default storage used for simulations"); return;}
+  TString defSimStore = ((TObjString*)pair->Value())->GetString();
+  RectifyCDBurl(defSimStore);
+  //
+  // get reconstruction CDB
+  const TMap *cdbMapRec = AliCDBManager::Instance()->GetStorageMap();   
+  const TList *cdbListRec = AliCDBManager::Instance()->GetRetrievedIds();
+  //
+  // get default path for reconstruction
+  pair = (TPair*)cdbMapRec->FindObject("default");
+  if (!pair) {AliFatal("Did not find default storage used for reconstruction"); return;}
+  TString defRecStore = ((TObjString*)pair->Value())->GetString();
+  RectifyCDBurl(defRecStore);
+  //
+  for (Int_t i=0;i<nent;i++) {
+    TNamed* cdbent = (TNamed*) fCheckRecoCDBvsSimuCDB[i];
+    if (!cdbent) continue;
+    //
+    AliInfo(Form("#%d Checking %s",i,cdbent->GetName()));
+    //
+    // find cdbID used for sim
+    TString idSim="",storSim="";
+    TIter nextSim(cdbListSim);
+    while ((stro=(TObjString*)nextSim())) {
+      if (stro->GetString().Contains(cdbent->GetName())) {
+       idSim = stro->GetString();
+       break;
+      }
+    }    
+    // find the storage used for sim
+    // check in the simuCDB special paths
+    pair = (TPair*)cdbMapSim->FindObject(cdbent->GetName());
+    if (pair) { // specific path is used
+      storSim = ((TObjString*)pair->Value())->GetString();
+      RectifyCDBurl(storSim);
+    }
+    else storSim = defSimStore;  // default storage list is used
+    //
+    if (!idSim.IsNull()) AliInfo(Form("Sim. used %s from %s",idSim.Data(), storSim.Data()));
+    else                 AliInfo("Sim. did not use this object");
+    //
+    // find cdbID used for rec
+    TString idRec="",storRec="";
+    TIter nextRec(cdbListRec);
+    AliCDBId* id=0;
+    while ((id=(AliCDBId*)nextRec())) {
+      idRec = id->ToString();
+      if (idRec.Contains(cdbent->GetName())) break;
+      idRec="";
+    }
+    //
+    // find storage used for the rec
+    pair = (TPair*)cdbMapRec->FindObject(cdbent->GetName());
+    if (pair) {  // specific path is used
+      storRec = ((TObjString*)pair->Value())->GetString();
+      RectifyCDBurl(storRec);
+    }
+    else storRec = defRecStore; // default storage list is used
+    //
+    if (!idRec.IsNull()) AliInfo(Form("Rec. used %s from %s",idRec.Data(), storRec.Data()));
+    else                 AliInfo("Rec. did not use this object");
+    //
+    if (!idSim.IsNull() && !idRec.IsNull() && ((idSim!=idRec) || (storSim!=storRec)) ) 
+      AliFatal("Different objects were used in sim and rec");
+  }
+  
+}
+
+//_________________________________________________________
+void AliReconstruction::RectifyCDBurl(TString& url)
+{
+  // TBD RS
+  // remove everything but the url
+  TString sbs;
+  if (!(sbs=url("\\?User=[^?]*")).IsNull())                url.ReplaceAll(sbs,"");
+  if (!(sbs=url("\\?DBFolder=[^?]*")).IsNull())            url.ReplaceAll("?DB","");
+  if (!(sbs=url("\\?SE=[^?]*")).IsNull())                  url.ReplaceAll(sbs,"");
+  if (!(sbs=url("\\?CacheFolder=[^?]*")).IsNull())         url.ReplaceAll(sbs,"");
+  if (!(sbs=url("\\?OperateDisconnected=[^?]*")).IsNull()) url.ReplaceAll(sbs,"");
+  if (!(sbs=url("\\?CacheSize=[^?]*")).IsNull())           url.ReplaceAll(sbs,"");  
+  if (!(sbs=url("\\?CleanupInterval=[^?]*")).IsNull())     url.ReplaceAll(sbs,"");  
+  Bool_t slash=kFALSE,space=kFALSE;
+  while ( (slash=url.EndsWith("/")) || (space=url.EndsWith(" ")) ) {
+    if (slash) url = url.Strip(TString::kTrailing,'/');
+    if (space) url = url.Strip(TString::kTrailing,' ');
+  }
+  //url.ToLower();
+  //
+}