]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/corrs/RunCopyMergeEff.C
Fixed references from PWG2 -> PWGLF - very efficiently done using ETags.
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / corrs / RunCopyMergeEff.C
1 /** 
2  * 
3  * 
4  * @param d 
5  * @param r 
6  * 
7  * @return 
8  * @ingroup pwglf_forward_scripts_corr
9  */
10 Color_t Color(UShort_t d, Char_t r ) const 
11
12   return ((d == 1 ? kRed : (d == 2 ? kGreen : kBlue))
13           + ((r == 'I' || r == 'i') ? 2 : -2));
14 }
15
16 /** 
17  * 
18  * 
19  * @param sys 
20  * @param cms 
21  * @param field 
22  * @param path 
23  * @ingroup pwglf_forward_scripts_corr
24  */
25 void
26 RunCopyMergeEff(UShort_t sys, UShort_t cms, Short_t field, const Char_t* path=0)
27 {
28   RunCopyMergeEff(sys == 1 ? "pp" : "PbPb", 
29                   cms, 
30                   field, 
31                   path);
32 }
33
34 /** 
35  * 
36  * @param sys       Collision system 
37  * @param cms       Center of mass energy per nucleon in GeV
38  * @param field     Magnetic field 
39  * @param path      File path
40  * 
41  * @ingroup pwglf_forward_scripts_corr
42  */
43 void
44 RunCopyMergeEff(const char* sys, UShort_t cms, 
45                 Short_t field, const char* path=0)
46 {
47   gROOT->Macro("$ALICE_ROOT/PWGLF/FORWARD/analysis2/scripts/LoadLibs.C");
48   gSystem->Load("libPWGLFforward.so");
49
50   AliFMDAnaParameters* p = AliFMDAnaParameters::Instance();
51   p->SetEnergy(Float_t(cms));
52   p->SetMagField(Float_t(field));
53   p->SetCollisionSystem(sys);
54   if (path) {
55     p->SetBackgroundPath(path);
56     p->SetEnergyPath(path);
57     p->SetEventSelectionPath(path);
58     p->SetSharingEfficiencyPath(path);
59   }
60   p->Init(true, AliFMDAnaParameters::kBackgroundCorrection|
61           AliFMDAnaParameters::kSharingEfficiency);
62   
63   Int_t    nVtx   = p->GetNvtxBins();
64   Double_t minVtx = -p->GetVtxCutZ();
65   Double_t maxVtx = -p->GetVtxCutZ();
66   Int_t    nEta   = p->GetNetaBins();
67   Double_t minEta = p->GetEtaMin();
68   Double_t maxEta = p->GetEtaMax();
69   
70   TAxis vtxAxis(nVtx, minVtx, maxVtx);
71   AliFMDCorrMergingEfficiency* m = new AliFMDCorrMergingEfficiency;
72   m->SetVertexAxis(nVtx, minVtx, maxVtx);
73   
74   for (UShort_t d = 1; d <= 3; d++) { 
75     UShort_t nQ = (d == 1 ? 1 : 2);
76     for (UShort_t q = 0; q < nQ; q++) { 
77       Char_t r = (q == 0 ? 'I' : 'O');
78       
79       for (UShort_t b = 1; b <= nVtx; b++) { 
80         TH1F* oldmap = p->GetSharingEfficiency(d, r, b-1);
81         if (!oldmap) {
82           Warning("RunCopyMergeEff",
83                   "Didn't find secondary map correction "
84                   "for FMD%d%c, vertex bin %3d", d, r, b);
85           continue;
86         }
87         
88         TH1D* newmap = new TH1D(Form("FMD%d%c_vtxbin%03d", d, r, b), 
89                                 Form("Merging efficiency for FMD%d%c "
90                                      "in vertex bin %d [%+8.4f,%+8.4f]", 
91                                      d, r, b, vtxAxis.GetBinLowEdge(b), 
92                                      vtxAxis.GetBinUpEdge(b)), 
93                                 nEta, minEta, maxEta);
94         newmap->SetXTitle("#eta");
95         newmap->SetYTitle("dN_{ch,i,incl}/d#eta / #sum_{i} N_{ch,i,FMD}");
96         newmap->SetDirectory(0);
97         newmap->SetStats(0);
98         newmap->Sumw2();
99         newmap->Add(oldmap);
100         
101         Info("RunCopyMergeEff",
102              "Copying %s to %s", oldmap->GetName(), newmap->GetName());
103         
104         m->SetCorrection(d, r, b, newmap);
105       }
106     }
107   }
108   UShort_t isys = AliForwardUtil::ParseCollisionSystem(sys);
109   AliForwardCorrectionManager& mgr = AliForwardCorrectionManager::Instance();
110   TString fname(mgr.GetFileName(AliForwardCorrectionManager::kMergingEfficiency,
111                                 isys, cms, field, false));
112   TFile* output = TFile::Open(fname.Data(), "RECREATE");
113   if (!output) { 
114     Warning("Run", "Failed to open output file %s", fname.Data());
115     return kFALSE;
116   }
117   m->Write(mgr.GetObjectName(AliForwardCorrectionManager::kMergingEfficiency));
118   output->Write();
119   output->Close();
120   Info("Run", "File %s created.  It should be copied to %s and stored in SVN",
121        fname.Data(),
122        mgr.GetFileDir(AliForwardCorrectionManager::kMergingEfficiency));
123   
124 }
125 //
126 // EOF
127 //