]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/testITSU/compClusHits.C
Removing printout (Davide)
[u/mrichter/AliRoot.git] / ITS / UPGRADE / testITSU / compClusHits.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include "TObjArray.h"
3 #include "TFile.h"
4 #include "TTree.h"
5 #include "TH1F.h"
6 #include "../ITS/UPGRADE/AliITSUClusterPix.h"
7 #include "../ITS/UPGRADE/AliITSURecoLayer.h"
8 #include "../ITS/UPGRADE/AliITSURecoDet.h"
9 #include "../ITS/UPGRADE/AliITSUHit.h"
10 #include "../ITS/UPGRADE/AliITSUGeomTGeo.h"
11 #include "AliITSsegmentation.h"
12 #include "AliGeomManager.h"
13 #include "AliStack.h"
14 #include "AliLoader.h"
15 #include "AliCDBManager.h"
16
17 #include "TROOT.h"
18 #include "TStyle.h"
19 #include "TGeoMatrix.h"
20 #include "TParticle.h"
21 #include "TCanvas.h"
22 #include "TPaveStats.h"
23 #include "TClonesArray.h"
24
25 #endif
26
27 TObjArray histoArr;
28 enum {kNPixAll=0,kNPixSPL=1,kDR=0,kDTXodd,kDTXeven,kDTZ, kDTXoddSPL,kDTXevenSPL,kDTZSPL};
29
30 TPaveStats* GetStPad(TH1* hst);
31 TPaveStats* SetStPadPos(TH1* hst,float x1,float x2,float y1,float y2, Int_t stl=-1,Int_t col=-1);
32 TCanvas* DrawNP(int np, TObjArray* harr=0, TCanvas* cnv=0);
33 TH1* GetHistoClSize(int npix,int id,TObjArray* harr=0);
34 void DrawReport(const char* psname, TObjArray* harr=0);
35
36
37 typedef struct {
38   Int_t evID;
39   Int_t volID;
40   Int_t lrID;
41   Int_t clID;
42   Int_t nPix;
43   Int_t nX;
44   Int_t nZ;
45   Int_t q;
46   Float_t pt;
47   Float_t eta;
48   Float_t phi;
49   Float_t xyz[3];
50   Float_t dX;
51   Float_t dZ;  
52   Bool_t split;  
53   Bool_t prim;
54   Int_t  pdg;
55   Int_t  ntr;
56 } clSumm;
57
58 void compClusHits(int nev=-1)
59 {
60   const int kSplit=0x1<<22;
61   const int kSplCheck=0x1<<23;
62   //
63   gSystem->Load("libITSUpgradeBase");
64   gSystem->Load("libITSUpgradeSim");
65   gSystem->Load("libITSUpgradeRec");
66   gROOT->SetStyle("Plain");
67
68   AliCDBManager* man = AliCDBManager::Instance();
69   man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
70   man->SetSpecificStorage("GRP/GRP/Data",
71                           Form("local://%s",gSystem->pwd()));
72   man->SetSpecificStorage("ITS/Align/Data",
73                           Form("local://%s",gSystem->pwd()));
74   man->SetSpecificStorage("ITS/Calib/RecoParam",
75                           Form("local://%s",gSystem->pwd()));
76   man->SetRun(0);
77
78
79   gAlice=NULL;
80   AliRunLoader* runLoader = AliRunLoader::Open("galice.root");
81   runLoader->LoadgAlice();
82
83   gAlice = runLoader->GetAliRun();
84
85   runLoader->LoadHeader();
86   runLoader->LoadKinematics();
87   runLoader->LoadRecPoints();
88   runLoader->LoadSDigits();
89   runLoader->LoadHits();
90
91   AliLoader *dl = runLoader->GetDetectorLoader("ITS");
92
93   AliGeomManager::LoadGeometry("geometry.root");
94   TObjArray algITS;
95   AliGeomManager::LoadAlignObjsFromCDBSingleDet("ITS",algITS);
96   AliGeomManager::ApplyAlignObjsToGeom(algITS);
97   //
98   AliITSUGeomTGeo* gm = new AliITSUGeomTGeo(kTRUE);
99   AliITSUClusterPix::SetGeom(gm);
100   //
101   AliITSURecoDet *its = new AliITSURecoDet(gm, "ITSinterface");
102   its->CreateClusterArrays();
103   //
104   Double_t xg1,yg1,zg1=0.,xg0,yg0,zg0=0.,tg0;
105   //
106   TTree * cluTree = 0x0;
107   TTree *hitTree = 0x0;
108   TClonesArray *hitList=new TClonesArray("AliITSUHit");
109   //
110   TObjArray arrMCTracks; // array of hit arrays for each particle
111   //
112   Float_t xyzClGloF[3];
113   Double_t xyzClGlo[3],xyzClTr[3];
114   Int_t labels[3];
115   int nLab = 0;
116   int nlr=its->GetNLayersActive();
117   int ntotev = (Int_t)runLoader->GetNumberOfEvents();
118   printf("N Events : %i \n",ntotev);
119   if (nev>0) ntotev = TMath::Min(nev,ntotev);
120   //
121   // output tree
122   TFile* flOut = TFile::Open("clInfo.root","recreate");
123   TTree* trOut = new TTree("clitsu","clitsu");
124   clSumm cSum;
125   trOut->Branch("evID", &cSum.evID ,"evID/I");
126   trOut->Branch("volID",&cSum.volID,"volID/I");
127   trOut->Branch("lrID", &cSum.lrID ,"lrID/I");  
128   trOut->Branch("clID", &cSum.clID ,"clID/I");  
129   trOut->Branch("nPix", &cSum.nPix ,"nPix/I");
130   trOut->Branch("nX"  , &cSum.nX   ,"nX/I");
131   trOut->Branch("nZ"  , &cSum.nZ   ,"nZ/I");
132   trOut->Branch("q"   , &cSum.q    ,"q/I");
133   trOut->Branch("pt"  , &cSum.pt   ,"pt/F");  
134   trOut->Branch("eta"  ,&cSum.eta  ,"eta/F");  
135   trOut->Branch("phi"  , &cSum.phi  ,"phi/F");  
136   trOut->Branch("xyz",   cSum.xyz,  "xyz[3]/F");  
137   trOut->Branch("dX"  , &cSum.dX   ,"dX/F");
138   trOut->Branch("dZ"  , &cSum.dZ   ,"dZ/F");  
139   trOut->Branch("split",&cSum.split,"split/O");
140   trOut->Branch("prim", &cSum.prim, "prim/O");
141   trOut->Branch("pdg",  &cSum.pdg,  "pdg/I");
142   trOut->Branch("ntr",  &cSum.ntr,  "ntr/I");
143   //
144   for (Int_t iEvent = 0; iEvent < ntotev; iEvent++) {
145     printf("\n Event %i \n",iEvent);
146     runLoader->GetEvent(iEvent);
147     AliStack *stack = runLoader->Stack();
148     cluTree=dl->TreeR();
149     hitTree=dl->TreeH();
150     hitTree->SetBranchAddress("ITS",&hitList);
151     // 
152     // read clusters
153     for (int ilr=nlr;ilr--;) {
154       TBranch* br = cluTree->GetBranch(Form("ITSRecPoints%d",ilr));
155       if (!br) {printf("Did not find cluster branch for lr %d\n",ilr); exit(1);}
156       br->SetAddress(its->GetLayerActive(ilr)->GetClustersAddress());
157     }
158     cluTree->GetEntry(0); 
159     its->ProcessClusters();
160     //
161     // read hits
162     for(Int_t iEnt=0;iEnt<hitTree->GetEntries();iEnt++){//entries loop degli hits
163       hitTree->GetEntry(iEnt);
164       int nh = hitList->GetEntries();      
165       for(Int_t iHit=0; iHit<nh;iHit++){
166         AliITSUHit *pHit = (AliITSUHit*)hitList->At(iHit);
167         int mcID = pHit->GetTrack();
168         TClonesArray* harr = arrMCTracks.GetEntriesFast()>mcID ? (TClonesArray*)arrMCTracks.At(mcID) : 0;
169         if (!harr) {
170           harr = new TClonesArray("AliITSUHit"); // 1st encounter of the MC track
171           arrMCTracks.AddAtAndExpand(harr,mcID);
172         }
173         //
174         new ( (*harr)[harr->GetEntriesFast()] ) AliITSUHit(*pHit);
175       }
176     }
177     //
178     // compare clusters and hits
179     //
180     printf(" tree entries: %lld\n",cluTree->GetEntries());
181     //
182     for (int ilr=0;ilr<nlr;ilr++) {
183       AliITSURecoLayer* lr = its->GetLayerActive(ilr);
184       TClonesArray* clr = lr->GetClusters();
185       int nClu = clr->GetEntries();
186       printf("Layer %d : %d clusters\n",ilr,nClu);
187       //
188       for (int icl=0;icl<nClu;icl++) {
189         AliITSUClusterPix *cl = (AliITSUClusterPix*)clr->At(icl);
190         int modID = cl->GetVolumeId();
191         
192         //------------ check if this is a split cluster
193         int sInL = modID - gm->GetFirstChipIndex(ilr);
194         if (!cl->TestBit(kSplCheck)) {
195           cl->SetBit(kSplCheck);
196           // check if there is no other cluster with same label on this module
197           AliITSURecoSens* sens = lr->GetSensor(sInL);
198           int nclSn = sens->GetNClusters();
199           int offs = sens->GetFirstClusterId();
200           //    printf("To check for %d (mod:%d) N=%d from %d\n",icl,modID,nclSn,offs);
201           for (int ics=0;ics<nclSn;ics++) {
202             AliITSUClusterPix* clusT = (AliITSUClusterPix*)lr->GetCluster(offs+ics); // access to clusters
203             if (clusT==cl) continue;
204             for (int ilb0=0;ilb0<3;ilb0++) {
205               int lb0 = cl->GetLabel(ilb0); if (lb0<=-1) break;
206               for (int ilb1=0;ilb1<3;ilb1++) {
207                 int lb1 = clusT->GetLabel(ilb1); if (lb1<=-1) break;
208                 if (lb1==lb0) {
209                   cl->SetBit(kSplit);
210                   clusT->SetBit(kSplit);
211                   /*
212                     printf("Discard clusters of module %d:\n",modID);
213                     cl->Print();
214                     clusT->Print();
215                   */
216                   break;
217                 }
218               }
219             }
220           }
221         }
222         //------------
223         const AliITSsegmentation* segm = gm->GetSegmentation(ilr);
224         //
225         cl->GetGlobalXYZ(xyzClGloF);
226         int clsize = cl->GetNPix();
227         for (int i=3;i--;) xyzClGlo[i] = xyzClGloF[i];
228         const TGeoHMatrix* mat = gm->GetMatrixSens(modID);
229         if (!mat) {printf("failed to get matrix for module %d\n",cl->GetVolumeId());}
230         mat->MasterToLocal(xyzClGlo,xyzClTr);
231         //
232         int col,row;
233         segm->LocalToDet(xyzClTr[0],xyzClTr[2],row,col); // effective col/row
234         nLab = 0;
235         for (int il=0;il<3;il++) {
236           if (cl->GetLabel(il)>=0) labels[nLab++] = cl->GetLabel(il);
237           else break;
238         }
239         // find hit info
240         for (int il=0;il<nLab;il++) {
241           TClonesArray* htArr = (TClonesArray*)arrMCTracks.At(labels[il]);
242           if (!htArr) {printf("did not find MChits for label %d ",labels[il]); cl->Print(); continue;}
243           //
244           int nh = htArr->GetEntriesFast();
245           AliITSUHit *pHit=0;
246           for (int ih=nh;ih--;) {
247             AliITSUHit* tHit = (AliITSUHit*)htArr->At(ih);
248             if (tHit->GetChip()!=modID) continue;
249             pHit = tHit;
250             break;
251           }
252           if (!pHit) {
253             printf("did not find MChit for label %d on module %d ",il,modID); 
254             cl->Print(); 
255             htArr->Print();
256             continue;
257           }
258           //
259           pHit->GetPositionG(xg1,yg1,zg1);
260           pHit->GetPositionG0(xg0,yg0,zg0,tg0);
261           //
262           double txyzH[3],gxyzH[3] = { (xg1+xg0)/2, (yg1+yg0)/2, (zg1+zg0)/2 };
263           mat->MasterToLocal(gxyzH,txyzH);
264           double rcl = TMath::Sqrt(xyzClTr[0]*xyzClTr[0]+xyzClTr[1]*xyzClTr[1]);
265           double rht = TMath::Sqrt(txyzH[0]*txyzH[0]+txyzH[1]*txyzH[1]);
266           //
267           GetHistoClSize(clsize,kDR,&histoArr)->Fill((rht-rcl)*1e4);
268           if (cl->TestBit(kSplit)) {
269             if (col%2) GetHistoClSize(clsize,kDTXoddSPL,&histoArr)->Fill((txyzH[0]-xyzClTr[0])*1e4);
270             else       GetHistoClSize(clsize,kDTXevenSPL,&histoArr)->Fill((txyzH[0]-xyzClTr[0])*1e4);
271             GetHistoClSize(clsize,kDTZSPL,&histoArr)->Fill((txyzH[2]-xyzClTr[2])*1e4);
272             GetHistoClSize(0,kNPixSPL,&histoArr)->Fill(clsize);
273           }
274           if (col%2) GetHistoClSize(clsize,kDTXodd,&histoArr)->Fill((txyzH[0]-xyzClTr[0])*1e4);
275           else       GetHistoClSize(clsize,kDTXeven,&histoArr)->Fill((txyzH[0]-xyzClTr[0])*1e4);
276           GetHistoClSize(clsize,kDTZ,&histoArr)->Fill((txyzH[2]-xyzClTr[2])*1e4);
277           GetHistoClSize(0,kNPixAll,&histoArr)->Fill(clsize);
278           //
279           cSum.evID = iEvent;
280           cSum.volID = cl->GetVolumeId();
281           cSum.lrID = ilr;
282           cSum.clID = icl;
283           cSum.nPix = cl->GetNPix();
284           cSum.nX   = cl->GetNx();
285           cSum.nZ   = cl->GetNz();
286           cSum.q    = cl->GetQ();
287           cSum.split = cl->TestBit(kSplit);
288           cSum.dX = (txyzH[0]-xyzClTr[0])*1e4;
289           cSum.dZ = (txyzH[2]-xyzClTr[2])*1e4;
290           int label = cl->GetLabel(0);
291           TParticle* part = 0;
292           if (label>=0 && (part=stack->Particle(label)) ) {
293             cSum.pdg = part->GetPdgCode();
294             cSum.eta = part->Eta();
295             cSum.pt  = part->Pt();
296             cSum.phi = part->Phi();
297             cSum.prim = stack->IsPhysicalPrimary(label);
298           }
299           cSum.ntr = 0;
300           for (int ilb=0;ilb<3;ilb++) if (cl->GetLabel(ilb)>=0) cSum.ntr++;
301           for (int i=0;i<3;i++) cSum.xyz[i] = xyzClGloF[i];
302           //
303           trOut->Fill();
304           /*
305           if (clsize==5) {
306             printf("\nL%d(%c) Mod%d, Cl:%d | %+5.1f %+5.1f (%d/%d)|H:%e %e %e | C:%e %e %e\n",ilr,cl->TestBit(kSplit) ? 'S':'N',
307                    modID,icl,(txyzH[0]-xyzClTr[0])*1e4,(txyzH[2]-xyzClTr[2])*1e4, row,col,
308                    gxyzH[0],gxyzH[1],gxyzH[2],xyzClGlo[0],xyzClGlo[1],xyzClGlo[2]);
309             cl->Print();
310             pHit->Print();
311             //
312             double a0,b0,c0,a1,b1,c1,e0;
313             pHit->GetPositionL0(a0,b0,c0,e0);
314             pHit->GetPositionL(a1,b1,c1);
315             float cloc[3];
316             cl->GetLocalXYZ(cloc);
317             printf("LocH: %e %e %e | %e %e %e\n",a0,b0,c0,a1,b1,c1);
318             printf("LocC: %e %e %e | %e %e %e\n",cloc[0],cloc[1],cloc[2],xyzClTr[0],xyzClTr[1],xyzClTr[2]);
319           }
320           */
321           //
322         }
323       }
324     }
325     
326     //    layerClus.Clear();
327     //
328     arrMCTracks.Delete();
329   }//event loop
330   //
331   flOut->cd();
332   trOut->Write();
333   delete trOut;
334   flOut->Close();
335   flOut->Delete();
336   DrawReport("clinfo.ps",&histoArr);
337   //
338 }
339
340 void DrawReport(const char* psname, TObjArray* harr) 
341 {
342   gStyle->SetOptFit(1);
343   if (!harr) harr = &histoArr;
344   TCanvas* cnv = new TCanvas("cl","cl",900,600);
345   //
346   TString psnm1 = psname;
347   if (psnm1.IsNull()) psnm1 = "clusters.ps";
348   TString psnm0 = psnm1.Data(); 
349   psnm0 += "[";
350   TString psnm2 = psnm1.Data(); 
351   psnm2 += "]";
352   cnv->Print(psnm0.Data());
353   //
354   TH1* clall = GetHistoClSize(0,kNPixAll,harr);
355   clall->SetLineColor(kRed);
356   clall->Draw();
357   TH1* clSpl = GetHistoClSize(0,kNPixSPL,harr);
358   clSpl->SetLineColor(kBlue);
359   clSpl->Draw("sames");
360   gPad->Modified();
361   gPad->Update();
362   SetStPadPos(clall,0.75,0.97,0.8,1.,-1,clall->GetLineColor());
363   SetStPadPos(clSpl,0.75,0.97,0.6,0.8,-1,clSpl->GetLineColor());
364   gPad->Modified();
365   gPad->Update();
366   gPad->SetLogy(1);
367   //
368   cnv->cd();
369   cnv->Print(psnm1.Data());
370   //
371   // plot cluster sized from 1 to 10
372   for (int i=1;i<=10;i++) {
373     if (clall->GetBinContent(clall->FindBin(i))<100) continue;
374     DrawNP(i,harr,cnv);
375     cnv->Print(psnm1.Data());
376   }
377   cnv->Print(psnm2.Data());
378 }
379
380 //------------------------------------------
381 TH1* GetHistoClSize(int npix,int id,TObjArray* harr)
382 {
383   // book histos
384   TH1* h = 0;
385   if (!harr) harr = &histoArr;
386   //
387   if (npix<1) {
388     if (harr->GetEntriesFast()>=id && (h=(TH1*)harr->At(id))) return h;
389     h = new TH1F("npixAll","npixAll",150,0.5,54.5); 
390     h->SetDirectory(0);
391     h->SetLineColor(kRed);
392     harr->AddAtAndExpand(h, kNPixAll);
393     //
394     h = new TH1F("npixSpl","npixSpl",150,0.5,54.5);
395     h->SetLineColor(kBlue);
396     h->SetDirectory(0);
397     harr->AddAtAndExpand(h, kNPixSPL);
398     //
399     h = (TH1*)harr->At(id);
400     if (!h) {printf("Unknown histo id=%d\n",id); exit(1);}
401     return h;
402   }
403   //
404   int idh = npix*10+id;
405   if (harr->GetEntriesFast()>=idh && (h=(TH1*)harr->At(idh))) return h;
406   //
407   const int nbin=100;
408   const double kdiff=80;
409   // need to create set of histos
410   //
411   h = new TH1F(Form("dxy_npix%d",npix),Form("dr_npix%d",npix),nbin,-kdiff,kdiff);
412   h->SetDirectory(0);
413   harr->AddAtAndExpand(h, npix*10 + kDR);
414   //
415   h  = new TH1F(Form("dtxODD_npix%d",npix),Form("dtxODD_npix%d",npix),nbin,-kdiff,kdiff);
416   h->SetDirectory(0);
417   h->SetLineColor(kRed);
418   harr->AddAtAndExpand(h, npix*10 + kDTXodd);
419   h  = new TH1F(Form("dtxEVN_npix%d",npix),Form("dtxEVN_npix%d",npix),nbin,-kdiff,kdiff);
420   h->SetDirectory(0);
421   h->SetLineColor(kBlue);
422   harr->AddAtAndExpand(h, npix*10 + kDTXeven);
423   //
424   h  = new TH1F(Form("dtz_npix%d",npix),Form("dtz_npix%d",npix),nbin,-kdiff,kdiff);
425   h->SetLineColor(kGreen);
426   h->SetDirectory(0);
427   harr->AddAtAndExpand(h, npix*10 + kDTZ);
428   //
429   //
430   h  = new TH1F(Form("SPL_dtxODD_npix%d",npix),Form("SPL_dtxODD_npix%d",npix),nbin,-kdiff,kdiff);
431   h->SetLineColor(kMagenta);
432   h->SetFillColor(kMagenta);
433   h->SetFillStyle(3001);  
434   h->SetLineStyle(2);
435   h->SetDirectory(0);
436
437   harr->AddAtAndExpand(h, npix*10 + kDTXoddSPL);
438   h  = new TH1F(Form("SPL_dtxEVN_npix%d",npix),Form("SPL_dtxEVN_npix%d",npix),nbin,-kdiff,kdiff);
439   h->SetLineColor(kCyan);
440   h->SetFillColor(kCyan);
441   h->SetFillStyle(3006);  
442   h->SetLineStyle(2);
443   h->SetDirectory(0);
444   harr->AddAtAndExpand(h, npix*10 + kDTXevenSPL);
445   //
446   h  = new TH1F(Form("SPL_dtz_npix%d",npix),Form("SPLdtz_npix%d",npix),nbin,-kdiff,kdiff);
447   harr->AddAtAndExpand(h, npix*10 + kDTZSPL);
448   h->SetDirectory(0);
449   //
450   h->SetLineColor(kGreen+2);
451   h->SetFillColor(kGreen+2);
452   h->SetLineStyle(2);
453   h->SetFillStyle(3001);
454   h = (TH1*)harr->At(idh);
455   if (!h) {printf("Unknown histo id=%d\n",idh); exit(1);}
456   return h;
457 }
458
459
460
461 TCanvas* DrawNP(int np, TObjArray* harr, TCanvas* cnv)
462 {
463   if (!harr) harr = &histoArr;
464   if (!cnv) cnv = new TCanvas(Form("cnv%d",np),Form("cnv%d",np),900,700);
465   cnv->Clear();
466   cnv->Divide(2,1);
467   cnv->cd(1);
468   //
469   TH1* dxodd = (TH1*)harr->At(np*10+kDTXodd);
470   TH1* dxevn = (TH1*)harr->At(np*10+kDTXeven);
471   TH1* dxoddS =(TH1*)harr->At(np*10+kDTXoddSPL);
472   TH1* dxevnS =(TH1*)harr->At(np*10+kDTXevenSPL);
473   double max = TMath::Max(dxodd->GetMaximum(),dxevn->GetMaximum());
474   dxodd->SetMaximum(1.1*max);
475   dxodd->GetXaxis()->SetTitle("#DeltaX, #mum");
476   dxodd->SetTitle(Form("#DeltaX for clSize=%d",np));
477   dxodd->Fit("gaus","","");
478   dxevn->Fit("gaus","","sames");
479   //
480   dxoddS->Draw("sames");
481   dxevnS->Draw("sames");
482   //
483   gPad->Modified();
484   gPad->Update();
485   SetStPadPos(dxodd,0.75,0.97,0.8,1., -1,dxodd->GetLineColor());
486   SetStPadPos(dxevn,0.75,0.97,0.6,0.8, -1,dxevn->GetLineColor());
487   SetStPadPos(dxoddS,0.75,0.97,0.4,0.6, -1,dxoddS->GetLineColor());
488   SetStPadPos(dxevnS,0.75,0.97,0.2,0.4, -1,dxevnS->GetLineColor());
489   //
490   cnv->cd(2);
491   TH1* dz  = (TH1*)harr->At(np*10+kDTZ);
492   dz->SetTitle(Form("#DeltaZ for clSize=%d",np));
493   dz->GetXaxis()->SetTitle("#DeltaZ, #mum");
494   dz->Fit("gaus");
495   TH1* dzS = (TH1*)harr->At(np*10+kDTZSPL);
496   dz->Draw("sames");
497   gPad->Modified();
498   gPad->Update();
499   SetStPadPos(dz,0.75,0.97,0.8,1., -1, dz->GetLineColor());
500   SetStPadPos(dzS,0.75,0.97,0.5,0.7, -1, dzS->GetLineColor());
501   gPad->Modified();
502   gPad->Update();
503   //
504   cnv->cd();
505   return cnv;
506 }
507
508
509 TPaveStats* GetStPad(TH1* hst)
510 {
511   TList *lst = hst->GetListOfFunctions();
512   if (!lst) return 0;
513   int nf = lst->GetSize();
514   for (int i=0;i<nf;i++) {
515     TPaveStats *fnc = (TPaveStats*) lst->At(i);
516     if (fnc->InheritsFrom("TPaveStats")) return fnc;
517   }
518   return 0;
519   //
520 }
521
522
523 TPaveStats* SetStPadPos(TH1* hst,float x1,float x2,float y1,float y2, Int_t stl, Int_t col)
524 {
525   TPaveStats* pad = GetStPad(hst);
526   if (!pad) return 0;
527   pad->SetX1NDC( x1 );
528   pad->SetX2NDC( x2 );
529   pad->SetY1NDC( y1 );
530   pad->SetY2NDC( y2 );
531   if (stl>=0) pad->SetFillStyle(stl);
532   if (col>=0) pad->SetTextColor(col);
533   pad->SetFillColor(0);
534   //
535   gPad->Modified();
536   return pad;
537 }
538