]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - ACORDE/AliACORDEQAChecker.cxx
changed order of adding histograms to the TList
[u/mrichter/AliRoot.git] / ACORDE / AliACORDEQAChecker.cxx
index b8486de2dc44a18307762fc221d036912a510148..3333119d1981ee82eff3ac9bbc6b02d8fd5ad719 100755 (executable)
@@ -19,6 +19,7 @@
 //     Mario Rodriguez Cahuantzi <mrodrigu@mail.cern.ch> (FCFM-BUAP) 
 //     Luciano Diaz Gonzalez <luciano.diaz@nucleares.unam.mx> (ICN-UNAM)
 //     Arturo Fernandez <afernan@mail.cern.ch> (FCFM-BUAP)
+//  Last update: Nov. 14t 2009 --> MRC <mrodrigu@mail.cern.ch> (FCFM-BUAP) 
 //...
 
 // --- ROOT system ---
@@ -28,7 +29,8 @@
 #include <TIterator.h> 
 #include <TKey.h> 
 #include <TFile.h> 
-
+#include <TPaveText.h>
+#include <TLine.h>
 // --- Standard library ---
 
 // --- AliRoot header files ---
 #include "AliCDBEntry.h"
 #include "AliQAManager.h"
 
+/*************************************************************************
+   Last update: Oct. 10th 2011 from Mario RC, <mrodrigu@mail.cern.ch>
+       |-> Adding the checker class for raw and esd index
+       |-> Setting the threshold lines and box for DQM shifter
+
+*************************************************************************/
+
 ClassImp(AliACORDEQAChecker)
 
 //____________________________________________________________________________
-Double_t * AliACORDEQAChecker::Check(AliQAv1::ALITASK_t /*index*/)
+
+AliACORDEQAChecker::AliACORDEQAChecker()  :
+AliQACheckerBase("ACORDE","ACORDE Quality Assurance Data Maker"),
+fTextDQMShifterInfo(new TPaveText(2.2,0.90,2.9,0.99,"T")),
+fMax(new TLine(1,0.5,3,0.5))
+{
+       // default constructor
+       fMax->SetLineColor(kGreen);
+       fMax->SetLineWidth(3);
+}
+//____________________________________________________________________________
+AliACORDEQAChecker::~AliACORDEQAChecker()
+{
+       // destructor
+       delete fTextDQMShifterInfo;
+       delete fMax;
+}
+//____________________________________________________________________________
+AliACORDEQAChecker::AliACORDEQAChecker(const AliACORDEQAChecker& qac) :
+AliQACheckerBase(qac.GetName(), qac.GetTitle()),
+fTextDQMShifterInfo(new TPaveText(2.2,0.90,2.9,0.99,"T")),
+fMax(static_cast<TLine*>(qac.fMax->Clone()))
+{
+       //
+}
+//____________________________________________________________________________
+AliACORDEQAChecker& AliACORDEQAChecker::operator = (const AliACORDEQAChecker &qac)
+{
+       // use cpy option from constructor
+       if (&qac==this) return *this;
+       new (this) AliACORDEQAChecker(qac);
+       return *this;
+}
+//____________________________________________________________________________
+void AliACORDEQAChecker::Check(Double_t * test, AliQAv1::ALITASK_t /*index*/)
 {
-  Double_t * rv = new Double_t[AliRecoParam::kNSpecies] ; 
   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) 
-    rv[specie] = 0.0 ; 
-  return rv ;  
+    test[specie] = 0.0 ; 
 }
-
-//__________________________________________________________________
-Double_t * AliACORDEQAChecker::Check(AliQAv1::ALITASK_t /*index*/, TObjArray ** list)
+//____________________________________________________________________________
+void AliACORDEQAChecker::Check(Double_t * test, AliQAv1::ALITASK_t /*index*/, TObjArray ** list, const AliDetectorRecoParam * /*recoParam*/)
 {
-
-       // We added one check to the ACORDE's QA histograms:
-       // 1.- We check if they are empty
-       // we check for the reference histogram to start the QAChecker. If not QAref object
-       // is found, we check that the number of hits per channel is not so far from
-       // the maximum number of hits.
-  Double_t * test = new Double_t[AliRecoParam::kNSpecies] ; 
-  Int_t * count   = new Int_t[AliRecoParam::kNSpecies] ; 
-  Double_t * acoTest = new Double_t[AliRecoParam::kNSpecies];
-  Double_t acoHitsNorm = 0;
- Double_t * acoRefTest = new Double_t[AliRecoParam::kNSpecies];
-
-       // Look at the QAref data for ACORDE
-
-       char * acoOCDBDir = Form("ACORDE/%s/%s",AliQAv1::GetRefOCDBDirName(),AliQAv1::GetRefDataDirName());
-       AliCDBEntry *acoQARefDir = AliQAManager::QAManager()->Get(acoOCDBDir);
-
-
-  for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
-    test[specie]    = 0.0 ; 
-    count[specie] = 0 ; 
-       acoTest[specie] = 0.0;
-  }
+// Close version to the final one for the ACORDE QA Checker
   
-  for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
-    if (list[specie]->GetEntries() == 0){  
-      test[specie] = 1. ; // nothing to check
-       acoTest[specie] = 1.;
-    }
-    else {
-      TIter next(list[specie]) ; 
-      TH1 * hdata ;
-      while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
-        if (hdata) { 
-          Double_t rv = 0.0 ; 
-          if(hdata->GetEntries()>0)rv=1; 
-          AliInfo(Form("%s -> %f", hdata->GetName(), rv)) ; 
-          count[specie]++ ; 
-          test[specie] += rv ; 
-       
-
-       // here we implement the second version for ACORDEQAChecker
-       // by the moment we only compare that the hits in every ACORDE's channel
-       // are close and > 0 
-       for (Int_t i=0;i<60;i++)
+// Loop over the run species (for specie!= cosmic by now we set QA to INFO) 
+  
+  for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) 
+  {
+       test[specie] = 1.0;
+       if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) continue ; 
+       if (list[specie]->GetEntries() == 0) test[specie] = 1.; // Nothing to check
+       else 
        {
-               acoHitsNorm =  hdata->GetBinContent(i)/hdata->GetMaximum();
-               if  (acoQARefDir)
+               TIter next(list[specie]) ; 
+               TH1 * hdata ; // Data created by the AliACORDEQADataMakerXXX (Sim/Rec)
+               while ( (hdata = dynamic_cast<TH1 *>(next())) ) 
                {
-       //              AliWarning("Using the QA Reference data for ACORDE !!!");
-                       test[specie] = CheckAcordeRefHits(list[specie],(TObjArray *)acoQARefDir->GetObject());
-                       if ((test[specie] = 0.86) || (acoHitsNorm>0.50)) 
-                       {
-                               acoRefTest[specie]=0.78;//printf("testMario: %f\n",acoRefTest[specie]);printf("histo:%f\n",hdata->GetMaximum());
-                       }
-               }else{
-       //      AliWarning("Using the inner ACORDE QA Checker !!!");
-               if ( (acoHitsNorm>0.40) && (acoHitsNorm<=1) ) acoTest[specie] = 0.75;
-               if ( (acoHitsNorm>0.0020) && (acoHitsNorm<=0.40) ) acoTest[specie] = 0.251;
-               if ( (acoHitsNorm>0.0) && (acoHitsNorm<=0.0020) ) acoTest[specie] = 0.0010;
-               if ( (acoHitsNorm>-1.0) && (acoHitsNorm<=0.0) ) acoTest[specie] = -0.5;
-               }
-       }
-        }
-        else{
-          AliError("Data type cannot be processed") ;
-        }
-      }
-      if (count[specie] != 0) { 
-        if (test[specie]==0) {
-         // AliWarning("Histograms are there, but they are all empty: setting flag to kWARNING");
-          test[specie] = 0.5;  //upper limit value to set kWARNING flag for a task
-        }
-        else {
-       if (acoQARefDir) test[specie] = acoRefTest[specie];
-       else{
-       test[specie] = acoTest[specie];//printf("testDyMa: %f\n",test[specie]);
-       }
-        }
-      }
-    }
-   // AliInfo(Form("Test Result = %f", test[specie])) ; 
-  }
-  return test ; 
+                       //if (hdata) 
+                       if (hdata->GetEntries()!=0)
+                       { 
+                               Double_t rv = 0.0 ; 
+                               if(hdata->GetEntries()>0) rv=1; 
+                               AliDebug(AliQAv1::GetQADebugLevel(), Form("%s -> %f", hdata->GetName(), rv)) ; 
+                               TString hdataName = hdata->GetName();
+                               //if (hdata->GetListOfFunctions()->GetEntries() == 0) continue;  
+                                               
+                               if (hdataName.Contains("fhACORDEStatusAMU_DQM"))
+                               {
+                                       if (!hdata->GetListOfFunctions()->Contains(fTextDQMShifterInfo))
+                                               hdata->GetListOfFunctions()->Add(fTextDQMShifterInfo);
+                                       if (!hdata->GetListOfFunctions()->Contains(fMax))
+                                               hdata->GetListOfFunctions()->Add(fMax);
+                               }
+                               // check with OCDB ref. data
+                               if ( (fRefOCDBSubDir[specie]) && (hdataName.Contains("ACORDEBitPattern")) ) 
+                               {
+                                       TH1 * href = NULL;
+                                       if (fRefSubDir) href = static_cast<TH1*>(fRefSubDir->Get(hdata->GetName()));
+                                       else if (fRefOCDBSubDir[specie]) href = static_cast<TH1*>(fRefOCDBSubDir[specie]->FindObject(hdata->GetName()));
+                                       test[specie] = CheckAcordeRefHits(href,hdata);
+                               }else
+                               {
+                                       if (hdataName.Contains("ACORDEBitPattern"))
+                                               // Here we use an inner QA Checher without the QAref data
+                                       {
+                                               Float_t acoDataMax = hdata->GetMaximum();
+                                               Int_t flagAcoQAChecker = 0;
+                                               if (acoDataMax == 0) continue;
+                                               for(Int_t i=0;i<60;i++)
+                                               {
+                                                       if ((hdata->GetBinContent(i)/acoDataMax) < 0.75) flagAcoQAChecker++; 
+                                               }
+                                               Double_t simpleFlag = 1.-flagAcoQAChecker/60.;
+                                               if ( (simpleFlag >= 0.90) && (simpleFlag <= 1.0) ) test[specie] = 0.75; // INFO
+                                               if ( (simpleFlag >= 0.70) && (simpleFlag < 0.90) ) test[specie] = 0.50; // WARNING
+                                               if ( (simpleFlag >= 0.25) && (simpleFlag < 0.70) ) test[specie] = 0.25; // ERROR
+                                               if ( (simpleFlag >= 0.0) && (simpleFlag < 0.25) )  test[specie] = -1.0; // FATAL
+                                               if (test[specie]>0) rv = test[specie];
+                                               else rv = 1.0;
+                                       }
+                                       if (hdataName.Contains("fhACORDEStatusAMU_DQM"))
+                                       {
+                                               Float_t goodModules = hdata->GetBinContent(1);
+                                               Float_t badModules = hdata->GetBinContent(2);
+                                               if (goodModules >= badModules)
+                                                       test[specie] = 0.75;
+                                               else test[specie] = 0.30;
+                                               
+                                       }        // check for DQM condition
+                                       rv = test[specie];
+                                       if (fTextDQMShifterInfo) // DQM-flag asignment
+                                       {
+                                               fTextDQMShifterInfo->Clear();
+                                               if (rv > 0.4){
+                                                       fTextDQMShifterInfo->SetFillColor(kGreen);
+                                                       fTextDQMShifterInfo->AddText("ACORDE: O.K.");
+                                               }else{
+                                                       fTextDQMShifterInfo->SetFillColor(kRed);
+                                                       fTextDQMShifterInfo->AddText("ACORDE: Not, O.K.");
+                                               }
+                                       } // end DQM flag asignment
+                               }// NO OCDB checking
+                       } // hdata with entries
+                       else AliError("Data type cannot be processed") ;
+               } // loop over histograms from QADataMaker
+       }  // end checker
+    if ( (specie == AliRecoParam::kHighMult) || (specie == AliRecoParam::kLowMult) || (specie == AliRecoParam::kCalib) ) test[specie] = 0.75;
+  } // end of species
 }
-Double_t AliACORDEQAChecker::CheckAcordeRefHits(TObjArray *AcordeList, TObjArray * /*AcordeRef */) const
+//____________________________________________________________________________
+Double_t AliACORDEQAChecker::CheckAcordeRefHits(const TH1 * href, const TH1 * hdata) const
 {
-       Double_t acoTest = 0;
-       TIter next(AcordeList);
-       TH1 *histo;
+       Double_t test = 0.;
+       Int_t flag=0;
        for (Int_t i=0;i<60;i++)
        {
-               while ( (histo = dynamic_cast<TH1 *>(next())) )
-               {       
-                 if (histo->GetMaximum() && ((histo->GetBinContent(i)/histo->GetMaximum())<1.0) ) acoTest = 0.86;
-//             if( histo->KolmogorovTest((TH1F *)AcordeRef->At(0))<0.8)  acoTest = 0.86;
-                       //printf("href:%f\n",histo->GetMaximum());
-               }
-       }       
-       return acoTest;
+               if (TMath::Abs(href->GetBinContent(i)-hdata->GetBinContent(i))>10) flag++;
+       }
+       if ((flag>50)&&(flag<=60)) test = -1.;
+       if ((flag>30)&&(flag<=50)) test = 0.25;
+       if ((flag>10)&&(flag<=30)) test = 0.5;
+       if ((flag>0)&&(flag<=10)) test = 0.75;
+       return test;
 }