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