6 #include "TDirectory.h"
29 float dz=-1;
float dr=-1;
31 if(i > 100 && i < 200) { dz=3.33;dr=0.4;}
33 if(i > 200 && i < 1000 && ( i%10 == 1 || i%10 == 7)) { dz=0.8;dr=0.4;}
34 if(i > 200 && i < 1000 && !( i%10 == 1 || i%10 == 7)) { dz=0.8;dr=0.8;}
36 if(i > 1000 && i < 2000) { dz=5.948;dr=0.4;}
38 if(i > 3000 && i < 4000) { dz=9.440;dr=0.4;}
40 if(i > 2000 && i < 3000 && (i%1000)/100 == 1) { dz=0.8;dr=5.647;}
41 if(i > 2000 && i < 3000 && (i%1000)/100 == 2) { dz=0.8;dr=4.512;}
42 if(i > 2000 && i < 3000 && (i%1000)/100 == 3) { dz=0.8;dr=5.637;}
44 if(i > 4000 && i < 6000 && (i%1000)/100 == 1) { dz=0.8;dr=4.362;}
45 if(i > 4000 && i < 6000 && (i%1000)/100 == 2) { dz=0.8;dr=4.512;}
46 if(i > 4000 && i < 6000 && (i%1000)/100 == 3) { dz=0.8;dr=5.637;}
47 if(i > 4000 && i < 6000 && (i%1000)/100 == 4) { dz=0.8;dr=5.862;}
48 if(i > 4000 && i < 6000 && (i%1000)/100 == 5) { dz=0.8;dr=7.501;}
49 if(i > 4000 && i < 6000 && (i%1000)/100 == 6) { dz=0.8;dr=9.336;}
50 if(i > 4000 && i < 6000 && (i%1000)/100 == 7) { dz=0.8;dr=10.373;}
52 std::pair<float,float> res(dz,dr);
58 float dz=-1;
float dr=-1;
60 if(i > 1000 && i < 1040) { dz=3.33;dr=0.4;}
61 if(i > 1040 && i < 1130) { dr=3.33;dz=0.4;}
63 if(i > 2000 && i < 2550) { dz=2.5;dr=0.1;}
64 if(i > 2550 && i < 3000) { dz=5.0;dr=0.1;}
66 if(i > 3000 && i < 5000) { dz=0.2;dr=5.0;}
68 if(i > 3100 && i < 3119) { dz=0.2;dr=2.5;}
69 if(i > 3140 && i < 3159) { dz=0.2;dr=2.5;}
70 if(i > 3180 && i < 3199) { dz=0.2;dr=2.5;}
71 if(i > 3220 && i < 3239) { dz=0.2;dr=2.5;}
72 if(i > 3260 && i < 3279) { dz=0.2;dr=2.5;}
74 if(i > 4100 && i < 4119) { dz=0.2;dr=2.5;}
75 if(i > 4140 && i < 4159) { dz=0.2;dr=2.5;}
76 if(i > 4180 && i < 4199) { dz=0.2;dr=2.5;}
77 if(i > 4220 && i < 4239) { dz=0.2;dr=2.5;}
78 if(i > 4260 && i < 4279) { dz=0.2;dr=2.5;}
80 std::pair<float,float> res(dz,dr);
86 float dz=-1;
float dr=-1;
88 if(i > 1000 && i < 1040) { dz=3.33;dr=0.4;}
89 if(i > 1040 && i < 1080) { dr=3.33;dz=0.4;}
91 if(i > 1100 && i < 2000) { dz=5.948;dr=0.4;}
93 if(i > 3000 && i < 4000) { dz=9.440;dr=0.4;}
95 if(i > 2000 && i < 3000 && (i%1000)/100 == 1) { dz=0.8;dr=5.647;}
96 if(i > 2000 && i < 3000 && (i%1000)/100 == 2) { dz=0.8;dr=4.512;}
97 if(i > 2000 && i < 3000 && (i%1000)/100 == 3) { dz=0.8;dr=5.637;}
99 if(i > 4000 && i < 6000 && (i%1000)/100 == 1) { dz=0.8;dr=4.362;}
100 if(i > 4000 && i < 6000 && (i%1000)/100 == 2) { dz=0.8;dr=4.512;}
101 if(i > 4000 && i < 6000 && (i%1000)/100 == 3) { dz=0.8;dr=5.637;}
102 if(i > 4000 && i < 6000 && (i%1000)/100 == 4) { dz=0.8;dr=5.862;}
103 if(i > 4000 && i < 6000 && (i%1000)/100 == 5) { dz=0.8;dr=7.501;}
104 if(i > 4000 && i < 6000 && (i%1000)/100 == 6) { dz=0.8;dr=9.336;}
105 if(i > 4000 && i < 6000 && (i%1000)/100 == 7) { dz=0.8;dr=10.373;}
107 std::pair<float,float> res(dz,dr);
113 h->SetAxisRange(min,max);
118 TText*
t =
new TText((max+min)/2,h->GetMaximum(),
label); t->SetTextAlign(22);
124 std::pair<float,float>(*
size)(
int), std::vector<SubDetParams>& vsub) {
126 gROOT->SetStyle(
"Plain");
130 TProfile* aveoccu= (TProfile*)gDirectory->Get(
"aveoccu");
131 TProfile* avemult= (TProfile*)gDirectory->Get(
"avemult");
132 TH1F* nchannels = (TH1F*)gDirectory->Get(
"nchannels_real");
134 TProfile* averadius = (TProfile*)gDirectory->Get(
"averadius");
135 TProfile* avez = (TProfile*)gDirectory->Get(
"avez");
137 std::cout <<
"pointers " << aveoccu <<
" " << avemult <<
" " << nchannels <<
" " << averadius <<
" " << avez << std::endl;
139 if(aveoccu && avemult && nchannels && averadius && avez) {
142 for(
int i=1;
i<nchannels->GetNbinsX()+1;++
i) {
143 nchannels->SetBinError(
i,0.);
146 TH1D* haveoccu = aveoccu->ProjectionX(
"haveoccu");
147 haveoccu->SetDirectory(0);
148 haveoccu->Divide(nchannels);
150 TH1D* havemult = avemult->ProjectionX(
"havemult");
151 havemult->SetDirectory(0);
152 havemult->Divide(nchannels);
154 TH1D* havewidth = (TH1D*)haveoccu->Clone(
"havewidth");
155 havewidth->SetDirectory(0);
156 havewidth->SetTitle(
"Average Cluster Size");
157 havewidth->Divide(havemult);
160 new TCanvas(
"occupancy",
"occupancy",1200,500);
162 haveoccu->SetStats(0);
163 haveoccu->SetLineColor(kRed);
164 haveoccu->SetMarkerColor(kRed);
165 haveoccu->DrawCopy();
167 new TCanvas(
"multiplicity",
"multiplicity",1200,500);
169 havemult->SetStats(0);
170 havemult->SetLineColor(kRed);
171 havemult->SetMarkerColor(kRed);
172 havemult->DrawCopy();
174 new TCanvas(
"width",
"width",1200,500);
175 havewidth->SetStats(0);
176 havewidth->SetLineColor(kRed);
177 havewidth->SetMarkerColor(kRed);
178 havewidth->DrawCopy();
181 TCanvas* o2 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(
"occupancy2");
184 haveoccu->SetLineColor(kBlue);
185 haveoccu->SetMarkerColor(kBlue);
188 o2 =
new TCanvas(
"occupancy2",
"occupancy2",1200,800);
191 for(
unsigned int isub=0;isub<vsub.size();++isub) {
192 printFrame(o2,haveoccu,vsub[isub].
label.c_str(),isub+1,vsub[isub].min,vsub[isub].max,same);
196 TCanvas* m2 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(
"multiplicity2");
199 havemult->SetLineColor(kBlue);
200 havemult->SetMarkerColor(kBlue);
203 m2 =
new TCanvas(
"multiplicity2",
"multiplicity2",1200,800);
206 for(
unsigned int isub=0;isub<vsub.size();++isub) {
207 printFrame(m2,havemult,vsub[isub].
label.c_str(),isub+1,vsub[isub].min,vsub[isub].max,same);
211 TCanvas*
w2 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(
"width2");
214 havewidth->SetLineColor(kBlue);
215 havewidth->SetMarkerColor(kBlue);
218 w2 =
new TCanvas(
"width2",
"width2",1200,800);
221 for(
unsigned int isub=0;isub<vsub.size();++isub) {
222 printFrame(w2,havewidth,vsub[isub].
label.c_str(),isub+1,vsub[isub].min,vsub[isub].max,same);
225 float (*
scale)(float);
230 drawMap(
"occumap",haveoccu,averadius,avez,min,max,
size,
scale,color,
"channel occupancy");
238 float combinedOccupancy(TFile* ff,
const char* module,
const int lowerbin,
const int upperbin) {
245 TProfile* aveoccu= (TProfile*)gDirectory->Get(
"aveoccu");
247 TH1F* nchannels = (TH1F*)gDirectory->Get(
"nchannels_real");
250 float sumnchannels=0;
253 for(
int i=lowerbin;
i<upperbin+1; ++
i) {
254 std::cout <<
"processing bin " <<
i <<
" " << aveoccu->GetBinContent(
i) <<
"+/-" << aveoccu->GetBinError(
i) << std::endl;
255 sumoccu += aveoccu->GetBinContent(
i);
256 sumnchannels += nchannels->GetBinContent(
i);
257 sumerrsq += aveoccu->GetBinError(
i)*aveoccu->GetBinError(
i);
259 cumoccu = sumnchannels!=0 ? sumoccu/sumnchannels : -1;
260 cumerr = sumnchannels!=0 ?
sqrt(sumerrsq)/sumnchannels : -1;
261 std::cout <<
"Cumulative occupancy: " << sumoccu <<
" " << sumnchannels <<
" " << cumoccu <<
"+/-" << cumerr;
268 void PlotOccupancyMap(TFile* ff,
const char* module,
const float min,
const float max,
const float mmin,
const float mmax,
const int color) {
270 std::pair<float,float> (*size)(int);
273 std::vector<SubDetParams> vsub;
274 SubDetParams ppix={
"BPIX+FPIX",100,270}; vsub.push_back(ppix);
275 SubDetParams ptib={
"TIB",1050,1450}; vsub.push_back(ptib);
276 SubDetParams ptid={
"TID",2070,2400}; vsub.push_back(ptid);
277 SubDetParams ptob={
"TOB",3000,3700}; vsub.push_back(ptob);
278 SubDetParams ptecm={
"TEC-",4000,4850}; vsub.push_back(ptecm);
279 SubDetParams ptecp={
"TEC+",5000,5850}; vsub.push_back(ptecp);
284 void PlotOccupancyMapPhase1(TFile* ff,
const char* module,
const float min,
const float max,
const float mmin,
const float mmax,
const int color) {
286 std::pair<float,float> (*size)(int);
289 std::vector<SubDetParams> vsub;
290 SubDetParams ppix={
"BPIX+FPIX",1000,1080}; vsub.push_back(ppix);
291 SubDetParams ptib={
"TIB",1090,1450}; vsub.push_back(ptib);
292 SubDetParams ptid={
"TID",2070,2400}; vsub.push_back(ptid);
293 SubDetParams ptob={
"TOB",3000,3700}; vsub.push_back(ptob);
294 SubDetParams ptecm={
"TEC-",4000,4850}; vsub.push_back(ptecm);
295 SubDetParams ptecp={
"TEC+",5000,5850}; vsub.push_back(ptecp);
300 void PlotOccupancyMapPhase2(TFile* ff,
const char* module,
const float min,
const float max,
const float mmin,
const float mmax,
const int color) {
302 std::pair<float,float> (*size)(int);
305 std::vector<SubDetParams> vsub;
306 SubDetParams ppix={
"BPIX+FPIX",1000,1090}; vsub.push_back(ppix);
307 SubDetParams ptob={
"TOB",2000,2900}; vsub.push_back(ptob);
308 SubDetParams ptecm={
"TEC-",3100,3300}; vsub.push_back(ptecm);
309 SubDetParams ptecp={
"TEC+",4100,4300}; vsub.push_back(ptecp);
314 void PlotOnTrackOccupancy(TFile* ff,
const char* module,
const char* ontrkmod,
const float mmin,
const float mmax,
const int color) {
316 std::pair<float,float> (*size)(int);
319 std::vector<SubDetParams> vsub;
320 SubDetParams ppix={
"BPIX+FPIX",100,270}; vsub.push_back(ppix);
321 SubDetParams ptib={
"TIB",1050,1450}; vsub.push_back(ptib);
322 SubDetParams ptid={
"TID",2070,2400}; vsub.push_back(ptid);
323 SubDetParams ptob={
"TOB",3000,3700}; vsub.push_back(ptob);
324 SubDetParams ptecm={
"TEC-",4000,4850}; vsub.push_back(ptecm);
325 SubDetParams ptecp={
"TEC+",5000,5850}; vsub.push_back(ptecp);
332 std::pair<float,float> (*size)(int);
335 std::vector<SubDetParams> vsub;
336 SubDetParams ppix={
"BPIX+FPIX",1000,1080}; vsub.push_back(ppix);
337 SubDetParams ptib={
"TIB",1090,1450}; vsub.push_back(ptib);
338 SubDetParams ptid={
"TID",2070,2400}; vsub.push_back(ptid);
339 SubDetParams ptob={
"TOB",3000,3700}; vsub.push_back(ptob);
340 SubDetParams ptecm={
"TEC-",4000,4850}; vsub.push_back(ptecm);
341 SubDetParams ptecp={
"TEC+",5000,5850}; vsub.push_back(ptecp);
348 std::pair<float,float> (*size)(int);
351 std::vector<SubDetParams> vsub;
352 SubDetParams ppix={
"BPIX+FPIX",1000,1090}; vsub.push_back(ppix);
353 SubDetParams ptob={
"TOB",2000,2900}; vsub.push_back(ptob);
354 SubDetParams ptecm={
"TEC-",3100,3300}; vsub.push_back(ptecm);
355 SubDetParams ptecp={
"TEC+",4100,4300}; vsub.push_back(ptecp);
361 std::pair<float,float>(*
size)(
int),
const std::vector<SubDetParams>& vsub) {
364 gROOT->SetStyle(
"Plain");
367 TProfile* aveontrkmult=0;
368 TProfile* averadius =0;
372 avemult= (TProfile*)gDirectory->Get(
"avemult");
373 averadius = (TProfile*)gDirectory->Get(
"averadius");
374 avez = (TProfile*)gDirectory->Get(
"avez");
376 if(ff->cd(ontrkmod)) aveontrkmult= (TProfile*)gDirectory->Get(
"avemult");
378 std::cout <<
"pointers " << avemult <<
" " << aveontrkmult <<
" " << averadius <<
" " << avez << std::endl;
380 if( averadius && avez && avemult && aveontrkmult) {
382 TH1D* havemult = avemult->ProjectionX(
"havemult");
383 TH1D* haveontrkmult = aveontrkmult->ProjectionX(
"haveontrkmult");
384 havemult->SetDirectory(0);
385 haveontrkmult->SetDirectory(0);
386 haveontrkmult->Divide(havemult);
388 new TCanvas(
"ontrkmult",
"ontrkmult",1200,500);
390 haveontrkmult->SetStats(0);
391 haveontrkmult->SetLineColor(kRed);
392 haveontrkmult->SetMarkerColor(kRed);
393 haveontrkmult->SetMarkerSize(.5);
394 haveontrkmult->SetMarkerStyle(20);
395 haveontrkmult->DrawCopy();
398 TCanvas* o2 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(
"ontrkmult2");
401 haveontrkmult->SetLineColor(kBlue);
402 haveontrkmult->SetMarkerColor(kBlue);
405 o2 =
new TCanvas(
"ontrkmult2",
"ontrkmult2",1200,800);
408 for(
unsigned int isub=0;isub<vsub.size();++isub) {
409 printFrame(o2,haveontrkmult,vsub[isub].
label.c_str(),isub+1,vsub[isub].min,vsub[isub].max,same);
412 float (*
scale)(float);
415 drawMap(
"ontrkmultmap",haveontrkmult,averadius,avez,mmin,mmax,
size,
scale,color);
419 TCanvas*
drawMap(
const char* cname,
const TH1* hval,
const TProfile* averadius,
const TProfile* avez,
const float mmin,
const float mmax,
420 std::pair<float,float>(*
size)(
int),
float(*
scale)(
float),
const int color,
const char* ptitle) {
424 const Int_t NRGBs = 5;
425 const Int_t NCont = 255;
426 Double_t stops[NRGBs] = { 0.00, 0.25, 0.50, 0.75, 1.00 };
427 Double_t
red[NRGBs] = { 0.00, 0.00, 0.40, 1.00, 1.00 };
428 Double_t
green[NRGBs] = { 0.00, 0.40, 0.70, 0.60, 1.00 };
429 Double_t blue[NRGBs] = { 0.30, 0.60, 0.00, 0.00, 0.20 };
430 TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
431 gStyle->SetNumberContours(NCont);
435 const Int_t NRGBs = 3;
436 const Int_t NCont = 255;
437 Double_t stops[NRGBs] = { 0.00, 0.50, 1.00 };
438 Double_t
red[NRGBs] = { 0.90, 0.50, 0.00};
439 Double_t
green[NRGBs] = { 0.90, 0.50, 0.00};
440 Double_t blue[NRGBs] = { 0.90, 0.50, 0.00};
441 TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
442 gStyle->SetNumberContours(NCont);
446 const Int_t NRGBs = 7;
447 const Int_t NCont = 255;
448 Double_t stops[NRGBs] = { 0.00, 0.15, 0.30, 0.45, 0.65, 0.85, 1.00 };
449 Double_t
red[NRGBs] = { 0.60, 0.30, 0.00, 0.00, 0.60, 0.40, 0.00 };
450 Double_t
green[NRGBs] = { 1.00, 0.90, 0.80, 0.75, 0.20, 0.00, 0.00 };
451 Double_t blue[NRGBs] = { 1.00, 1.00, 1.00, 0.30, 0.00, 0.00, 0.00 };
452 TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
453 gStyle->SetNumberContours(NCont);
456 int ncol = gStyle->GetNumberOfColors();
457 std::cout <<
"Number of colors " << ncol << std::endl;
462 for(
int i=1;
i<hval->GetNbinsX();++
i) {
464 if(averadius->GetBinEntries(
i)*avez->GetBinEntries(
i)) {
473 if(dz<0 && dr<0)
continue;
476 TBox* modmult =
new TBox(avez->GetBinContent(
i)-dz,averadius->GetBinContent(
i)-dr,avez->GetBinContent(
i)+dz,averadius->GetBinContent(
i)+dr);
477 modmult->SetFillStyle(1001);
479 modmult->SetFillColor(kBlack);
484 if(icol > (ncol-1)) icol=(ncol-1);
485 std::cout <<
i <<
" " << icol <<
" " << hval->GetBinContent(
i) << std::endl;
486 modmult->SetFillColor(gStyle->GetColorPalette(icol));
488 modulesmult.Add(modmult);
495 double etavalext[] = {4.,3.5,3.,2.8,2.6,2.4,2.2,2.0,1.8,1.6};
496 double etavalint[] = {-1.4,-1.2,-1.0,-0.8,-0.6,-0.4,-0.2,0.,0.2,0.4,0.6,0.8,1.0,1.2,1.4};
500 for(
int i=0;
i<10;++
i) {
502 double eta = etavalext[
i];
503 TLine* lin =
new TLine(295,2*295/(
exp(eta)-
exp(-eta)),305,2*305/(
exp(eta)-
exp(-eta)));
506 sprintf(lab,
"%3.1f",eta);
507 TText*
label =
new TText(285,2*285/(
exp(eta)-
exp(-eta)),lab);
508 label->SetTextSize(.03);
509 label->SetTextAlign(22);
510 etalabels.Add(label);
512 for(
int i=0;
i<10;++
i) {
514 double eta = -1*etavalext[
i];
515 TLine* lin =
new TLine(-295,-2*295/(
exp(eta)-
exp(-eta)),-305,-2*305/(
exp(eta)-
exp(-eta)));
518 sprintf(lab,
"%3.1f",eta);
519 TText*
label =
new TText(-285,-2*285/(
exp(eta)-
exp(-eta)),lab);
520 label->SetTextSize(.03);
521 label->SetTextAlign(22);
522 etalabels.Add(label);
524 for(
int i=0;
i<15;++
i) {
526 double eta = etavalint[
i];
527 TLine* lin =
new TLine(130.*(
exp(eta)-
exp(-eta))/2.,130,138.*(
exp(eta)-
exp(-eta))/2.,138);
530 sprintf(lab,
"%3.1f",eta);
531 TText*
label =
new TText(125.*(
exp(eta)-
exp(-eta))/2.,125,lab);
532 label->SetTextSize(.03);
533 label->SetTextAlign(22);
534 etalabels.Add(label);
536 TLatex* etalab =
new TLatex(0,115,
"#eta");
537 etalab->SetTextSize(.03);
538 etalab->SetTextAlign(22);
539 etalabels.Add(etalab);
542 TLatex *cmslab =
new TLatex(0.15,0.965,
"CMS");
544 cmslab->SetTextSize(0.04);
545 cmslab->SetTextAlign(31);
546 paperlabels.Add(cmslab);
547 TLatex *enelab =
new TLatex(0.92,0.965,
"#sqrt{s} = 7 TeV");
549 enelab->SetTextSize(0.04);
550 enelab->SetTextAlign(31);
551 paperlabels.Add(enelab);
561 TGaxis *raxis =
new TGaxis(-310,0,-310,140,0,140,10,
"S");
562 TGaxis *zaxis =
new TGaxis(-310,0,310,0,-310,310,10,
"S");
563 raxis->SetTickSize(.01); zaxis->SetTickSize(.01);
564 raxis->SetTitle(
"r (cm)"); zaxis->SetTitle(
"z (cm)");
568 for(
int i = 0;
i< ncol ; ++
i) {
569 TBox* box=
new TBox(315,0+140./ncol*
i,330,0+140./ncol*(i+1));
570 box->SetFillStyle(1001);
571 box->SetFillColor(gStyle->GetColorPalette(i));
578 mpaxis =
new TGaxis(330,0,330,140,mmin,mmax,510,
"SLG+");
581 mpaxis =
new TGaxis(330,0,330,140,mmin,mmax,510,
"SL+");
583 mpaxis->SetTickSize(.02);
584 mpaxis->SetLabelOffset(mpaxis->GetLabelOffset()*0.5);
585 mpaxis->SetTitle(ptitle);
586 mpalette.Add(mpaxis);
588 TCanvas* cc2 =
new TCanvas(cname,cname,1000,500);
589 cc2->Range(-370.,-20.,390.,150.);
590 TFrame* fr2 =
new TFrame(-310,0,310,140);
591 fr2->UseCurrentStyle();
593 raxis->Draw(); zaxis->Draw();
594 std::cout << modulesmult.GetSize() << std::endl;
597 if(color>=0) mpalette.Draw();
607 gROOT->SetStyle(
"Plain");
609 TCanvas* cc =
new TCanvas(name,name,750,750);
610 cc->Range(-25,-25,25,25);
611 TFrame* fr1 =
new TFrame(-20,-20,20,20);
612 fr1->UseCurrentStyle();
616 TProfile* avex = (TProfile*)gDirectory->Get(
"avex");
617 TProfile* avey = (TProfile*)gDirectory->Get(
"avey");
618 TProfile* avez = (TProfile*)gDirectory->Get(
"avez");
620 if(avex && avey && avez) {
621 TText* tittext =
new TText(0,0,name);
622 tittext->SetTextSize(.04); tittext->SetTextAlign(22);
624 for(
unsigned int mod=ioffset+1;
mod<ioffset+57;++
mod) {
625 double x = avex->GetBinContent(
mod);
626 double y = avey->GetBinContent(
mod);
629 sprintf(modstring,
"%d",
mod%100);
630 TText* modtext =
new TText(x,y,modstring);
631 modtext->SetTextAngle(atan(y/x)*180/3.14159);
632 modtext->SetTextSize(.02); modtext->SetTextAlign(22); modtext->SetTextColor(kRed);
637 for(
unsigned int mod=ioffset+101;
mod<ioffset+157;++
mod) {
638 double x = avex->GetBinContent(
mod);
639 double y = avey->GetBinContent(
mod);
642 sprintf(modstring,
"%d",
mod%100);
643 TText* modtext =
new TText(x,y,modstring);
644 modtext->SetTextAngle(atan(y/x)*180/3.14159);
645 modtext->SetTextSize(.02); modtext->SetTextAlign(22); modtext->SetTextColor(kBlue);
646 std::cout <<
mod <<
" " << x <<
" " << y <<
" " << atan(y/x) << std::endl;
656 gROOT->SetStyle(
"Plain");
660 TProfile* averadius = (TProfile*)gDirectory->Get(
"averadius");
661 TProfile* avez = (TProfile*)gDirectory->Get(
"avez");
663 std::cout <<
"pointers " << averadius <<
" " << avez << std::endl;
665 if(averadius && avez) {
667 std::pair<float,float> (*size)(int);
669 float (*
scale)(float);
673 drawMap(
"trackermap",averadius,averadius,avez,0,0,size,
scale,-1);
void PlotDebugFPIX_XYMap(TFile *ff, const char *module, const unsigned int ioffset, const char *name)
void PlotOnTrackOccupancy(TFile *ff, const char *module, const char *ontrkmod, const float mmin, const float mmax, const int color)
void PlotOccupancyMapPhase1(TFile *ff, const char *module, const float min, const float max, const float mmin, const float mmax, const int color)
common ppss p3p6s2 common epss epspn46 common const1 w2
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
std::pair< float, float > phase2bin(int i)
void PlotOccupancyMap(TFile *ff, const char *module, const float min, const float max, const float mmin, const float mmax, const int color)
std::pair< float, float > presentbin(int i)
void PlotOccupancyMapGeneric(TFile *ff, const char *module, const float min, const float max, const float mmin, const float mmax, const int color, std::pair< float, float >(*size)(int), std::vector< SubDetParams > &vsub)
void printFrame(TCanvas *c, TH1D *h, const char *label, const int frame, const int min, const int max, const bool same)
U second(std::pair< T, U > const &p)
void PlotTrackerXsect(TFile *ff, const char *module)
T x() const
Cartesian x coordinate.
std::pair< float, float > phase1bin(int i)
float combinedOccupancy(TFile *ff, const char *module, const int lowerbin, const int upperbin)
TCanvas * drawMap(const char *cname, const TH1 *hval, const TProfile *averadius, const TProfile *avez, const float mmin, const float mmax, std::pair< float, float >(*size)(int), float(*scale)(float), const int color, const char *ptitle)
void PlotOnTrackOccupancyPhase2(TFile *ff, const char *module, const char *ontrkmod, const float mmin, const float mmax, const int color)
void PlotOccupancyMapPhase2(TFile *ff, const char *module, const float min, const float max, const float mmin, const float mmax, const int color)
void PlotOnTrackOccupancyPhase1(TFile *ff, const char *module, const char *ontrkmod, const float mmin, const float mmax, const int color)
T mod(const T &a, const T &b)
tuple size
Write out results.
void PlotOnTrackOccupancyGeneric(TFile *ff, const char *module, const char *ontrkmod, const float mmin, const float mmax, const int color, std::pair< float, float >(*size)(int), const std::vector< SubDetParams > &vsub)