]>
Commit | Line | Data |
---|---|---|
7e24fa8c | 1 | /**************************************************************************** |
2 | * Very important, delicate and rather obscure macro. * | |
3 | * ("a la" AliTPCComparison.C) * | |
4 | * * | |
5 | * Creates list of "trackable" tracks, * | |
6 | * calculates efficiency, resolutions etc. * | |
7 | * (To get the list of the "trackable" tracks one should * | |
8 | * first run the AliTPCComparison.C macro) * | |
9 | * * | |
10 | * Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch * | |
11 | ****************************************************************************/ | |
12 | ||
13 | #if !defined(__CINT__) || defined(__MAKECINT__) | |
14 | #include <TMath.h> | |
15 | #include <TError.h> | |
16 | #include <Riostream.h> | |
17 | #include <TH1F.h> | |
18 | #include <TH2F.h> | |
19 | #include <TTree.h> | |
20 | #include <TParticle.h> | |
21 | #include <TCanvas.h> | |
22 | #include <TLine.h> | |
23 | #include <TText.h> | |
24 | #include <TBenchmark.h> | |
25 | #include <TStyle.h> | |
26 | #include <TFile.h> | |
27 | #include <TROOT.h> | |
28 | ||
29 | #include "AliHeader.h" | |
30 | #include "AliTrackReference.h" | |
31 | #include "AliRunLoader.h" | |
32 | #include "AliRun.h" | |
33 | #include "AliESD.h" | |
34 | #endif | |
35 | ||
36 | Int_t GoodTracksTRD(const Char_t *dir="."); | |
37 | ||
38 | extern AliRun *gAlice; | |
39 | extern TBenchmark *gBenchmark; | |
40 | extern TROOT *gROOT; | |
41 | ||
42 | static Int_t allgood=0; | |
43 | static Int_t allselected=0; | |
44 | static Int_t allfound=0; | |
45 | ||
46 | Int_t AliTRDComparisonV2 | |
47 | (Float_t ptcutl=0.2, Float_t ptcuth=10., const Char_t *dir=".") { | |
48 | gBenchmark->Start("AliTRDComparisonV2"); | |
49 | ||
50 | ::Info("AliTRDComparisonV2.C","Doing comparison..."); | |
51 | ||
52 | ||
53 | TH1F *hp=(TH1F*)gROOT->FindObject("hp"); | |
54 | if (!hp) { | |
55 | hp=new TH1F("hp","PHI resolution",50,-70.,70.); | |
56 | hp->SetFillColor(4); | |
57 | hp->SetXTitle("(mrad)"); | |
58 | } | |
59 | TH1F *hl=(TH1F*)gROOT->FindObject("hl"); | |
60 | if (!hl) { | |
61 | hl=new TH1F("hl","LAMBDA resolution",50,-70,70); | |
62 | hl->SetFillColor(4); | |
63 | hl->SetXTitle("(mrad)"); | |
64 | } | |
65 | TH1F *hc=(TH1F*)gROOT->FindObject("hc"); | |
66 | if (!hc) { | |
67 | hc=new TH1F("hc","Number of the assigned clusters",25,110,135); | |
68 | hc->SetLineColor(2); | |
69 | } | |
70 | TH1F *hpt=(TH1F*)gROOT->FindObject("hpt"); | |
71 | if (!hpt) { | |
72 | hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.); | |
73 | hpt->SetFillColor(2); | |
74 | hpt->SetXTitle("(%)"); | |
75 | } | |
76 | TH1F *hmpt=(TH1F*)gROOT->FindObject("hmpt"); | |
77 | if (!hmpt) { | |
78 | hmpt=new TH1F("hmpt","Y and Z resolution",30,-30,30); | |
79 | hmpt->SetFillColor(6); | |
80 | hmpt->SetXTitle("(mm)"); | |
81 | } | |
82 | TH1F *hz=(TH1F*)gROOT->FindObject("hz"); | |
83 | if (!hz) hz=new TH1F("hz","Z resolution",30,-30,30); | |
84 | ||
85 | ||
86 | ||
87 | TH1F *hgood=(TH1F*)gROOT->FindObject("hgood"); | |
88 | if (!hgood) hgood=new TH1F("hgood","Good tracks",30,0.2,6.1); | |
89 | ||
90 | TH1F *hfound=(TH1F*)gROOT->FindObject("hfound"); | |
91 | if (!hfound) hfound=new TH1F("hfound","Found tracks",30,0.2,6.1); | |
92 | ||
93 | TH1F *hfake=(TH1F*)gROOT->FindObject("hfake"); | |
94 | if (!hfake) hfake=new TH1F("hfake","Fake tracks",30,0.2,6.1); | |
95 | ||
96 | TH1F *hg=(TH1F*)gROOT->FindObject("hg"); | |
97 | if (!hg) hg=new TH1F("hg","Efficiency for good tracks",30,0.2,6.1); | |
98 | hg->SetLineColor(4); hg->SetLineWidth(2); | |
99 | ||
100 | TH1F *hf=(TH1F*)gROOT->FindObject("hf"); | |
101 | if (!hf) hf=new TH1F("hf","Efficiency for fake tracks",30,0.2,6.1); | |
102 | hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2); | |
103 | ||
104 | TH1F *he=(TH1F*)gROOT->FindObject("he"); | |
105 | if (!he) | |
106 | he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,1000.); | |
107 | ||
108 | TH2F *hep=(TH2F*)gROOT->FindObject("hep"); | |
109 | if (!hep) hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,2000.); | |
110 | hep->SetMarkerStyle(8); | |
111 | hep->SetMarkerSize(0.4); | |
112 | ||
113 | ||
114 | Char_t fname[100]; | |
115 | sprintf(fname,"%s/GoodTracksTRD.root",dir); | |
116 | ||
117 | TFile *refFile=TFile::Open(fname,"old"); | |
118 | if (!refFile || !refFile->IsOpen()) { | |
119 | ::Info("AliTRDComparisonV2.C","Marking good tracks (will take a while)..."); | |
120 | if (GoodTracksTRD(dir)) { | |
121 | ::Error("AliTRDComparisonV2.C","Can't generate the reference file !"); | |
122 | return 1; | |
123 | } | |
124 | } | |
125 | refFile=TFile::Open(fname,"old"); | |
126 | if (!refFile || !refFile->IsOpen()) { | |
127 | ::Error("AliTRDComparisonV2.C","Can't open the reference file !"); | |
128 | return 1; | |
129 | } | |
130 | ||
131 | TTree *trdTree=(TTree*)refFile->Get("trdTree"); | |
132 | if (!trdTree) { | |
133 | ::Error("AliTRDComparisonV2.C","Can't get the reference tree !"); | |
134 | return 2; | |
135 | } | |
136 | TBranch *branch=trdTree->GetBranch("TRD"); | |
137 | if (!branch) { | |
138 | ::Error("AliTRDComparisonV2.C","Can't get the TRD branch !"); | |
139 | return 3; | |
140 | } | |
141 | TClonesArray dummy("AliTrackReference",1000), *refs=&dummy; | |
142 | branch->SetAddress(&refs); | |
143 | ||
144 | ||
145 | sprintf(fname,"%s/AliESDs.root",dir); | |
146 | TFile *ef=TFile::Open(fname); | |
147 | if ((!ef)||(!ef->IsOpen())) { | |
148 | sprintf(fname,"%s/AliESDtrd.root",dir); | |
149 | ef=TFile::Open(fname); | |
150 | if ((!ef)||(!ef->IsOpen())) { | |
151 | ::Error("AliTRDComparisonV2.C","Can't open AliESDtrd.root !"); | |
152 | return 4; | |
153 | } | |
154 | } | |
155 | AliESD* event = new AliESD; | |
156 | TTree* esdTree = (TTree*) ef->Get("esdTree"); | |
157 | if (!esdTree) { | |
158 | ::Error("AliTRDComparisonV2.C", "no ESD tree found"); | |
159 | return 6; | |
160 | } | |
161 | esdTree->SetBranchAddress("ESD", &event); | |
162 | ||
163 | ||
164 | //******* Loop over events ********* | |
165 | Int_t e=0; | |
166 | while (esdTree->GetEvent(e)) { | |
167 | cout<<endl<<endl<<"********* Processing event number: "<<e<<"*******\n"; | |
168 | ||
169 | Int_t nentr=event->GetNumberOfTracks(); | |
170 | allfound+=nentr; | |
171 | ||
172 | if (trdTree->GetEvent(e++)==0) { | |
173 | cerr<<"No reconstructable tracks !\n"; | |
174 | continue; | |
175 | } | |
176 | ||
177 | Int_t ngood=refs->GetEntriesFast(); | |
178 | allgood+=ngood; | |
179 | ||
180 | const Int_t MAX=15000; | |
181 | Int_t notf[MAX], nnotf=0; | |
182 | Int_t fake[MAX], nfake=0; | |
183 | Int_t mult[MAX], numb[MAX], nmult=0; | |
184 | Int_t k; | |
185 | for (k=0; k<ngood; k++) { | |
186 | AliTrackReference *ref=(AliTrackReference*)refs->UncheckedAt(k); | |
187 | Int_t lab=ref->Label(), tlab=-1; | |
188 | Float_t ptg=ref->Pt(); | |
189 | ||
190 | if (ptg<ptcutl) continue; | |
191 | if (ptg>ptcuth) continue; | |
192 | ||
193 | allselected++; | |
194 | ||
195 | hgood->Fill(ptg); | |
196 | ||
197 | AliESDtrack *esd=0; | |
198 | Int_t cnt=0; | |
199 | for (Int_t i=0; i<nentr; i++) { | |
200 | AliESDtrack *t=event->GetTrack(i); | |
201 | UInt_t status=t->GetStatus(); | |
202 | ||
203 | if ((status&AliESDtrack::kTRDout)==0) continue; | |
204 | ||
205 | Int_t lbl=t->GetTRDLabel(); | |
206 | if (lab==TMath::Abs(lbl)) { | |
207 | if (cnt==0) {esd=t; tlab=lbl;} | |
208 | cnt++; | |
209 | } | |
210 | } | |
211 | if (cnt==0) { | |
212 | notf[nnotf++]=lab; | |
213 | continue; | |
214 | } else if (cnt>1){ | |
215 | mult[nmult]=lab; | |
216 | numb[nmult]=cnt; nmult++; | |
217 | } | |
218 | ||
219 | if (lab==tlab) hfound->Fill(ptg); | |
220 | else { | |
221 | fake[nfake++]=lab; | |
222 | hfake->Fill(ptg); | |
223 | } | |
224 | ||
225 | AliExternalTrackParam out(*esd->GetOuterParam()); | |
226 | ||
227 | Double_t bz=event->GetMagneticField(); | |
228 | Double_t xg = ref->LocalX(); | |
229 | out.PropagateTo(xg,bz); | |
230 | ||
231 | Float_t phi=TMath::ASin(out.GetSnp()) + out.GetAlpha(); | |
232 | if (phi < 0) phi+=2*TMath::Pi(); | |
233 | if (phi >=2*TMath::Pi()) phi-=2*TMath::Pi(); | |
234 | Float_t phig=ref->Phi(); | |
235 | hp->Fill((phi - phig)*1000.); | |
236 | ||
237 | Float_t lam=TMath::ATan(out.GetTgl()); | |
238 | Float_t lamg=TMath::ATan2(ref->Pz(),ptg); | |
239 | hl->Fill((lam - lamg)*1000.); | |
240 | ||
241 | hc->Fill(esd->GetTRDclusters(0)); | |
242 | ||
243 | Float_t pt_1=TMath::Abs(out.Get1Pt()); | |
244 | hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.); | |
245 | ||
246 | Float_t y=out.GetY(); | |
247 | Float_t yg=ref->LocalY(); | |
248 | hmpt->Fill(10*(y - yg)); | |
249 | ||
250 | Float_t z=out.GetZ(); | |
251 | Float_t zg=ref->Z(); | |
252 | hz->Fill(10*(z - zg)); | |
253 | ||
254 | Float_t mom=1./(pt_1*TMath::Cos(lam)); | |
255 | Float_t dedx=esd->GetTRDsignal(); | |
256 | hep->Fill(mom,dedx,1.); | |
257 | ||
258 | Int_t pdg=(Int_t)ref->GetLength(); //this is particle's PDG ! | |
259 | ||
260 | if (TMath::Abs(pdg)==211) //pions | |
261 | if (mom>0.4 && mom<0.5) he->Fill(dedx,1.); | |
262 | ||
263 | } | |
264 | ||
265 | cout<<"\nList of Not found tracks :\n"; | |
266 | for (k=0; k<nnotf; k++){ | |
267 | cout<<notf[k]<<"\t"; | |
268 | if ((k%9)==8) cout<<"\n"; | |
269 | } | |
270 | cout<<"\n\nList of fake tracks :\n"; | |
271 | for (k=0; k<nfake; k++){ | |
272 | cout<<fake[k]<<"\t"; | |
273 | if ((k%9)==8) cout<<"\n"; | |
274 | } | |
275 | cout<<"\n\nList of multiple found tracks :\n"; | |
276 | for (k=0; k<nmult; k++) { | |
277 | cout<<"id. "<<mult[k] | |
278 | <<" found - "<<numb[k]<<"times\n"; | |
279 | } | |
280 | cout<<endl; | |
281 | ||
282 | cout<<"Number of found tracks : "<<nentr<<endl; | |
283 | cout<<"Number of \"good\" tracks : "<<ngood<<endl; | |
284 | ||
285 | refs->Clear(); | |
286 | } //***** End of the loop over events | |
287 | ||
288 | delete event; | |
289 | ef->Close(); | |
290 | ||
291 | delete trdTree; | |
292 | refFile->Close(); | |
293 | ||
294 | Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries(); | |
295 | if (ng!=0) cout<<"\n\nIntegral efficiency is about "<<nf/ng*100.<<" %\n"; | |
296 | cout<<"Total number selected of \"good\" tracks ="<<allselected<<endl<<endl; | |
297 | cout<<"Total number of found tracks ="<<allfound<<endl; | |
298 | cout<<"Total number of \"good\" tracks ="<<allgood<<endl; | |
299 | cout<<endl; | |
300 | ||
301 | gStyle->SetOptStat(111110); | |
302 | gStyle->SetOptFit(1); | |
303 | ||
304 | TCanvas *c1=new TCanvas("c1","",0,0,700,850); | |
305 | ||
306 | Int_t minc=33; | |
307 | ||
308 | TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw(); | |
309 | p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10); | |
310 | if (hp->GetEntries()<minc) hp->Draw(); else hp->Fit("gaus"); c1->cd(); | |
311 | ||
312 | TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw(); | |
313 | p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10); | |
314 | if (hl->GetEntries()<minc) hl->Draw(); else hl->Fit("gaus"); c1->cd(); | |
315 | ||
316 | TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw(); | |
317 | p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10); | |
318 | if (hpt->GetEntries()<minc) hpt->Draw(); else hpt->Fit("gaus"); c1->cd(); | |
319 | ||
320 | TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw(); | |
321 | p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10); | |
322 | if (hmpt->GetEntries()<minc) hmpt->Draw(); else hmpt->Fit("gaus"); | |
323 | hz->Draw("same"); c1->cd(); | |
324 | ||
325 | ||
326 | TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd(); | |
327 | p5->SetFillColor(41); p5->SetFrameFillColor(10); | |
328 | hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2(); | |
329 | hg->Divide(hfound,hgood,1,1.,"b"); | |
330 | hf->Divide(hfake,hgood,1,1.,"b"); | |
331 | hg->SetMaximum(1.4); | |
332 | hg->SetYTitle("Tracking efficiency"); | |
333 | hg->SetXTitle("Pt (GeV/c)"); | |
334 | hg->Draw(); | |
335 | ||
336 | TLine *line1 = new TLine(0.2,1.0,6.1,1.0); line1->SetLineStyle(4); | |
337 | line1->Draw("same"); | |
338 | TLine *line2 = new TLine(0.2,0.9,6.1,0.9); line2->SetLineStyle(4); | |
339 | line2->Draw("same"); | |
340 | ||
341 | hf->SetFillColor(1); | |
342 | hf->SetFillStyle(3013); | |
343 | hf->SetLineColor(2); | |
344 | hf->SetLineWidth(2); | |
345 | hf->Draw("histsame"); | |
346 | TText *text = new TText(0.461176,0.248448,"Fake tracks"); | |
347 | text->SetTextSize(0.05); | |
348 | text->Draw(); | |
349 | text = new TText(0.453919,1.11408,"Good tracks"); | |
350 | text->SetTextSize(0.05); | |
351 | text->Draw(); | |
352 | ||
353 | TCanvas *c2=new TCanvas("c2","",320,32,530,590); | |
354 | TPad *p6=new TPad("p6","",0.,0.,1.,.5); p6->Draw(); | |
355 | p6->cd(); p6->SetFillColor(42); p6->SetFrameFillColor(10); | |
356 | he->SetFillColor(2); he->SetFillStyle(3005); | |
357 | he->SetXTitle("Arbitrary Units"); | |
358 | if (he->GetEntries()<minc) he->Draw(); else he->Fit("gaus"); c2->cd(); | |
359 | ||
360 | TPad *p7=new TPad("p7","",0.,0.5,1.,1.); p7->Draw(); | |
361 | p7->cd(); p7->SetFillColor(42); p7->SetFrameFillColor(10); | |
362 | hep->SetXTitle("p (Gev/c)"); hep->SetYTitle("dE/dX (Arb. Units)"); | |
363 | hep->Draw(); c1->cd(); | |
364 | ||
365 | TFile fc("AliTRDComparisonV2.root","RECREATE"); | |
366 | c1->Write(); | |
367 | c2->Write(); | |
368 | fc.Close(); | |
369 | ||
370 | gBenchmark->Stop("AliTRDComparisonV2"); | |
371 | gBenchmark->Show("AliTRDComparisonV2"); | |
372 | ||
373 | return 0; | |
374 | } | |
375 | ||
376 | ||
377 | ||
378 | Int_t GoodTracksTRD(const Char_t *dir) { | |
379 | if (gAlice) { | |
380 | delete gAlice->GetRunLoader(); | |
381 | delete gAlice;//if everything was OK here it is already NULL | |
382 | gAlice = 0x0; | |
383 | } | |
384 | ||
385 | Char_t fname[100]; | |
386 | sprintf(fname,"%s/galice.root",dir); | |
387 | ||
388 | AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON"); | |
389 | if (!rl) { | |
390 | ::Error("GoodTracksTRD","Can't start session !"); | |
391 | return 1; | |
392 | } | |
393 | ||
394 | rl->LoadgAlice(); | |
395 | rl->LoadHeader(); | |
396 | rl->LoadTrackRefs(); | |
397 | ||
398 | Int_t nev=rl->GetNumberOfEvents(); | |
399 | ::Info("GoodTracksTRD","Number of events : %d\n",nev); | |
400 | ||
401 | sprintf(fname,"%s/GoodTracksTPC.root",dir); | |
402 | TFile *tpcFile=TFile::Open(fname); | |
403 | if ((!tpcFile)||(!tpcFile->IsOpen())) { | |
404 | ::Error("GoodTracksTRD","Can't open the GoodTracksTPC.root !"); | |
405 | delete rl; | |
406 | return 5; | |
407 | } | |
408 | TClonesArray dum("AliTrackReference",1000), *tpcRefs=&dum; | |
409 | TTree *tpcTree=(TTree*)tpcFile->Get("tpcTree"); | |
410 | if (!tpcTree) { | |
411 | ::Error("GoodTracksTRD","Can't get the TPC reference tree !"); | |
412 | delete rl; | |
413 | return 6; | |
414 | } | |
415 | TBranch *tpcBranch=tpcTree->GetBranch("TPC"); | |
416 | if (!tpcBranch) { | |
417 | ::Error("GoodTracksTRD","Can't get the TPC reference branch !"); | |
418 | delete rl; | |
419 | return 7; | |
420 | } | |
421 | tpcBranch->SetAddress(&tpcRefs); | |
422 | ||
423 | sprintf(fname,"%s/GoodTracksTRD.root",dir); | |
424 | TFile *trdFile=TFile::Open(fname,"recreate"); | |
425 | TClonesArray dummy2("AliTrackReference",1000), *trdRefs=&dummy2; | |
426 | TTree trdTree("trdTree","Info about the reconstructable TRD tracks"); | |
427 | trdTree.Branch("TRD",&trdRefs); | |
428 | ||
429 | //******** Loop over generated events | |
430 | for (Int_t e=0; e<nev; e++) { | |
431 | Int_t k; | |
432 | ||
433 | rl->GetEvent(e); trdFile->cd(); | |
434 | ||
435 | Int_t np = rl->GetHeader()->GetNtrack(); | |
436 | cout<<"Event "<<e<<" Number of particles: "<<np<<endl; | |
437 | ||
438 | ||
439 | TTree *TR=rl->TreeTR(); | |
440 | TBranch *branch=TR->GetBranch("TRD"); | |
441 | if (branch==0) { | |
442 | ::Error("GoodTracksTRD","No TRD track references !"); | |
443 | delete rl; | |
444 | return 5; | |
445 | } | |
446 | TClonesArray dummy("AliTrackReference",1000), *refs=&dummy; | |
447 | branch->SetAddress(&refs); | |
448 | Int_t nr=TR->GetEntries(); | |
449 | ||
450 | //Preselect the "good" track candidates (TRD info only) | |
451 | Int_t nt=0; | |
452 | for (Int_t r=0; r<nr; r++) { | |
453 | refs->Clear(); | |
454 | TR->GetEntry(r); | |
455 | Int_t n=refs->GetEntriesFast(); | |
456 | ||
457 | if (n<=1) continue; | |
458 | ||
459 | AliTrackReference *ref0=(AliTrackReference *)refs->UncheckedAt(0); | |
460 | if (ref0->LocalX() > 300.) continue; | |
461 | ||
462 | AliTrackReference *refn=(AliTrackReference *)refs->UncheckedAt(n-1); | |
463 | if (refn->LocalX() < 363.) continue; | |
464 | ||
465 | if (TMath::Abs(ref0->Alpha() - refn->Alpha()) > 1e-5) continue; | |
466 | ||
467 | new((*trdRefs)[nt++]) AliTrackReference(*refn); | |
468 | } | |
469 | ||
470 | //Check if the candidates are "good" from the TPC point of view | |
471 | tpcTree->GetEvent(e); | |
472 | Int_t ntpc=tpcRefs->GetEntriesFast(); | |
473 | for (Int_t t=0; t<nt; t++) { | |
474 | AliTrackReference *ref=(AliTrackReference *)trdRefs->UncheckedAt(t); | |
475 | Int_t lab=ref->Label(); | |
476 | for (k=0; k<ntpc; k++) { | |
477 | AliTrackReference | |
478 | *tpcRef=(AliTrackReference *)tpcRefs->UncheckedAt(k); | |
479 | if (tpcRef->Label()==lab) { | |
480 | ref->SetLength(tpcRef->GetLength()); | |
481 | break; | |
482 | } | |
483 | } | |
484 | if (k==ntpc) delete trdRefs->RemoveAt(t); | |
485 | } | |
486 | trdRefs->Compress(); | |
487 | trdTree.Fill(); | |
488 | ||
489 | trdRefs->Clear(); | |
490 | tpcRefs->Clear(); | |
491 | ||
492 | } //*** end of the loop over generated events | |
493 | ||
494 | trdTree.Write(); | |
495 | trdFile->Close(); | |
496 | ||
497 | delete tpcTree; | |
498 | tpcFile->Close(); | |
499 | ||
500 | delete rl; | |
501 | return 0; | |
502 | } | |
503 | ||
504 |