]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/corrs/ExtractAcceptance.C
This commit has two major parts:
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / corrs / ExtractAcceptance.C
1 //_____________________________________________________________________
2 /** 
3  * 
4  * 
5  * @param d 
6  * @param r 
7  * @param vz 
8  * @param nDead 
9  * 
10  * @return 
11  *
12  * @ingroup pwglf_forward_scripts_corr
13  */
14 TH2D* MakeOneRing(UShort_t      d, 
15                   Char_t        r, 
16                   Double_t      vz, 
17                   Int_t&        nDead, 
18                   std::ostream* deadScript)
19 {
20   AliFMDGeometry*   geom = AliFMDGeometry::Instance();
21   AliFMDParameters* pars = AliFMDParameters::Instance();
22
23   UShort_t nS = (r == 'I' || r == 'i' ?  20 :  40);
24   UShort_t nT = (r == 'I' || r == 'i' ? 512 : 256);
25   
26   // Make our two histograms 
27   TH2D* hAll = new TH2D("all","All",200,-4,6,nS,0,2*TMath::Pi());
28   hAll->SetXTitle("#eta");
29   hAll->SetYTitle("#phi");
30   hAll->Sumw2();
31   hAll->SetDirectory(0);
32   TH2D* hOK  = static_cast<TH2D*>(hAll->Clone());
33   hOK->SetDirectory(0);
34   
35   if (deadScript)
36     *deadScript << "\n // FMD" << d << r << std::endl;
37   // Loop over all sectors and strips in this ring 
38   Int_t nOK  = 0;
39   Int_t nAll = 0;
40   Int_t nPhi = hAll->GetNbinsY();
41   for (UShort_t s = 0; s < nS; s++) { 
42     for (UShort_t t = 0; t < nT; t++) { 
43       // Get eta,phi by quering the geometry (first for (x,y,z), then
44       // correcting for the vertex position, and then calculating 
45       // (eta, phi))
46       Double_t x, y, z;
47       geom->Detector2XYZ(d, r, s, t, x, y, z);
48       z -= vz;
49       Double_t q, eta, phi, theta;
50       AliFMDGeometry::XYZ2REtaPhiTheta(x, y, z, q, eta, phi, theta);
51       if (phi < 0) phi += 2*TMath::Pi();
52
53       // Check if this is a dead channel or not 
54       Bool_t isDead = pars->IsDead(d, r, s, t);
55
56       // Special check for FMD2i - upper part of sectors 16/17 have 
57       // have anomalous gains and/or noise - common sources are 
58       // power regulartors for bias currents and the like 
59       Int_t VA = t/128;
60       if(d==2 && r=='I' && VA>1 && (s==16 || s==17)) isDead =true;
61
62       // Find the eta bin number and corresponding overflow bin
63       Int_t etaBin = hAll->GetXaxis()->FindBin(eta);
64       Int_t ovrBin = hAll->GetBin(etaBin, nPhi+1); 
65       
66       // Increment all histogram 
67       hAll->Fill(eta, phi);
68       hAll->AddBinContent(ovrBin);
69       nAll++;
70
71       // If not dead, increment OK histogram 
72       if (!isDead) {
73         hOK->Fill(eta, phi);
74         hOK->AddBinContent(ovrBin);
75         nOK++;
76       }
77       else {
78         nDead++;
79         if (deadScript)
80           *deadScript << "  filter->AddDead(" << d << ",'" << r << "'," 
81                       << s << ',' << t << ");" << std::endl;
82       }
83     }
84   }
85   // Divide out the efficiency. 
86   // Note, that the overflow bins along eta now contains the ratio 
87   // nOK/nAll Strips for a given eta bin. 
88   hOK->Divide(hOK,hAll,1,1,"B");
89
90   // Invert overflow bin
91   for (Int_t etaBin = 1; etaBin <= hOK->GetNbinsX(); etaBin++) { 
92     Double_t ovr    = hOK->GetBinContent(etaBin, nPhi+1);
93     Double_t novr   = (ovr < 1e-12 ? 0 : 1./ovr);
94     hOK->SetBinContent(etaBin, nPhi+1, novr);
95 #if 0
96     if (ovr > 0 && ovr != 1)
97       Info("", "Setting overflow bin (%3d,%3d) to 1/%f=%f", etaBin, nPhi+1, 
98            ovr, hOK->GetBinContent(etaBin, nPhi+1));
99 #endif
100   }
101   // Clean up
102   delete hAll;
103
104   Printf("=== FMD%d%c at vz=%+5.1f - %d/%d=%3d%% OK (w/overflow)", 
105          d, r, vz, nOK, nAll, (100*nOK)/nAll);
106
107   // Return result 
108   return hOK;
109 }
110
111 //_____________________________________________________________________
112 /** 
113  * 
114  * 
115  * @param runNo 
116  * @param system 
117  * @param energy 
118  * @param field 
119  * @param nVtxBins 
120  * @param vtxLow 
121  * @param vtxHigh 
122  *
123  * @ingroup pwglf_forward_scripts_corr
124  */
125 void ExtractAcceptance(Int_t   runNo=121526, 
126                        Int_t   nVtxBins=10, 
127                        Float_t vtxLow=-10, 
128                        Float_t vtxHigh=10)
129 {  
130   const char* fwd = "$ALICE_ROOT/../trunk/PWGLF/FORWARD/analysis2";
131   gSystem->AddIncludePath(Form("-I%s", fwd));
132   gROOT->Macro(Form("%s/scripts/LoadLibs.C", fwd));
133
134   // gSystem->Load("libANALYSIS");
135   // gSystem->Load("libANALYSISalice");
136   // gSystem->Load("libPWGLFforward2");
137   // Float_t delta = (vtxHigh - vtxLow) / (Float_t)nVtxBins;
138   
139   Bool_t gridOnline = kTRUE; 
140   if(!(TGrid::Connect("alien://",0,0,"t")))
141     gridOnline = kFALSE;
142
143   // --- Figure out the year --------------------------------------
144   UShort_t year = 0;
145   if      (runNo <= 99999)  year = 2009;
146   else if (runNo <= 139667) year = 2010;
147   else if (runNo <= 170718) year = 2011;
148   else if (runNo <= 194306) year = 2012;
149   else if (runNo <= 197709) year = 2013;
150   if (year <= 0) { 
151     Error("", "Couldn't deduce the year from the run number");
152     return;
153   }
154
155   // --- Initialisations ------------------------------------------
156   //Set up CDB manager
157   Printf("=== Setting up OCDB");
158   
159   AliCDBManager* cdb = AliCDBManager::Instance();
160   if(gridOnline) {
161     cdb->SetDefaultStorageFromRun(runNo);
162     // cdb->SetDefaultStorage(Form("alien://Folder=/alice/data/%4d/OCDB", year));
163   }
164   else
165     cdb->SetDefaultStorage("local://$(ALICE_ROOT)/OCDB");
166   cdb->SetRun(runNo);
167   
168   // Get the geometry 
169   Printf("=== Loading geometry");
170   AliGeomManager::LoadGeometry();
171
172   // Get an initialize parameters 
173   Printf("=== Intialising parameters");
174   AliFMDParameters* pars = AliFMDParameters::Instance();
175   pars->Init();
176
177   // Get an initialise geometry 
178   Printf("=== Initialising geomtry");
179   AliFMDGeometry* geom = AliFMDGeometry::Instance();
180   geom->Init();
181   geom->InitTransformations();
182
183   // --- Get the general run parameters ------------------------------
184   AliCDBEntry* grpE = cdb->Get("GRP/GRP/Data");
185   if (!grpE) { 
186     AliWarningF("No GRP entry found for run %d", runNo);
187     return;
188   }
189   AliGRPObject* grp = static_cast<AliGRPObject*>(grpE->GetObject());
190   if (!grp) { 
191     AliWarningF("No GRP object found for run %d", runNo);
192     return;
193   }
194   Float_t  beamE = grp->GetBeamEnergy();
195   TString  beamT = grp->GetBeamType();
196 # if 0 
197   // This isn't really needed as the acceptance map is indifferent to
198   // the field settings.
199   Float_t  l3cur = grp->GetL3Current(AliGRPObject::kMean);
200   Char_t   l3pol = grp->GetL3Polarity();
201   Bool_t   l3lhc = grp->IsPolarityConventionLHC();
202   Bool_t   l3uni = grp->IsUniformBMap();
203   AliMagF* fldM  = 
204     AliMagF::CreateFieldMap(TMath::Abs(l3cur) * (l3pol ? -1:1), 0, 
205                             (l3lhc ? 0 : 1), l3uni, beamE, beamT.Data());
206   Float_t  l3fld = fldM->SolenoidField();
207 #endif
208   
209   UShort_t sys = AliForwardUtil::ParseCollisionSystem(beamT);
210   UShort_t sNN = AliForwardUtil::ParseCenterOfMassEnergy(sys, 2 * beamE);
211   Short_t  fld = 0; // AliForwardUtil::ParseMagneticField(l3fld);
212   Printf("=== Run=%d, year=%d, sys=%d, sNN=%d, fld=%d", 
213          runNo, year, sys, sNN, fld);
214
215   // --- Output object -----------------------------------------------
216   // Make our correction object 
217   AliFMDCorrAcceptance* corr = new AliFMDCorrAcceptance();
218   corr->SetVertexAxis(nVtxBins, vtxLow, vtxHigh);
219
220   // --- Output script -----------------------------------------------
221   std::ofstream deadScript("deadstrips.C");
222   deadScript << "// Automatically generaeted by ExtractAcceptance.C\n"
223              << "// Add additional dead strips to sharing filter\n"
224              << "// Information taken from OCDB entry for\n" 
225              << "//\n"
226              << "//    run       = " << runNo << "\n"
227              << "//    year      = " << year << "\n"
228              << "//    system    = " << sys << "\n"
229              << "//    sqrt{sNN} = " << sNN << "GeV\n"
230              << "//    L3 field  = " << fld << "kG\n"
231              << "void deadstrips(AliFMDSharingFilter* filter)\n"
232              << "{" << std::endl;
233   
234   // --- Loop over verticies and rings -------------------------------
235   Int_t  nDead   = 0;
236   Bool_t gotDead = false;
237   Float_t dV     = (vtxHigh - vtxLow) / nVtxBins;
238   Printf("=== Looping over vertices: %d bins from %+6.2f to %+6.2f "
239          "in steps of %+6.2f", 
240          nVtxBins, vtxLow, vtxHigh, dV);
241   for (Double_t v = vtxLow+dV/2; v < vtxHigh; v += dV) { 
242     for(UShort_t d = 1; d <= 3;d++) { 
243       UShort_t nR = (d == 1 ? 1 : 2);
244       for (UShort_t q = 0; q < nR; q++) { 
245         Char_t   r  = (q == 0 ? 'I' : 'O');
246
247         // Delegate to other function 
248         Int_t nLocal = 0;
249         TH2D* ratio = MakeOneRing(d, r, v, nLocal, 
250                                   !gotDead ? &deadScript : 0);
251         nDead += nLocal;
252         if (!ratio) {
253           Warning("ExtractAcceptance", "Didn't get correction from "
254                   "FMD%d%c @ vz=%+6.2fcm", d, r, v);
255           continue;
256         }
257         // Printf("v=%+6.2f FMD%d%c, got %d dead strips", v, d, r, nLocal);
258
259         // Set the correction 
260         corr->SetCorrection(d, r, v, ratio);
261       }
262     }
263     gotDead = true;
264   }
265   corr->SetHasOverflow();
266   corr->Print();
267   // corr->ls();
268
269   deadScript << "}\n"
270              << "// EOF" << std::endl;
271   deadScript.close();
272
273   // Write to a file 
274   Printf("=== Writing to disk");
275   AliForwardCorrectionManager& cm = AliForwardCorrectionManager::Instance();
276   if (!cm.Store(corr, runNo, sys, sNN, fld, false, false, 
277                 "fmd_corrections.root")) { 
278     Error("", "Failed to store acceptance correction in local file");
279     return;
280   }
281
282   std::ofstream f("Upload.C");
283   f << "// Generated by ExtractELoss.C\n"
284     << "TString MakeDest(const TString& dest, const TString& fname)\n"
285     << "{\n"
286     << "  TString tmp(dest);\n"
287     << "  if (!tmp.IsNull()) {\n"
288     << "    if (!tmp.EndsWith(\"/\")) tmp.Append(\"/\");\n"
289     << "    tmp.Append(fname);\n"
290     << "  }\n"
291     << "  return tmp;\n"
292     << "}\n\n"
293     << "void Upload(const TString& dest=\"\")\n"
294     << "{\n"
295     << "  gROOT->Macro(\"" << fwd << "/scripts/LoadLibs.C\");\n"
296     << "  \n"
297     << "  const char* fmdFile = \"fmd_corrections.root\";\n"
298     << "  TString fdest = MakeDest(dest, fmdFile);\n"
299     << "  \n"
300     << "  AliForwardCorrectionManager::Instance().Append(fmdFile, fdest);\n"
301     << "}\n"
302     << "// EOF\n"
303     << std::endl;
304   f.close();
305   
306   Info("ExtracAcceptance", 
307        "Run generated Upload.C(DEST) script to copy files in place");
308 }
309
310 //
311 // EOF
312 //
313
314