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