]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FORWARD/analysis2/corrs/RunCopySecMap.C
2dc75f122a4c551f929d720e3657f3d1b849375b
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis2 / corrs / RunCopySecMap.C
1 /** 
2  * 
3  * 
4  * @param d 
5  * @param r 
6  * 
7  * @return 
8  * 
9  * @ingroup pwg2_forward_analysis_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  * @param sys 
20  * @param cms 
21  * @param field 
22  * @param path 
23  * 
24  * @ingroup pwg2_forward_analysis_scripts_corr
25  */
26 void
27 RunCopySecMap(UShort_t sys, UShort_t cms, Short_t field, const Char_t* path=0)
28 {
29   RunCopySecMap(sys == 1 ? "pp" : "PbPb", 
30                 cms, 
31                 field, 
32                 path);
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  * 
42  * @ingroup pwg2_forward_analysis_scripts_corr
43  */
44 void
45 RunCopySecMap(const char* sys, UShort_t cms, Short_t field,
46               const Char_t* path=0)
47 {
48   gROOT->Macro("$ALICE_ROOT/PWG2/FORWARD/analysis2/scripts/LoadLibs.C");
49   gSystem->Load("libPWG2forward.so");
50
51   AliFMDAnaParameters* p = AliFMDAnaParameters::Instance();
52   p->SetEnergy(Float_t(cms));
53   p->SetMagField(Float_t(field));
54   p->SetCollisionSystem(sys);
55   if (path) {
56     std::cout<<"Setting path to "<<path<<std::endl;
57     p->SetBackgroundPath(path);
58     p->SetEnergyPath(path);
59     p->SetEventSelectionPath(path);
60     p->SetSharingEfficiencyPath(path);
61   }
62   p->Init(true, AliFMDAnaParameters::kBackgroundCorrection);
63  
64   Int_t    nVtx   = p->GetNvtxBins();
65   Double_t minVtx = -p->GetVtxCutZ();
66   Double_t maxVtx = p->GetVtxCutZ();
67   Int_t    nEta   = p->GetNetaBins();
68   Double_t minEta = p->GetEtaMin();
69   Double_t maxEta = p->GetEtaMax();
70
71   TAxis vtxAxis(nVtx, minVtx, maxVtx);
72   AliFMDCorrSecondaryMap* m = new AliFMDCorrSecondaryMap;
73   m->SetVertexAxis(nVtx, minVtx, maxVtx);
74   m->SetEtaAxis(nEta,minEta,maxEta);
75   AliFMDCorrDoubleHit* dh = new AliFMDCorrDoubleHit;
76
77   for (UShort_t d = 1; d <= 3; d++) { 
78     UShort_t nQ = (d == 1 ? 1 : 2);
79     for (UShort_t q = 0; q < nQ; q++) { 
80       Char_t r = (q == 0 ? 'I' : 'O');
81      
82       TH1F* old2hit = p->GetDoubleHitCorrection(d, r);
83       if (!old2hit) 
84         Warning("RunCopySecMap",
85                 "Didn't find double hit correction for FMD%d%c", d, r);
86       else { 
87         TH1D* new2hit = new TH1D(Form("FMD%d%c", d, r), 
88                                  Form("Double hit correction for FMD%d%c",d,r), 
89                                  nEta, minEta, maxEta);
90         new2hit->SetXTitle("#eta");
91         new2hit->SetYTitle("#sum_{i} N_{i,strips hit}(#eta)/"
92                            "#sum_{i} N_{i,total hits}(#eta)");
93         new2hit->SetFillColor(Color(d,r));
94         new2hit->SetFillStyle(3001);
95         new2hit->SetDirectory(0);
96         new2hit->SetStats(0);
97         new2hit->Sumw2();
98         new2hit->Add(old2hit);
99
100         Info("RunCopySecMap", 
101              "Copying %s to %s", old2hit->GetName(), new2hit->GetName());
102
103         dh->SetCorrection(d, r, new2hit);
104       }
105
106       for (UShort_t b = 1; b <= nVtx; b++) { 
107         TH2F* oldmap = p->GetBackgroundCorrection(d, r, b-1);
108         if (!oldmap) {
109           Warning("RunCopySecMap",
110                   "Didn't find secondary map correction "
111                   "for FMD%d%c, vertex bin %3d", d, r, b);
112           continue;
113         }
114        
115         TH2D* newmap = new TH2D(Form("FMD%d%c_vtxbin%03d", d, r, b), 
116                                 Form("Secondary map correction for FMD%d%c "
117                                      "in vertex bin %d [%+8.4f,%+8.4f]", 
118                                      d, r, b, vtxAxis.GetBinLowEdge(b), 
119                                      vtxAxis.GetBinUpEdge(b)), 
120                                 nEta, minEta, maxEta, 
121                                 oldmap->GetYaxis()->GetNbins(), 
122                                 oldmap->GetYaxis()->GetXmin(), 
123                                 oldmap->GetYaxis()->GetXmax());
124         newmap->SetXTitle("#eta");
125         newmap->SetYTitle("#phi [radians]");
126         newmap->SetZTitle("#sum_{i} N_{ch,i,primary} / #sum_{i} N_{ch,i,FMD}");
127         newmap->SetDirectory(0);
128         newmap->SetStats(0);
129         newmap->Sumw2();
130         newmap->Add(oldmap);
131
132         Info("RunCopySecMap",
133              "Copying %s to %s", oldmap->GetName(), newmap->GetName());
134
135         m->SetCorrection(d, r, b, newmap);
136       }
137     }
138   }
139   UShort_t isys = AliForwardUtil::ParseCollisionSystem(sys);
140   AliForwardCorrectionManager& mgr = AliForwardCorrectionManager::Instance();
141   TString fname(mgr.GetFileName(AliForwardCorrectionManager::kSecondaryMap,
142                                 isys, cms, field, false));
143   TFile* output = TFile::Open(fname.Data(), "RECREATE");
144   if (!output) { 
145     Warning("Run", "Failed to open output file %s", fname.Data());
146     return kFALSE;
147   }
148   m->Write(mgr.GetObjectName(AliForwardCorrectionManager::kSecondaryMap));
149   output->Write();
150   output->Close();
151   Info("Run", "File %s created.  It should be copied to %s and stored in SVN",
152        fname.Data(),mgr.GetFileDir(AliForwardCorrectionManager::kSecondaryMap));
153
154   fname = mgr.GetFileName(AliForwardCorrectionManager::kDoubleHit,
155                           isys, cms, field, false);
156   output = TFile::Open(fname.Data(), "RECREATE");
157   if (!output) { 
158     Warning("Run", "Failed to open output file %s", fname.Data());
159     return kFALSE;
160   }
161   dh->Write(mgr.GetObjectName(AliForwardCorrectionManager::kDoubleHit));
162   output->Write();
163   output->Close();
164   Info("Run", "File %s created.  It should be copied to %s and stored in SVN",
165        fname.Data(),mgr.GetFileDir(AliForwardCorrectionManager::kDoubleHit));
166 }
167 //
168 // EOF
169 //