]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - STEER/STEER/AliReconstruction.cxx
Implement comparison of sim and rec CDB's, by default for the TPC/RecoParam only
[u/mrichter/AliRoot.git] / STEER / STEER / AliReconstruction.cxx
index ae3c39646afb9d22d9a86b82f5acb22cc5958df2..453687716252a18bb86dd7bbfd72049162a77ea6 100644 (file)
@@ -211,7 +211,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","AD","MFT", "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) :
@@ -219,6 +221,7 @@ AliReconstruction::AliReconstruction(const char* gAliceFilename) :
   fRunVertexFinder(kTRUE),
   fRunVertexFinderTracks(kTRUE),
   fRunMuonTracking(kFALSE),
+  fRunMFTTrackingMU(kFALSE),
   fRunV0Finder(kTRUE),
   fRunCascadeFinder(kTRUE),
   fRunMultFinder(kTRUE),
@@ -282,6 +285,7 @@ AliReconstruction::AliReconstruction(const char* gAliceFilename) :
   fCDBUri(),
   fQARefUri(),
   fSpecCDBUri(), 
+  fCheckRecoCDBvsSimuCDB(),
   fInitCDBCalled(kFALSE),
   fCDBSnapshotMode(kFALSE),
   fSetRunNumberFromDataCalled(kFALSE),
@@ -342,7 +346,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;
 }
 
@@ -352,6 +358,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),
@@ -415,6 +422,7 @@ AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
   fCDBUri(rec.fCDBUri),
   fQARefUri(rec.fQARefUri),
   fSpecCDBUri(), 
+  fCheckRecoCDBvsSimuCDB(),
   fInitCDBCalled(rec.fInitCDBCalled),
   fCDBSnapshotMode(rec.fCDBSnapshotMode),
   fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
@@ -482,6 +490,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];
 
 }
@@ -501,6 +513,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;
@@ -596,6 +609,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;
@@ -661,7 +680,7 @@ AliReconstruction::~AliReconstruction()
     delete fAlignObjArray;
   }
   fSpecCDBUri.Delete();
-
+  fCheckRecoCDBvsSimuCDB.Delete();
   AliCodeTimer::Instance()->Print();
 }
 
@@ -824,6 +843,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
@@ -866,6 +886,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()
 {
@@ -1745,6 +1807,8 @@ void AliReconstruction::SlaveBegin(TTree*)
   }
   AliSysInfo::AddStamp("LoadLoader");
  
+  CheckRecoCDBvsSimuCDB();
+
   ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
 
   // get trackers
@@ -2120,6 +2184,20 @@ 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,*fESDpid)) {
@@ -2920,6 +2998,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)
 {
@@ -3519,6 +3644,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)) {
@@ -4516,3 +4645,147 @@ 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()));
+    //
+    // check in the simuCDB special params
+    pair = (TPair*)cdbMapSim->FindObject(cdbent->GetName());
+    TString idSimD = "";
+    TString idSimS = "";
+    if (pair) { // specific path is used
+      idSimS = ((TObjString*)pair->Value())->GetString();
+      RectifyCDBurl(idSimS);
+    }
+    else { // check in default storage list
+      TIter nextSim(cdbListSim);
+      while ((stro=(TObjString*)nextSim())) {
+       if (stro->GetString().Contains(cdbent->GetName())) {
+         idSimD = stro->GetString();
+         break;
+       }
+      }
+    }
+    //
+    // check in the recoCDB special params
+    pair = (TPair*)cdbMapRec->FindObject(cdbent->GetName());
+    TString idRecD = "";
+    TString idRecS = "";
+    if (pair) {  // specific path is used
+      idRecS = ((TObjString*)pair->Value())->GetString();
+      RectifyCDBurl(idRecS);
+    }
+    else { // check in default storage list
+      TIter nextRec(cdbListRec);
+      while ((stro=(TObjString*)nextRec())) {
+       if (stro->GetString().Contains(cdbent->GetName())) {
+         idRecD = stro->GetString();
+         break;
+       }
+      }
+    }
+    //-----------------------------
+    Bool_t ok = kTRUE;
+    if (!idSimD.IsNull()) {  // simulation used object from default storage
+      AliInfo(Form("Simulation used default storage %s\nentry %s",defSimStore.Data(),idSimD.Data()));
+      if (!idRecD.IsNull()) { // reco also
+       AliInfo(Form("Reconstruction used default storage %s\nentry %s",defRecStore.Data(),idRecD.Data()));
+       if ( (idSimD!=idRecD) || (defSimStore!=defRecStore) ) ok = kFALSE;
+      }
+      else if (!idRecS.IsNull()) { // reco used specific storage, strict check of version is not possible
+       AliInfo(Form("Reconstruction used specific storage %s",idRecS.Data()));
+       if (defSimStore!=idRecS) ok = kFALSE;
+      }
+      else {
+       AliInfo("Did not find object used in reconstruction");
+       ok = kFALSE;
+      }
+    }
+    else if (!idSimS.IsNull()) { // simulation used object from specific storage
+      AliInfo(Form("Simulation used specific storage %s",idSimS.Data()));
+      if (!idRecS.IsNull()) { // reco also     
+       AliInfo(Form("Reconstruction used specific storage %s",idRecS.Data()));
+       if (idSimS!=idRecS) ok = kFALSE;
+      }
+      else if (!idRecD.IsNull()) {
+       AliInfo(Form("Reconstruction used default storage %s\nentry",idRecD.Data()));
+       if (idSimS!=defRecStore) ok = kFALSE;
+      }
+      else {
+       AliInfo("Did not find object used in reconstruction");
+       ok = kFALSE;
+      }      
+    }
+    else {
+      AliInfo("Did not find object used in simulation");
+      ok = kFALSE;      
+    }
+    if (!ok) 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();
+  //
+}