ROOT logo
#if !defined(__CINT__) || defined(__MAKECINT__)
#include <TSystem.h>
#include <TFile.h>
#include <TStyle.h>
#include <TTree.h>
#include <TClonesArray.h>
#include <TObjArray.h>
#include <TF1.h>
#include <TH1F.h>
#include <TH2F.h>
#include <TProfile.h>
#include <TCanvas.h>

#include <AliHMPIDParam.h>
#include <AliHMPIDHit.h>
#include <AliHMPIDCluster.h>
#include <AliHMPIDDigit.h>
#endif

Int_t nEntries = 0;

TObjArray *CreateContainer(const char *classname,TTree *pTree);
void Hits(Int_t mode,TTree *pTree=0x0);
void Digs(Int_t mode, TTree *pTree=0x0);
void Clus(Int_t mode, TTree *pTree=0x0);
void Summary();
  
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
TObjArray *CreateContainer(const char *classname,TTree *pTree)
{
  TObjArray *pOA=new TObjArray(AliHMPIDParam::kMaxCh+1);
  for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
    TClonesArray *pCA=new TClonesArray(classname);
    pOA->AddAt(pCA,iCh);    
  }
  
  pOA->SetOwner(kTRUE);  
  
  for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
    pTree->SetBranchAddress(Form("HMPID%i",iCh),&(*pOA)[iCh]);
  }
  return pOA;
}
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
TH1F *hHitQdc,*hHitQdcCh[AliHMPIDParam::kMaxCh+1]; TH2F *hHitMap[AliHMPIDParam::kMaxCh+1];  
Double_t fHitMean[AliHMPIDParam::kMaxCh+1],fHitErr[AliHMPIDParam::kMaxCh+1];
void Hits(Int_t mode,TTree *pTree)
{
  switch(mode){
    case 1:
      hHitQdc=new TH1F("HitQdc","Hit Qdc all chamber;QDC",4000,0,4000);
  for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
        hHitMap[iCh]=new TH2F(Form("HitMap%i",iCh),Form("Ch%i;x_{Hit};y_{Hit}",iCh),
            Int_t(AliHMPIDParam::SizeAllX()/AliHMPIDParam::SizePadX()),0,AliHMPIDParam::SizeAllX(),
            Int_t(AliHMPIDParam::SizeAllY()/AliHMPIDParam::SizePadY()),0,AliHMPIDParam::SizeAllY());
        hHitQdcCh[iCh]=new TH1F(Form("HMPID%i",iCh),Form("Charge for HMPID%i",iCh),4000,0,4000);
      }
      break;
    case 2:
      if(pTree==0) return;
      TClonesArray *pHits=new TClonesArray("AliHMPIDHit");  pTree->SetBranchAddress("HMPID",&pHits);  
      for(Int_t iEnt=0;iEnt<pTree->GetEntriesFast();iEnt++){//entries loop
        pTree->GetEntry(iEnt);
        for(Int_t iHit=0;iHit<pHits->GetEntriesFast();iHit++){//hits loop
          AliHMPIDHit *pHit = (AliHMPIDHit*)pHits->UncheckedAt(iHit);
          hHitMap[pHit->Ch()]->Fill(pHit->LorsX(),pHit->LorsY());
          hHitQdc->Fill(pHit->Q());
          hHitQdcCh[pHit->Ch()]->Fill(pHit->Q());
        }//hits loop      
      }//entries loop
      delete pHits;
      break;
    case 3:
      TCanvas *c1=new TCanvas("HitCan","Hits",1280,800); c1->Divide(3,3);
  
      for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
        if(iCh==6) c1->cd(1); if(iCh==5) c1->cd(2);
        if(iCh==4) c1->cd(4); if(iCh==3) c1->cd(5); if(iCh==2) c1->cd(6);
                              if(iCh==1) c1->cd(8); if(iCh==0) c1->cd(9);
        gStyle->SetPalette(1);      
        hHitMap[iCh]->Draw("colz");
      }  
      for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
        fHitMean[iCh] = 0;
        fHitErr[iCh]  = 0;
        if((Int_t)hHitQdcCh[iCh]->GetEntries()<nEntries) continue;
        c1->cd(3);hHitQdcCh[iCh]->Fit("expo","Q");
        TF1 *funcfit = (TF1*)hHitQdcCh[iCh]->FindObject("expo");
        fHitMean[iCh] = funcfit->GetParameter(1);
        fHitErr[iCh]  = funcfit->GetParError(1);
      }
      TPad *pad = (TPad*)c1->cd(3); hHitQdc->SetFillColor(5); pad->SetLogy(); hHitQdc->Fit("expo","Q");
      break;
  }
}//Hits()
//+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
TH1F     *hDigQ;
TProfile *hDigHighQ,*hDigChEvt;
void Digs(Int_t mode, TTree *pTree)
{
  switch(mode){
    case 1:
      hDigHighQ=new TProfile("hDigHighQ","Highest charge in chamber  ",AliHMPIDParam::kMaxCh+1,AliHMPIDParam::kMinCh,AliHMPIDParam::kMaxCh+1);
      hDigChEvt=new TProfile("hDigChEvt","Chamber occupancy per event",AliHMPIDParam::kMaxCh+1,AliHMPIDParam::kMinCh,AliHMPIDParam::kMaxCh+1);
      hDigQ    =new TH1F("hDigQ        ","Charge of digits (ADC)     ",3000,0,3000);
      break;
    case 2:
      if(pTree==0) return;
      TObjArray *pLst=CreateContainer("AliHMPIDDigit",pTree); pTree->GetEntry(0);
      for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){//chambers loop
        TClonesArray *pDigs=(TClonesArray *)pLst->UncheckedAt(iCh);
        hDigChEvt->Fill(iCh,pDigs->GetEntriesFast()/(48.*80.*6.)*100.);
        Double_t highQ = 0;
        for(Int_t iDig=0;iDig<pDigs->GetEntriesFast();iDig++){//digits loop
          AliHMPIDDigit *pDig=(AliHMPIDDigit*)pDigs->UncheckedAt(iDig);
          hDigQ->Fill(pDig->Q());
          if(pDig->Q()>highQ) highQ = pDig->Q();
        }
        hDigHighQ->Fill(iCh,highQ);
      }
      delete pLst;
      break;
    case 3:
      TCanvas *c1=new TCanvas("DigQa","Digit Check",1280,800); c1->Divide(2,2);
      gStyle->SetOptFit(1);
      TPad *pad = (TPad*)c1->cd(1); pad->SetLogy(); hDigQ->Fit("expo","QN"); 
      c1->cd(2); hDigHighQ->Draw();
      c1->cd(3); hDigChEvt->Draw();
      break;
  }//switch
}//Dig()
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
TH1F *hCluEvt,*hCluChi2,*hCluFlg,*hCluSize,*hCluQ;
void Clus(Int_t mode, TTree *pTree)
{
  switch(mode){
    case 1:
      hCluEvt=new TH1F("CluPerEvt","Cluster multiplicity"   ,100,0,100);
      hCluChi2  =new TH1F("CluChi2"  ,"Chi2 "               ,1000,0,100);
      hCluFlg   =new TH1F("CluFlg"   ,"Cluster flag"        ,14,-1.5,12.5);                       hCluFlg->SetFillColor(5);
      hCluSize  =new TH1F("CluSize"  ,"Cluster size        ",100,0,100);
      hCluQ     =new TH1F("CluQ"     ,"Cluster charge (ADC)",1000,0,5000);
      break;
    case 2:      
      if(pTree==0) return;

      TObjArray *pLst=CreateContainer("AliHMPIDCluster",pTree); pTree->GetEntry(0);
      for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){//chambers loop
        TClonesArray *pClus=(TClonesArray *)pLst->UncheckedAt(iCh);
        hCluEvt->Fill(pClus->GetEntriesFast());
        for(Int_t iClu=0;iClu<pClus->GetEntriesFast();iClu++){//clusters loop
          AliHMPIDCluster *pClu=(AliHMPIDCluster*)pClus->UncheckedAt(iClu);
          hCluFlg->Fill(pClu->Status());
          hCluChi2->Fill(pClu->Chi2());
          hCluSize->Fill(pClu->Size());
          hCluQ->Fill(pClu->Q());
        }
      }
      delete pLst;           
      break;
    case 3:
      TCanvas *c1=new TCanvas("CluComCan","Clusters in common",1280,800); c1->Divide(3,3);
      c1->cd(1); hCluEvt->SetFillColor(5);      hCluEvt->Draw();
      c1->cd(2); hCluChi2->SetFillColor(5);     hCluChi2->Draw(); 
      c1->cd(3); hCluFlg->SetFillColor(5);      hCluFlg->Draw(); 
      c1->cd(4); hCluSize->SetFillColor(5);     hCluSize->Draw(); 
      TPad *pad = (TPad*)c1->cd(5); hCluQ->SetFillColor(5); pad->SetLogy(); hCluQ->Draw(); 
      break;
  }//switch
}//Clus()
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void Hqa()
{
  AliHMPIDParam::Instance();
  
  TFile *fh=0; if(gSystem->IsFileInIncludePath("HMPID.Hits.root"))      fh=TFile::Open("HMPID.Hits.root"     ,"read");if(fh->IsZombie()) fh=0;
  TFile *fd=0; if(gSystem->IsFileInIncludePath("HMPID.Digits.root"))    fd=TFile::Open("HMPID.Digits.root"   ,"read");if(fd->IsZombie()) fd=0;
  TFile *fc=0; if(gSystem->IsFileInIncludePath("HMPID.RecPoints.root")) fc=TFile::Open("HMPID.RecPoints.root","read");if(fc->IsZombie()) fc=0;
  if(fh==0 && fd==0 && fc==0){Printf("Nothing to do!"); return;}
  if(fh) Hits(1); if(fc) Clus(1);  if(fd) Digs(1);//book
  Int_t iEvt=0;
  while(1){
    TTree *th=0; if(fh) th=(TTree*)fh->Get(Form("Event%i/TreeH",iEvt));
    TTree *td=0; if(fd) td=(TTree*)fd->Get(Form("Event%i/TreeD",iEvt));
    TTree *tc=0; if(fc) tc=(TTree*)fc->Get(Form("Event%i/TreeR",iEvt));
    
    Hits(2,th); Digs(2,td); Clus(2,tc); //fill
    if(th==0 && td==0 && tc==0) break;
    iEvt++;
    if(!(iEvt%50)) Printf("Event %i processed",iEvt);
  }
  
  nEntries = iEvt;
  
  if(fd) Clus(3);//plot everything
  if(fc) Digs(3); 
  if(fh) Hits(3);
  Summary();
}
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void Summary()
{
  //info for hits...
  for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
    Printf(" #################### Summary for HMPID %i#################### ",iCh);
    //info for hits...
    Printf("-----Summary of Hits----- ");
    if(fHitMean[iCh]==0) {
      Printf("gain %5.2f +/- %5.2f",fHitMean[iCh],fHitErr[iCh]);
    } else {   
      Double_t gain = 1./TMath::Abs(fHitMean[iCh]);
      Double_t errgain = gain*gain*fHitErr[iCh];
      Printf("gain %5.2f +/- %5.2f",gain,errgain);
    }
    //info for digits...
    Printf("-----Summary of Digits-----");
    Printf(" Chamber %d with occupancy (%) %5.2f +/- %5.2f",iCh,hDigChEvt->GetBinContent(iCh+1),hDigChEvt->GetBinError(iCh+1));
    Printf(" Chamber %d with higest Q (ADC) %7.0f +/- %7.0f",iCh,hDigHighQ->GetBinContent(iCh+1),hDigHighQ->GetBinError(iCh+1));
  }
}
 Hqa.C:1
 Hqa.C:2
 Hqa.C:3
 Hqa.C:4
 Hqa.C:5
 Hqa.C:6
 Hqa.C:7
 Hqa.C:8
 Hqa.C:9
 Hqa.C:10
 Hqa.C:11
 Hqa.C:12
 Hqa.C:13
 Hqa.C:14
 Hqa.C:15
 Hqa.C:16
 Hqa.C:17
 Hqa.C:18
 Hqa.C:19
 Hqa.C:20
 Hqa.C:21
 Hqa.C:22
 Hqa.C:23
 Hqa.C:24
 Hqa.C:25
 Hqa.C:26
 Hqa.C:27
 Hqa.C:28
 Hqa.C:29
 Hqa.C:30
 Hqa.C:31
 Hqa.C:32
 Hqa.C:33
 Hqa.C:34
 Hqa.C:35
 Hqa.C:36
 Hqa.C:37
 Hqa.C:38
 Hqa.C:39
 Hqa.C:40
 Hqa.C:41
 Hqa.C:42
 Hqa.C:43
 Hqa.C:44
 Hqa.C:45
 Hqa.C:46
 Hqa.C:47
 Hqa.C:48
 Hqa.C:49
 Hqa.C:50
 Hqa.C:51
 Hqa.C:52
 Hqa.C:53
 Hqa.C:54
 Hqa.C:55
 Hqa.C:56
 Hqa.C:57
 Hqa.C:58
 Hqa.C:59
 Hqa.C:60
 Hqa.C:61
 Hqa.C:62
 Hqa.C:63
 Hqa.C:64
 Hqa.C:65
 Hqa.C:66
 Hqa.C:67
 Hqa.C:68
 Hqa.C:69
 Hqa.C:70
 Hqa.C:71
 Hqa.C:72
 Hqa.C:73
 Hqa.C:74
 Hqa.C:75
 Hqa.C:76
 Hqa.C:77
 Hqa.C:78
 Hqa.C:79
 Hqa.C:80
 Hqa.C:81
 Hqa.C:82
 Hqa.C:83
 Hqa.C:84
 Hqa.C:85
 Hqa.C:86
 Hqa.C:87
 Hqa.C:88
 Hqa.C:89
 Hqa.C:90
 Hqa.C:91
 Hqa.C:92
 Hqa.C:93
 Hqa.C:94
 Hqa.C:95
 Hqa.C:96
 Hqa.C:97
 Hqa.C:98
 Hqa.C:99
 Hqa.C:100
 Hqa.C:101
 Hqa.C:102
 Hqa.C:103
 Hqa.C:104
 Hqa.C:105
 Hqa.C:106
 Hqa.C:107
 Hqa.C:108
 Hqa.C:109
 Hqa.C:110
 Hqa.C:111
 Hqa.C:112
 Hqa.C:113
 Hqa.C:114
 Hqa.C:115
 Hqa.C:116
 Hqa.C:117
 Hqa.C:118
 Hqa.C:119
 Hqa.C:120
 Hqa.C:121
 Hqa.C:122
 Hqa.C:123
 Hqa.C:124
 Hqa.C:125
 Hqa.C:126
 Hqa.C:127
 Hqa.C:128
 Hqa.C:129
 Hqa.C:130
 Hqa.C:131
 Hqa.C:132
 Hqa.C:133
 Hqa.C:134
 Hqa.C:135
 Hqa.C:136
 Hqa.C:137
 Hqa.C:138
 Hqa.C:139
 Hqa.C:140
 Hqa.C:141
 Hqa.C:142
 Hqa.C:143
 Hqa.C:144
 Hqa.C:145
 Hqa.C:146
 Hqa.C:147
 Hqa.C:148
 Hqa.C:149
 Hqa.C:150
 Hqa.C:151
 Hqa.C:152
 Hqa.C:153
 Hqa.C:154
 Hqa.C:155
 Hqa.C:156
 Hqa.C:157
 Hqa.C:158
 Hqa.C:159
 Hqa.C:160
 Hqa.C:161
 Hqa.C:162
 Hqa.C:163
 Hqa.C:164
 Hqa.C:165
 Hqa.C:166
 Hqa.C:167
 Hqa.C:168
 Hqa.C:169
 Hqa.C:170
 Hqa.C:171
 Hqa.C:172
 Hqa.C:173
 Hqa.C:174
 Hqa.C:175
 Hqa.C:176
 Hqa.C:177
 Hqa.C:178
 Hqa.C:179
 Hqa.C:180
 Hqa.C:181
 Hqa.C:182
 Hqa.C:183
 Hqa.C:184
 Hqa.C:185
 Hqa.C:186
 Hqa.C:187
 Hqa.C:188
 Hqa.C:189
 Hqa.C:190
 Hqa.C:191
 Hqa.C:192
 Hqa.C:193
 Hqa.C:194
 Hqa.C:195
 Hqa.C:196
 Hqa.C:197
 Hqa.C:198
 Hqa.C:199
 Hqa.C:200
 Hqa.C:201
 Hqa.C:202
 Hqa.C:203
 Hqa.C:204
 Hqa.C:205
 Hqa.C:206
 Hqa.C:207
 Hqa.C:208
 Hqa.C:209
 Hqa.C:210
 Hqa.C:211
 Hqa.C:212
 Hqa.C:213
 Hqa.C:214
 Hqa.C:215
 Hqa.C:216
 Hqa.C:217
 Hqa.C:218
 Hqa.C:219
 Hqa.C:220
 Hqa.C:221