]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Adding correction for FMD dead channels
authorhdalsgaa <hdalsgaa@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 1 Oct 2010 10:31:54 +0000 (10:31 +0000)
committerhdalsgaa <hdalsgaa@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 1 Oct 2010 10:31:54 +0000 (10:31 +0000)
PWG2/FORWARD/analysis/AliFMDAnaCalibBackgroundCorrection.cxx
PWG2/FORWARD/analysis/AliFMDAnaCalibBackgroundCorrection.h
PWG2/FORWARD/analysis/AliFMDAnaParameters.cxx
PWG2/FORWARD/analysis/AliFMDAnaParameters.h
PWG2/FORWARD/analysis/AliFMDAnalysisTaskBackgroundCorrection.cxx
PWG2/FORWARD/analysis/AliFMDAnalysisTaskDensity.cxx
PWG2/FORWARD/analysis/AliFMDAnalysisTaskGenerateCorrection.cxx

index 0b6a4fb31c7820d2f393f139398b2176b339fd2b..52bc92e799b77c5010fe15783ab49ee0b8be2cc3 100644 (file)
@@ -107,6 +107,21 @@ void AliFMDAnaCalibBackgroundCorrection::SetSPDDeadCorrection(Int_t vtxbin,
   hCorrection->SetName(Form("hSPDDeadCorrection_vtx%d",vtxbin));
   fListOfDoubleHitCorrection.Add(hCorrection);    
 }
+
+//____________________________________________________________________
+TH1F* AliFMDAnaCalibBackgroundCorrection::GetFMDDeadCorrection(Int_t  vtxbin) {
+  
+  TH1F* hCorrection    = (TH1F*)fListOfDoubleHitCorrection.FindObject(Form("hFMDDeadCorrection_vtx%d",vtxbin));
+  return hCorrection;
+}
+
+//____________________________________________________________________
+void AliFMDAnaCalibBackgroundCorrection::SetFMDDeadCorrection(Int_t vtxbin, 
+                                                             TH1F* hCorrection) {
+  hCorrection->SetName(Form("hFMDDeadCorrection_vtx%d",vtxbin));
+  fListOfDoubleHitCorrection.Add(hCorrection);    
+}
+
 //____________________________________________________________________
 void AliFMDAnaCalibBackgroundCorrection::SetRefAxis(TAxis* axis) {
   
index a4664b1e49d35e98e11dab25faee2837df6067ba..34c0a65d316c64d5b6609db09c03fdbcca05efc0 100644 (file)
@@ -33,6 +33,8 @@ class AliFMDAnaCalibBackgroundCorrection : public TObject
   void    SetDoubleHitCorrection(Int_t det, Char_t ring, TH1F* hCorrection);
   TH1F*   GetSPDDeadCorrection(Int_t vtxbin);
   void    SetSPDDeadCorrection(Int_t vtxbin, TH1F* hCorrection);
+  TH1F*   GetFMDDeadCorrection(Int_t vtxbin);
+  void    SetFMDDeadCorrection(Int_t vtxbin, TH1F* hCorrection);
   void    SetRefAxis(TAxis* axis);
   Int_t   GetNvtxBins();
   Float_t GetVtxCutZ();
index 1f1586ac38483473cf29a4582cb5787f6e3d7595..bbc34137573c4702d687a0e0d910a1364f23cf11 100644 (file)
@@ -283,7 +283,7 @@ void AliFMDAnaParameters::SetParametersFromESD(AliESDEvent* esd) {
   
   if(TMath::Abs(energy - 450.) < 1)  fEnergy = k900;
   if(TMath::Abs(energy - 3500.) < 1) fEnergy = k7000;
-  if(TMath::Abs(energy - 1380.) < 20) fEnergy = k1380;
+  if(TMath::Abs(energy - 2750.) < 20) fEnergy = k2750;
   
   
   if(TMath::Abs(magfield - 30000.) < 10 ) fMagField = k5G;
@@ -312,6 +312,8 @@ void AliFMDAnaParameters::PrintStatus() const
     energystring.Form("10000 GeV"); break;
   case k14000:
     energystring.Form("14000 GeV"); break;
+  case k2750:
+    energystring.Form("2750 GeV"); break;
   default:
     energystring.Form("invalid energy"); break;
   }
@@ -334,8 +336,10 @@ void AliFMDAnaParameters::PrintStatus() const
     magstring.Form("5 kGaus");   break;
   case k0G:
     magstring.Form("0 kGaus");   break;
+  case k5Gnegative:
+    magstring.Form("-5 kGaus");   break;
   default:
-    magstring.Form("invalid mag field"); break;
+    magstring.Form("invalid mag field %d", fMagField); break;
   }
   TString collsystemstring;
   switch(fSpecies) {
@@ -598,6 +602,17 @@ TH1F* AliFMDAnaParameters::GetSPDDeadCorrection(Int_t vtxbin) {
   return fBackground->GetSPDDeadCorrection(vtxbin);
 }
 //_____________________________________________________________________
+TH1F* AliFMDAnaParameters::GetFMDDeadCorrection(Int_t vtxbin) {
+                                               
+  //Get correction for several hits in strips for p+p
+  if(!fIsInit) {
+    AliWarning("Not initialized yet. Call Init() to remedy");
+    return 0;
+  }
+  
+  return fBackground->GetFMDDeadCorrection(vtxbin);
+}
+//_____________________________________________________________________
 Float_t AliFMDAnaParameters::GetEventSelectionEfficiency(Int_t vtxbin) {
   //Get event selection efficiency object
   if(!fIsInit) {
index 7292504914ef6fa69411c90009418af4cc486933..1717d55214913c910db0918d71928bc7868d24b9 100644 (file)
@@ -66,7 +66,7 @@ public:
   
   enum Trigger { kMB1 = 0, kMB2, kSPDFASTOR, kNOCTP, kEMPTY , kNSD};
   
-  enum Energy { k900 , k10000, k14000 , k7000, k2400, k5500, k1380};
+  enum Energy { k900 , k10000, k14000 , k7000, k2400, k5500, k2750};
   
   enum MagField {k0G, k5G, k5Gnegative};
   
@@ -101,6 +101,7 @@ public:
   TH2F* GetBackgroundCorrectionNSD(Int_t det, Char_t ring, Int_t vtxbin);
   TH1F* GetDoubleHitCorrection(Int_t det, Char_t ring);
   TH1F* GetSPDDeadCorrection(Int_t vtxbin);
+  TH1F* GetFMDDeadCorrection(Int_t vtxbin);
   
   TH1F* GetSharingEfficiency(Int_t det, Char_t ring, Int_t vtxbin);
   TH1F* GetSharingEfficiencyTrVtx(Int_t det, Char_t ring, Int_t vtxbin);
index 5f28a4ca156cdef0b0b1c0dd5c5ec46e0827fd73..35b3528b5e56cd9adda65699b353966850544a85 100644 (file)
@@ -256,13 +256,13 @@ void AliFMDAnalysisTaskBackgroundCorrection::Exec(Option_t */*option*/)
        TH1F* hSharingEff = pars->GetSharingEfficiencyTrVtx(det,ringChar,vtxbin); 
        TH1F* hSharingEffTrVtx = pars->GetSharingEfficiencyTrVtx(det,ringChar,vtxbin);  
       
-       for(Int_t nx=1; nx<hMult->GetNbinsX(); nx++) {
+       for(Int_t nx=1; nx<=hMult->GetNbinsX(); nx++) {
          Float_t correction      = hSharingEff->GetBinContent(nx);
          Float_t correctionTrVtx = hSharingEffTrVtx->GetBinContent(nx);
          //FIXME : This should be for NSD events
          Float_t correctionNSD   = hSharingEff->GetBinContent(nx);
          
-         for(Int_t ny=1; ny<hMult->GetNbinsY(); ny++) {
+         for(Int_t ny=1; ny<=hMult->GetNbinsY(); ny++) {
            
            if(correction != 0){
              hMult->SetBinContent(nx,ny,hMult->GetBinContent(nx,ny)/correction);
index 45345504e8dff4545cef672098c29ad684a0249b..61867d720c6c2749e16f55cfea562334f2064b2c 100644 (file)
@@ -97,6 +97,7 @@ void AliFMDAnalysisTaskDensity::CreateOutputObjects()
                              hBg->GetXaxis()->GetXmin(),
                              hBg->GetXaxis()->GetXmax(),
                              nSec, 0, 2*TMath::Pi());
+           hMult->Sumw2();
            
            fOutputList->Add(hMult);
          }
@@ -235,6 +236,13 @@ void AliFMDAnalysisTaskDensity::Exec(Option_t */*option*/)
            
          }
          
+         //Dead channel correction:
+         TH1F* hFMDDeadCorrection = pars->GetFMDDeadCorrection(vtxbin);
+         if(hFMDDeadCorrection)
+           if(hFMDDeadCorrection->GetBinContent(hFMDDeadCorrection->FindBin(eta)) != 0)
+             correction = correction*hFMDDeadCorrection->GetBinContent(hFMDDeadCorrection->FindBin(eta));
+         
+         
          //std::cout<<det<<"   "<<ring<<"    "<<sec<<"    "<<strip<<"    "<<vertex[2]<<"   "<<eta<<std::endl;
          
          //Float_t signal = hSignalDist->Integral(hSignalDist->FindBin(0.5),hSignalDist->FindBin(2)) ;//pars->GetConstant(det,ring,eta);
index 0d3c80a753d43618206c29f92b716f1f07e0c9eb..bde27afb1411d2c9b0f80d0068b091dad7d584d9 100644 (file)
@@ -84,6 +84,19 @@ void AliFMDAnalysisTaskGenerateCorrection::UserCreateOutputObjects()
     hSPDhits->Sumw2();
     fListOfHits.Add(hSPDhits);
     
+    TH1F* hReadChannels  = new TH1F(Form("hFMDReadChannels_vtx%d",v),
+                                   Form("hFMDReadChannels_vtx%d",v),
+                                   fNbinsEta,-4,6);
+    hReadChannels->Sumw2();
+    fListOfHits.Add(hReadChannels);
+    TH1F* hAllChannels  = new TH1F(Form("hFMDAllChannels_vtx%d",v),
+                                  Form("hFMDAllChannels_vtx%d",v),
+                                  fNbinsEta,-4,6);
+    hAllChannels->Sumw2();
+    fListOfHits.Add(hAllChannels);
+    
+    
+    
     for(Int_t iring = 0; iring<2;iring++) {
       Char_t ringChar = (iring == 0 ? 'I' : 'O');
       Int_t nSec = (iring == 1 ? 40 : 20);
@@ -237,6 +250,10 @@ void AliFMDAnalysisTaskGenerateCorrection::UserExec(Option_t */*option*/)
   Double_t esdvertex[3];
   Bool_t vtxStatus =  pars->GetVertex(esdevent,esdvertex);
   
+  /* Double_t deltaEsd           = 2*fZvtxCut/fNvtxBins;
+  Double_t vertexBinDoubleEsd = (esdvertex[2] + fZvtxCut) / deltaEsd;
+  Int_t    vertexBinEsd       = (Int_t)vertexBinDoubleEsd;*/
+  
   AliMCParticle* particle = 0;
   AliStack* stack = mcevent->Stack();
   
@@ -326,8 +343,32 @@ void AliFMDAnalysisTaskGenerateCorrection::UserExec(Option_t */*option*/)
   hEventsAll->Fill(vertexBin);
   if(nsd) hEventsAllNSD->Fill(vertexBin);
   
-  //if(!vtxFound || !isTriggered) return; //Not important for FMD but crucial for SPD since maps are done from ESD
-  
+  // FMD dead channels
+  
+  TH1F* hFMDReadChannels  = (TH1F*)fListOfHits.FindObject(Form("hFMDReadChannels_vtx%d",vertexBin));
+  TH1F* hFMDAllChannels   = (TH1F*)fListOfHits.FindObject(Form("hFMDAllChannels_vtx%d",vertexBin));
+  
+  AliESDFMD* fmd = esdevent->GetFMDData();
+  
+  for(UShort_t d=1;d<=3;d++) {
+    Int_t nRings = (d==1 ? 1 : 2);
+    for (UShort_t ir = 0; ir < nRings; ir++) {
+      Char_t   r = (ir == 0 ? 'I' : 'O');
+      UShort_t nsec = (ir == 0 ? 20  : 40);
+      UShort_t nstr = (ir == 0 ? 512 : 256);
+      for(UShort_t s =0; s < nsec;  s++) {
+       for(UShort_t str = 0; str < nstr; str++) {
+         Float_t mult = fmd->Multiplicity(d,r,s,str);
+         
+         hFMDAllChannels->Fill(pars->GetEtaFromStrip(d,r,s,str,esdvertex[2]));
+         
+         if(mult != AliESDFMD::kInvalidMult)
+           hFMDReadChannels->Fill(pars->GetEtaFromStrip(d,r,s,str,esdvertex[2]));
+       }
+      }
+    }
+  }
+  // End, FMD dead channels
 
   for(Int_t i = 0 ;i<nTracks;i++) {
     particle = (AliMCParticle*) mcevent->GetTrack(i);
@@ -512,6 +553,17 @@ void AliFMDAnalysisTaskGenerateCorrection::GenerateCorrection() {
     fEventSelectionEff->SetCorrection("NSD",v,'I',hAnalysedNSDInner);
     fEventSelectionEff->SetCorrection("NSD",v,'O',hAnalysedNSDOuter);
     
+    //FMD dead channels
+    TH1F* hFMDReadChannels   = (TH1F*)fListOfHits.FindObject(Form("hFMDReadChannels_vtx%d",v));
+    
+    TH1F* hFMDAllChannels    = (TH1F*)fListOfHits.FindObject(Form("hFMDAllChannels_vtx%d",v));
+    if(hFMDReadChannels) {
+      TH1F* hFMDDeadCorrection = (TH1F*)hFMDReadChannels->Clone("hFMDDeadCorrection");
+      hFMDDeadCorrection->Divide(hFMDAllChannels);
+      fBackground->SetFMDDeadCorrection(v,hFMDDeadCorrection);
+      fListOfCorrection.Add(hFMDDeadCorrection);
+    }
+    else AliWarning("No Dead Channel Correction generated");
     
   }
   
@@ -600,7 +652,7 @@ void AliFMDAnalysisTaskGenerateCorrection::GenerateCorrection() {
     hCorrection->SetTitle(hCorrection->GetName());
     fListOfCorrection.Add(hCorrection);
     hCorrection->Divide(hPrimary);
-    fBackground->SetBgCorrection(0,'Q',vertexBin,hCorrection);
+
     
     TH1F* hAlive = new TH1F(Form("hAliveSPD_vtxbin%d",vertexBin),Form("hAliveSPD_vtxbin%d",vertexBin),hSPDMult->GetNbinsX(),hSPDMult->GetXaxis()->GetXmin(), hSPDMult->GetXaxis()->GetXmax());
     TH1F* hPresent = new TH1F(Form("hPresentSPD_vtxbin%d",vertexBin),Form("hPresentSPD_vtxbin%d",vertexBin),hSPDMult->GetNbinsX(),hSPDMult->GetXaxis()->GetXmin(), hSPDMult->GetXaxis()->GetXmax());
@@ -630,11 +682,9 @@ void AliFMDAnalysisTaskGenerateCorrection::GenerateCorrection() {
     hDeadCorrection->Divide(hPresent);
     fBackground->SetSPDDeadCorrection(vertexBin,hDeadCorrection);
     fListOfCorrection.Add(hDeadCorrection);
-  }
-  
-  
-  
   
+    fBackground->SetBgCorrection(0,'Q',vertexBin,hCorrection);
+  }
   
   
   TAxis refAxis(fNvtxBins,-1*fZvtxCut,fZvtxCut);
@@ -701,7 +751,13 @@ void AliFMDAnalysisTaskGenerateCorrection::ReadFromFile(const Char_t* filename,
   for(Int_t v=0; v<fNvtxBins;v++) {
     TH2F* hSPDHits          = (TH2F*)listOfHits->FindObject(Form("hSPDhits_vtx%d", v));   
     fListOfHits.Add(hSPDHits);
+    TH1F* hFMDReadChannels  = (TH1F*)listOfHits->FindObject(Form("hFMDReadChannels_vtx%d",v));
+    TH1F* hFMDAllChannels   = (TH1F*)listOfHits->FindObject(Form("hFMDAllChannels_vtx%d",v));
     
+    if(hFMDReadChannels && hFMDAllChannels) {
+      fListOfHits.Add(hFMDReadChannels);
+      fListOfHits.Add(hFMDAllChannels);
+    }
     for(Int_t iring = 0; iring<2;iring++) {
       Char_t ringChar = (iring == 0 ? 'I' : 'O');