]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FORWARD/analysis2/corrs/MakeAcceptanceCorrection.C
51c68edaca6a1619b87bfd0b19cb9e86ad34d641
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis2 / corrs / MakeAcceptanceCorrection.C
1 //_____________________________________________________________________
2 /** 
3  * 
4  * 
5  * @param d 
6  * @param r 
7  * @param vz 
8  * @param nDead 
9  * 
10  * @return 
11  *
12  * @ingroup pwg2_forward_analysis_scripts_corr
13  */
14 TH2D* MakeOneRing(UShort_t d, Char_t r, Double_t vz, Int_t& nDead)
15 {
16   AliFMDGeometry*   geom = AliFMDGeometry::Instance();
17   AliFMDParameters* pars = AliFMDParameters::Instance();
18
19   UShort_t nS = (r == 'I' || r == 'i' ?  20 :  40);
20   UShort_t nT = (r == 'I' || r == 'i' ? 512 : 256);
21   
22   // Make our two histograms 
23   TH2D* hAll = new TH2D("all","All",200,-4,6,nS,0,2*TMath::Pi());
24   hAll->SetXTitle("#eta");
25   hAll->SetYTitle("#phi");
26   hAll->Sumw2();
27   hAll->SetDirectory(0);
28   TH2D* hOK  = static_cast<TH2D*>(hAll->Clone());
29   hOK->SetDirectory(0);
30   
31   // Loop over all sectors and strips in this ring 
32   Int_t nOK  = 0;
33   Int_t nAll = 0;
34   for (UShort_t s = 0; s < nS; s++) { 
35     for (UShort_t t = 0; t < nT; t++) { 
36       // Get eta,phi by quering the geometry (first for (x,y,z), then
37       // correcting for the vertex position, and then calculating 
38       // (eta, phi))
39       Double_t x, y, z;
40       geom->Detector2XYZ(d, r, s, t, x, y, z);
41       z -= vz;
42       Double_t q, eta, phi, theta;
43       AliFMDGeometry::XYZ2REtaPhiTheta(x, y, z, q, eta, phi, theta);
44       if (phi < 0) phi += 2*TMath::Pi();
45
46       // Check if this is a dead channel or not 
47       Bool_t isDead = pars->IsDead(d, r, s, t);
48       hAll->Fill(eta, phi);
49       nAll++;
50       if (!isDead) {
51         hOK->Fill(eta, phi);
52         nOK++;
53       }
54       else         nDead++;
55     }
56   }
57   // Divide out the efficiency. 
58   hOK->Divide(hOK,hAll,1,1,"B");
59
60   // Clean up
61   delete hAll;
62
63   Info("MakeAcceptanceCorrections","Made correction for FMD%d%c at vz=%f - "
64        "%d strips out of %d OK", d, r, vz, nOK, nAll);
65
66   // Return result 
67   return hOK;
68 }
69
70 //_____________________________________________________________________
71 /** 
72  * 
73  * 
74  * @param runNo 
75  * @param system 
76  * @param energy 
77  * @param field 
78  * @param nVtxBins 
79  * @param vtxLow 
80  * @param vtxHigh 
81  *
82  * @ingroup pwg2_forward_analysis_scripts_corr
83  */
84 void MakeAcceptanceCorrection(Int_t   runNo=121526, 
85                               Int_t   system = 1,
86                               Float_t energy = 900,
87                               Float_t field  = 5,
88                               Int_t   nVtxBins=10, 
89                               Float_t vtxLow=-10, 
90                               Float_t vtxHigh=10){
91   
92   gSystem->Load("libANALYSIS");
93   gSystem->Load("libANALYSISalice");
94   gSystem->Load("libPWG2forward2");
95   
96   // Float_t delta = (vtxHigh - vtxLow) / (Float_t)nVtxBins;
97   
98   Bool_t kGridOnline = kTRUE; 
99   if(!(TGrid::Connect("alien://",0,0,"t")))
100     kGridOnline = kFALSE;
101   
102   // --- Initialisations ------------------------------------------
103   //Set up CDB manager
104   Info("MakeAcceptanceCorrections","Setting up OCDB");
105   
106   AliCDBManager* cdb = AliCDBManager::Instance();
107   if(kGridOnline)
108     cdb->SetDefaultStorage("alien://Folder=/alice/data/2010/OCDB");
109   else
110     cdb->SetDefaultStorage("local://$(ALICE_ROOT)/OCDB");
111   cdb->SetRun(runNo);
112   
113   // Get the geometry 
114   Info("MakeAcceptanceCorrections","Loading geometry");
115   AliGeomManager::LoadGeometry();
116
117   // Get an initialize parameters 
118   Info("MakeAcceptanceCorrections","Intialising parameters");
119   AliFMDParameters* pars = AliFMDParameters::Instance();
120   pars->Init();
121
122   // Get an initialise geometry 
123   Info("MakeAcceptanceCorrections","Initialising geomtry");
124   AliFMDGeometry* geom = AliFMDGeometry::Instance();
125   geom->Init();
126   geom->InitTransformations();
127
128   // --- Output object -----------------------------------------------
129   // Make our correction object 
130   AliFMDCorrAcceptance* corr = new AliFMDCorrAcceptance();
131   corr->SetVertexAxis(nVtxBins, vtxLow, vtxHigh);
132
133   // --- Loop over verticies and rings -------------------------------
134   Int_t nDead = 0;
135   Float_t dV = (vtxHigh - vtxLow) / nVtxBins;
136   for (Double_t v = vtxLow+dV/2; v < vtxHigh; v += dV) { 
137     for(UShort_t d = 1; d <= 3;d++) { 
138       UShort_t nR = (d == 1 ? 1 : 2);
139       for (UShort_t q = 0; q < nR; q++) { 
140         Char_t   r  = (q == 0 ? 'I' : 'O');
141
142         // Delegate to other function 
143         TH2D* ratio = MakeOneRing(d, r, v, nDead);
144         if (!ratio) continue;
145
146         // Set the correction 
147         corr->SetCorrection(d, r, v, ratio);
148       }
149     }
150   }
151
152   // Write to a file 
153 #if 0
154   TFile* out = TFile::Open(Form("acceptance_%d.root", runNo), "RECREATE");
155   corr->Write("acceptance");
156   out->Write();
157   out->Close();
158 #else 
159   Info("MakeAcceptanceCorrections","Writing to disk");
160   AliForwardCorrectionManager& cm = AliForwardCorrectionManager::Instance();
161   TString fname = cm.GetFileName(AliForwardCorrectionManager::kAcceptance, 
162                                  system, energy, field, false);
163   TFile* out = TFile::Open(fname.Data(), "RECREATE");
164   corr->Write(cm.GetObjectName(AliForwardCorrectionManager::kAcceptance));
165   out->Write();
166   out->Close();
167   Info("MakeAcceptanceCorrections","File %s generated.  "
168        "It should be copied to %s", fname.Data(), 
169        cm.GetFileDir(AliForwardCorrectionManager::kAcceptance));
170   out = TFile::Open(fname.Data(), "READ");
171   new TBrowser;
172 #endif
173   
174 }
175
176 //
177 // EOF
178 //
179