]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FORWARD/analysis2/corrs/RunCopyELossFit.C
Added 2012 geom
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis2 / corrs / RunCopyELossFit.C
1 /** 
2  * 
3  * 
4  * @param d 
5  * @param r 
6  * 
7  * @return 
8  * 
9  * @ingroup pwg2_forward_scripts_corr
10  */
11 Color_t Color(UShort_t d, Char_t r ) const 
12
13   return ((d == 1 ? kRed : (d == 2 ? kGreen : kBlue))
14           + ((r == 'I' || r == 'i') ? 2 : -2));
15 }
16
17 /** 
18  * 
19  * 
20  * @param sys 
21  * @param cms 
22  * @param field 
23  * @param mc 
24  * @param path 
25  * 
26  * @ingroup pwg2_forward_scripts_corr
27  */
28 void
29 RunCopyELossFit(UShort_t sys, UShort_t cms, Short_t field, Bool_t mc=false,
30                 const Char_t* path=0)
31 {
32   RunCopyELossFit(sys == 1 ? "pp" : "PbPb", 
33                   cms, 
34                   field, 
35                   mc,
36                   path);
37 }
38
39 /** 
40  * 
41  * @param sys       Collision system 
42  * @param cms       Center of mass energy per nucleon in GeV
43  * @param field     Magnetic field 
44  * @param mc        Monte-carlo flag
45  * @param path      File path
46  * 
47  * 
48  * @ingroup pwg2_forward_scripts_corr
49  */
50 void
51 RunCopyELossFit(const char* sys, UShort_t cms, Short_t field, bool mc=false,
52                 const char* path=0)
53 {
54   gROOT->Macro("$ALICE_ROOT/PWG2/FORWARD/analysis2/scripts/LoadLibs.C");
55   gSystem->Load("libPWG2forward.so");
56
57   AliFMDAnaParameters* p = AliFMDAnaParameters::Instance();
58   p->SetRealData(!mc);
59   p->SetEnergy(Float_t(cms));
60   p->SetMagField(Float_t(field));
61   p->SetCollisionSystem(sys);
62   if (path) {
63     p->SetBackgroundPath(path);
64     p->SetEnergyPath(path);
65     p->SetEventSelectionPath(path);
66     p->SetSharingEfficiencyPath(path);
67   }
68   p->Init(true, AliFMDAnaParameters::kBackgroundCorrection|
69           AliFMDAnaParameters::kEnergyDistributions);
70   
71   Int_t    nVtx   = p->GetNvtxBins();
72   Double_t minVtx = -p->GetVtxCutZ();
73   Double_t maxVtx = -p->GetVtxCutZ();
74   Int_t    nEta   = p->GetNetaBins();
75   Double_t minEta = p->GetEtaMin();
76   Double_t maxEta = p->GetEtaMax();
77   
78   TAxis vtxAxis(nVtx, minVtx, maxVtx);
79   TAxis etaAxis(nEta, minEta, maxEta);
80   AliFMDCorrELossFit* m = new AliFMDCorrELossFit;
81   m->SetEtaAxis(etaAxis);
82   
83   for (UShort_t d = 1; d <= 3; d++) { 
84     UShort_t nQ = (d == 1 ? 1 : 2);
85     for (UShort_t q = 0; q < nQ; q++) { 
86       Char_t r = (q == 0 ? 'I' : 'O');
87       
88       for (UShort_t b = 1; b < nEta; b++) { 
89         Double_t eta = etaAxis.GetBinCenter(b);
90         Info("RunCopyELossFit", "FMD%d%c, bin %3d, eta %+8.4f", d, r, b, eta);
91         if (eta < minEta || eta > maxEta) continue;
92         TH1F* oldhist = p->GetEnergyDistribution(d,r,eta);
93         if (!oldhist) {
94           Warning("RunCopyELossFit",
95                   "Didn't find energy distribution for FMD%d%c, eta bin %3d", 
96                   d, r, b);
97           continue;
98         }
99         TF1*  oldfunc = oldhist->GetFunction("FMDfitFunc");
100         if (!oldfunc) {
101           // Warning("RunCopyELossFit",
102           // "Didn't find energy fit for FMD%d%c, eta bin %3d", 
103           // d, r, b);
104           continue;
105         }
106         Double_t chi2   = oldfunc->GetChisquare();
107         UShort_t nu     = oldfunc->GetNDF();
108         Double_t c      = oldfunc->GetParameter(0);
109         Double_t delta  = oldfunc->GetParameter(1);
110         Double_t xi     = oldfunc->GetParameter(2);
111         Double_t a2 = 0, a3 = 0, ea2 = 0, ea3 = 0;
112         if (oldfunc->GetNpar() >= 4) {
113           a2     = oldfunc->GetParameter(3);
114           ea2    = oldfunc->GetParError(3);
115         }
116         if (oldfunc->GetNpar() >= 5) { 
117           a3     = oldfunc->GetParameter(4);
118           ea3    = oldfunc->GetParError(4);
119         }
120         Int_t    n      = (a3 < 1e-5 ? (a2 < 1e-5 ? 1 : 2) : 3);
121         Double_t ec     = oldfunc->GetParError(0);
122         Double_t edelta = oldfunc->GetParError(1);
123         Double_t exi    = oldfunc->GetParError(2);
124         delta           = delta - xi * -0.22278298;
125         edelta          = TMath::Sqrt(TMath::Power(edelta,2) - 
126                                       TMath::Power(0.22278298*exi,2));
127         Info("RunCopyELossFit", "FMD%d%c etabin=%3d: n=%d\n" 
128              "  C=%f+/-%f, Delta=%f+/-%f, xi=%f+/-%f, a2=%f+/-%f, a3=%f+/-%f",
129              d, r, b, n, c, ec, delta, edelta, xi, exi, a2, ea2, a3, ea3);
130         Double_t a[]  = { a2, a3, -9999 };
131         Double_t ea[] = { ea2, ea3, -9999 };
132
133         AliFMDCorrELossFit::ELossFit* fit = new 
134           AliFMDCorrELossFit::ELossFit(0, n, 
135                                        chi2,   nu, 
136                                        c,      ec, 
137                                        delta,  edelta, 
138                                        xi,     exi, 
139                                        0,      0, 
140                                        0,      0, 
141                                        a,      ea);
142         fit->CalculateQuality();
143         fit->fBin = b;
144         fit->fDet = d;
145         fit->fRing = r;
146         fit->Print();
147         
148         Info("RunCopyELossFit", "Setting fit FMD%d%c etabin=%3d: %p", 
149              d, r, b, fit);
150         m->SetFit(d, r, b, fit);
151       }
152     }
153   }
154   UShort_t isys = AliForwardUtil::ParseCollisionSystem(sys);
155   AliForwardCorrectionManager& mgr = AliForwardCorrectionManager::Instance();
156   TString fname(mgr.GetFileName(AliForwardCorrectionManager::kELossFits,
157                                 isys, cms, field, mc));
158   TFile* output = TFile::Open(fname.Data(), "RECREATE");
159   if (!output) { 
160     Warning("RunCopyELossFit", "Failed to open output file %s", fname.Data());
161     return kFALSE;
162   }
163   m->Write(mgr.GetObjectName(AliForwardCorrectionManager::kELossFits));
164   output->Write();
165   output->Close();
166   Info("RunCopyELossFit", 
167        "File %s created.  It should be copied to %s and stored in SVN",
168        fname.Data(), mgr.GetFileDir(AliForwardCorrectionManager::kELossFits));
169   
170 }
171 //
172 // EOF
173 //