CMS 3D CMS Logo

Functions | Variables
compareAlignments.cc File Reference
#include <string.h>
#include <cstring>
#include "TChain.h"
#include "TFile.h"
#include "TH1.h"
#include "TTree.h"
#include "TKey.h"
#include "TMath.h"
#include "Riostream.h"
#include <vector>
#include <sstream>
#include "TCanvas.h"
#include "TLegend.h"
#include "TROOT.h"
#include "TPaveStats.h"
#include "TObjArray.h"
#include "TObjString.h"
#include "TStyle.h"
#include "TEnv.h"
#include "Alignment/OfflineValidation/plugins/TkAlStyle.cc"

Go to the source code of this file.

Functions

void ColourStatsBoxes (TObjArray *hists)
 
void compareAlignments (TString namesandlabels="readFromFile", TString legendheader="", TString lefttitle="", TString righttitle="", bool bigtext=false)
 
void MergeRootfile (TDirectory *target, TList *sourcelist, TList *labellist, bool bigtext)
 
void nicePad (Int_t logx, Int_t logy)
 
void SetMinMaxRange (TObjArray *hists)
 

Variables

TList * FileList
 
TList * LabelList
 
std::vector< std::string > lowestlevels
 
std::vector< int > phases
 
TFile * Target
 
std::vector< int > theColors
 
std::vector< int > theStyles
 

Function Documentation

void ColourStatsBoxes ( TObjArray *  hists)

Definition at line 350 of file compareAlignments.cc.

References fftjetcommon_cfi::Class.

Referenced by MergeRootfile().

351 {
352 
353  Double_t fStatsX1 = 0.85, fStatsX2 = 1., fStatsY1 = 0.77, fStatsY2 = 0.92;
354  // colours stats boxes like hists' line colors and moves the next to each other
355  if (!hists) return;
356  Double_t x1 = fStatsX1, x2 = fStatsX2, y1 = fStatsY1, y2 = fStatsY2;
357  for (Int_t iH = 0; iH < hists->GetEntries(); ++iH) {
358  TH1 *h = static_cast<TH1*>(hists->At(iH));
359  if (!h) continue;
360  TObject *statObj = h->GetListOfFunctions()->FindObject("stats");
361  if (statObj && statObj->InheritsFrom(TPaveStats::Class())) {
362  TPaveStats *stats = static_cast<TPaveStats*>(statObj);
363  stats->SetLineColor(static_cast<TH1*>(hists->At(iH))->GetLineColor());
364  stats->SetTextColor(static_cast<TH1*>(hists->At(iH))->GetLineColor());
365  stats->SetFillColor(10);
366  stats->SetX1NDC(x1);
367  stats->SetX2NDC(x2);
368  stats->SetY1NDC(y1);
369  stats->SetY2NDC(y2);
370  y2 = y1 - 0.005; // shift down 2
371  y1 = y2 - (fStatsY2 - fStatsY1); // shift down 1
372  if (y1 < 0.) {
373  y1 = fStatsY1; y2 = fStatsY2; // restart y-positions
374  x2 = x1 - 0.005; // shift left 2
375  x1 = x2 - (fStatsX2 - fStatsX1); // shift left 1
376  if (x1 < 0.) { // give up, start again:
377  x1 = fStatsX1, x2 = fStatsX2, y1 = fStatsY1, y2 = fStatsY2;
378  }
379  }
380  stats->DrawClone();
381  //} else if (gStyle->GetOptStat() != 0) { // failure in case changed in list via TExec....
382  //this->Warning("ColourStatsBoxes", "No stats found for %s", hists->At(iH)->GetName());
383  }
384  }
385 }
void compareAlignments ( TString  namesandlabels = "readFromFile",
TString  legendheader = "",
TString  lefttitle = "",
TString  righttitle = "",
bool  bigtext = false 
)

Definition at line 38 of file compareAlignments.cc.

References gather_cfg::cout, CUSTOM, FileList, i, LabelList, TkAlStyle::legendheader, lowestlevels, MergeRootfile(), NONE, NULL, phases, TkAlStyle::set(), Target, theColors, and theStyles.

39 {
40  cout << "Comparing using: >"<<namesandlabels<<"<"<<endl;
41 
42  TkAlStyle::legendheader = legendheader;
43  TkAlStyle::set(CUSTOM, NONE, lefttitle, righttitle);
44  gStyle->SetOptStat(111110);
45  gStyle->SetOptTitle(0);
46 
47  Target = TFile::Open( "result.root", "RECREATE" );
48  FileList = new TList();
49  LabelList = new TList();
50 
51  int formatCounter = 1;
52  //TObjArray* stringarray = namesandlabels.Tokenize(",");
53  TObjArray *nameandlabelpairs = namesandlabels.Tokenize(",");
54  for (Int_t i = 0; i < nameandlabelpairs->GetEntries(); ++i) {
55  TObjArray *aFileLegPair = TString(nameandlabelpairs->At(i)->GetName()).Tokenize("=");
56 
57  if(aFileLegPair->GetEntries() == 2) {
58  TFile* currentFile = TFile::Open(aFileLegPair->At(0)->GetName());
59  if( currentFile != NULL && !currentFile->IsZombie() ){
60  FileList->Add( currentFile ); // 2
61  if ( currentFile->Get("TrackerOfflineValidationStandalone/Pixel/P1PXBBarrel_1") ) {
62  cout << "phase 1" << endl;
63  phases.push_back(1);
64  } else if ( currentFile->Get("TrackerOfflineValidationStandalone/Pixel/TPBBarrel_1") ) {
65  cout << "phase 0" << endl;
66  phases.push_back(0);
67  } else {
68  cout << "Unknown phase for file " << aFileLegPair->At(0)->GetName() << endl;
69  assert(false);
70  }
71  if(TString(aFileLegPair->At(1)->GetName()).Contains("|")){
72  TObjArray* formatedLegendEntry = TString(aFileLegPair->At(1)->GetName()).Tokenize("|");
73  LabelList->Add( formatedLegendEntry->At(0) );
74  if(formatedLegendEntry->GetEntries() > 1){
75  theColors.push_back(atoi(formatedLegendEntry->At(1)->GetName()));
76 
77  if(formatedLegendEntry->GetEntries() > 2)
78  theStyles.push_back(atoi(formatedLegendEntry->At(2)->GetName()));
79  else
80  theStyles.push_back( formatCounter );
81  }else{
82  std::cout <<"if you give a \"|\" in the legend name you will need to at least give a int for the color"<<std::endl;
83  }
84  formatCounter++;
85  }else{
86  LabelList->Add( aFileLegPair->At(1) );
87  theColors.push_back(formatCounter);
88  theStyles.push_back(formatCounter);
89  formatCounter++;
90  }
91  }else{
92  std::cout << "Could not open: "<<aFileLegPair->At(0)->GetName()<<std::endl;
93  }
94  }
95  else {
96  std::cout << "Please give file name and legend entry in the following form:\n"
97  << " filename1=legendentry1,filename2=legendentry2[|color[|style]]"<<std::endl;
98 
99  }
100 
101  }
102 
103  // ************************************************************
104  // List of Files
105  //FileList->Add( TFile::Open("../test/AlignmentValidation_Elliptical.root") );
106  //FileList->Add( TFile::Open("../test/AlignmentValidation_10pb.root") ); // 2
107  //FileList->Add( TFile::Open("../test/AlignmentValidation_custom.root") ); // 2
108  // ************************************************************
109 
110  // put here the lowest level up to which you want to combine the
111  // histograms
112  lowestlevels.push_back("TPBLadder");
113  lowestlevels.push_back("TPEPanel");
114  lowestlevels.push_back("TIBHalfShell");
115  lowestlevels.push_back("TIDRing");
116  lowestlevels.push_back("TOBRod");
117  lowestlevels.push_back("TECSide");
118  // phase 1
119  // it checks each one independently so no harm in having
120  // both phase 0 and 1 together in the vector
121  lowestlevels.push_back("P1PXBLadder");
122  lowestlevels.push_back("P1PXECPanel");
123 
124 
125  MergeRootfile( Target, FileList, LabelList, bigtext );
126 
127 }
std::vector< int > theStyles
int i
Definition: DBlmapReader.cc:9
std::vector< int > theColors
#define NULL
Definition: scimark2.h:8
static void set(const PublicationStatus status, const Era era=NONE, const TString customTitle="", const TString customRightTitle="")
Definition: TkAlStyle.cc:401
void MergeRootfile(TDirectory *target, TList *sourcelist, TList *labellist, bool bigtext)
static TString legendheader
Definition: TkAlStyle.cc:102
TList * LabelList
TFile * Target
std::vector< int > phases
std::vector< std::string > lowestlevels
TList * FileList
void MergeRootfile ( TDirectory *  target,
TList *  sourcelist,
TList *  labellist,
bool  bigtext 
)

Definition at line 129 of file compareAlignments.cc.

References EnergyCorrector::c, MessageLogger_cfi::cerr, fftjetcommon_cfi::Class, ColourStatsBoxes(), gather_cfg::cout, compare_using_db::dirname, TkAlStyle::drawStandardTitle(), createfilelist::int, crabWrapper::key, create_public_lumi_plots::leg, TkAlStyle::legendheader, lowestlevels, hpstanc_transforms::max, nicePad(), NULL, MuonAssociatorByHits_cfi::obj, callgraph::path, phases, lumiQueryAPI::q, Scenarios_cff::scale, SetMinMaxRange(), mps_update::status, AlCaHLTBitMon_QueryRunRegistry::string, TkAlStyle::textSize, theColors, and theStyles.

Referenced by compareAlignments().

129  {
130 
131  if( sourcelist->GetSize() == 0){
132  std::cout<< "Cowardly refuse to merge empty SourceList! " <<std::endl;
133  return;
134  }
135 
136  TString path( (char*)strstr( target->GetPath(), ":" ) );
137  path.Remove( 0, 2 );
138 
139  TFile *first_source = (TFile*)sourcelist->First();
140  TObjString *first_label = (TObjString*)labellist->First();
141 
142  first_source->cd( path );
143  TDirectory *current_sourcedir = gDirectory;
144  //gain time, do not add the objects in the list in memory
145  Bool_t status = TH1::AddDirectoryStatus();
146  TH1::AddDirectory(kFALSE);
147 
148  // loop over all keys in this directory
149 
150  TIter nextkey( current_sourcedir->GetListOfKeys() );
151  TKey *key, *oldkey=0;
152  while ( (key = (TKey*)nextkey())) {
153 
154  //keep only the highest cycle number for each key
155  if (oldkey && !strcmp(oldkey->GetName(),key->GetName())) continue;
156 
157  // read object from first source file
158  first_source->cd( path );
159  TObject *obj = key->ReadObj();
160 
161  auto itphase = phases.begin();
162  int firstfilephase = *itphase;
163 
164  if ( obj->IsA()->InheritsFrom( TH1::Class() ) ) {
165  // descendant of TH1 -> merge it
166  TCanvas c(obj->GetName(),obj->GetName());
167  c.SetFillColor(10);
168 
169  bool is2d = false;
170  if(strstr(obj->ClassName() ,"TH2") != NULL )
171  is2d = true;
172  TH1 *h1 = static_cast<TH1*>(obj);
173 
174  int q = 1;
175  TObjArray *histarray = new TObjArray;
176 
177  h1->SetLineStyle(theStyles.at(q-1));
178  h1->SetLineWidth(2);
179 
180  h1->SetTitle("");
181 
182  h1->SetLineColor(theColors.at(q-1));
183  h1->GetYaxis()->SetTitleOffset(1.5);
184  if(strstr(h1->GetName(),"summary") != NULL )
185  h1->Draw("x0e1*H");
186  else if(is2d)
187  h1->Draw();
188  else
189  h1->Draw();
190 
191  double max = h1->GetMaximum();
192  double scale = 1;
193  if (max > 1000) {
194  int power = (int)(TMath::Log10(max / 100));
195  scale = 1/TMath::Power(10, power);
196  h1->GetYaxis()->SetTitle(((TString(h1->GetYaxis()->GetTitle()) + " [#times 10^{") += (int)power) + "}]");
197  }
198  h1->Scale(scale);
199 
200  int nPlots = sourcelist->GetSize();
201  double legendY = 0.80;
202  if (nPlots > 3) { legendY -= 0.01 * (nPlots - 3); }
203  if (bigtext) { legendY -= 0.05; }
204  if (legendY < 0.6) {
205  std::cerr << "Warning: Huge legend!" << std::endl;
206  legendY = 0.6;
207  }
208  TLegend leg(0.17, legendY, 0.85, 0.88);
209  bool hasheader = (TkAlStyle::legendheader != "");
210  if (hasheader) leg.SetHeader(TkAlStyle::legendheader);
211  if (bigtext) leg.SetTextSize(TkAlStyle::textSize);
212  leg.AddEntry(h1,first_label->String().Data(),"L");
213  leg.SetBorderSize(0);
214  leg.SetFillColor(10);
215  // loop over all source files and add the content of the
216  // correspondant histogram to the one pointed to by "h1"
217  TFile *nextsource = (TFile*)sourcelist->After( first_source );
218  TObjString *nextlabel = (TObjString*)labellist->After( labellist->First() );
219 
220  histarray->Add(h1);
221  while ( nextsource ) {
222 
223  TKey *key2;
224  bool wrongphase = false;
225  ++itphase;
226 
227  if (firstfilephase != *itphase && path.Contains("TrackerOfflineValidationStandalone/Pixel")) {
228  //skip this one
229  key2 = 0;
230  wrongphase = true;
231  } else {
232  // make sure we are at the correct directory level by cd'ing to path
233  nextsource->cd( path );
234  key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(h1->GetName());
235  }
236  if (key2) {
237  ++q;
238  TH1 *h2 = (TH1*)key2->ReadObj();
239 
240  if(!is2d){
241  h2->SetLineStyle(theStyles.at(q-1));
242  h2->SetLineWidth(2);
243  }
244 
245  h2->SetLineColor(theColors.at(q-1));
246  std::stringstream newname;
247  newname << h2->GetName() << q;
248  h2->Scale(scale);
249 
250  h2->SetName(newname.str().c_str());
251  if(strstr(newname.str().c_str(),"summary") != NULL )
252  h2->DrawClone("x0*He1sames");
253  else if(is2d)
254  h2->DrawClone("sames");
255  else
256  h2->DrawClone("sames");
257  leg.AddEntry(c.FindObject(h2->GetName()),nextlabel->String().Data(),"L");
258  histarray->Add(c.FindObject(h2->GetName()));
259  delete h2;
260  } else if (wrongphase) {
261  //nothing
262  } else {
263  std::cerr << "Histogram "<< key2->GetTitle() << " is not present in file " << nextsource->GetName() << std::endl;
264  }
265 
266  nextsource = (TFile*)sourcelist->After( nextsource );
267  nextlabel = (TObjString*)labellist->After(nextlabel);
268  }
269  nicePad(0,0);
270  leg.Draw();
272  c.Update();
273  if(strstr(h1->GetName(),"summary") == NULL )
274  SetMinMaxRange(histarray);
275  ColourStatsBoxes(histarray);
276  target->cd();
277  c.Update();
278  c.Write();
279  histarray->Delete();
280 
281 
282 
283 
284  } else if ( obj->IsA()->InheritsFrom( TDirectory::Class() ) ) {
285  // it's a subdirectory
286 
287  std::string dirname = obj->GetName();
288  for( std::vector< std::string >::const_iterator lowlevelit = lowestlevels.begin(),
289  lowlevelitend = lowestlevels.end(); lowlevelit != lowlevelitend; ++lowlevelit)
290  if( dirname.find(*lowlevelit) != std::string::npos )
291  return;
292 
293  // create a new subdir of same name and title in the target file
294  target->cd();
295  TDirectory *newdir = target->mkdir( obj->GetName(), obj->GetTitle() );
296 
297  // newdir is now the starting point of another round of merging
298  // newdir still knows its depth within the target file via
299  // GetPath(), so we can still figure out where we are in the recursion
300  MergeRootfile( newdir, sourcelist, labellist, bigtext );
301 
302 
303  } else {
304 
305  // object is of no type that we know or can handle
306  cout << "Unknown object type, name: "
307  << obj->GetName() << " title: " << obj->GetTitle() << endl;
308  }
309 
310 
311  } // while ( ( TKey *key = (TKey*)nextkey() ) )
312 
313  // save modifications to target file
314  target->SaveSelf(kTRUE);
315  TH1::AddDirectory(status);
316 }
std::vector< int > theStyles
std::vector< int > theColors
#define NULL
Definition: scimark2.h:8
static double textSize
Definition: TkAlStyle.cc:104
void ColourStatsBoxes(TObjArray *hists)
static void drawStandardTitle()
Definition: TkAlStyle.cc:75
void MergeRootfile(TDirectory *target, TList *sourcelist, TList *labellist, bool bigtext)
static TString legendheader
Definition: TkAlStyle.cc:102
void SetMinMaxRange(TObjArray *hists)
std::vector< int > phases
std::vector< std::string > lowestlevels
void nicePad(Int_t logx, Int_t logy)
void nicePad ( Int_t  logx,
Int_t  logy 
)

Definition at line 321 of file compareAlignments.cc.

Referenced by MergeRootfile().

322 {
323 /*
324  gPad->SetBottomMargin(0.10);
325  gPad->SetRightMargin(0.1);
326  gPad->SetLeftMargin(0.15);
327  gPad->SetTopMargin(0.15);
328  gPad->SetTickx(1);
329  gPad->SetTicky(1);
330 */
331  if(logy==1)
332  {
333  gPad->SetLogy();
334  }
335  else
336  {
337  gPad->SetLogy(0);
338  }
339  if(logx==1)
340  {
341  gPad->SetLogx();
342  }
343  else
344  {
345  gPad->SetLogx(0);
346  }
347 }
void SetMinMaxRange ( TObjArray *  hists)

Definition at line 388 of file compareAlignments.cc.

References i, hpstanc_transforms::max, and min().

Referenced by MergeRootfile().

389 {
390  Double_t min = 100000;
391  Double_t max = -100000;
392  for (Int_t iH = 0; iH < hists->GetEntries(); ++iH) {
393  TH1 *h = static_cast<TH1*>(hists->At(iH));
394  if (!h) continue;
395  for(int i = 1; i <= h->GetNbinsX(); ++i) {
396  if(h->GetBinContent(i) + h->GetBinError(i) > max ) max = h->GetBinContent(i) + h->GetBinError(i);
397  if(h->GetBinContent(i) - h->GetBinError(i) < min ) min = h->GetBinContent(i) - h->GetBinError(i);
398  }
399  }
400 
401  TH1 *h_first = static_cast<TH1*>(hists->At(0));
402  h_first->SetMaximum(max*1.3);
403  if(min == 0.) {
404  min = -1111;
405  h_first->SetMinimum(min);
406  } else {
407  h_first->SetMinimum(min-min*0.1);
408  }
409 }
int i
Definition: DBlmapReader.cc:9
T min(T a, T b)
Definition: MathUtil.h:58

Variable Documentation

TList* FileList

Definition at line 24 of file compareAlignments.cc.

Referenced by compareAlignments().

TList* LabelList

Definition at line 25 of file compareAlignments.cc.

Referenced by compareAlignments().

std::vector< std::string > lowestlevels

Definition at line 27 of file compareAlignments.cc.

Referenced by compareAlignments(), and MergeRootfile().

std::vector< int > phases
TFile* Target

Definition at line 26 of file compareAlignments.cc.

Referenced by compareAlignments(), and HistoData::setResultTarget().

std::vector< int > theColors

Definition at line 28 of file compareAlignments.cc.

Referenced by compareAlignments(), and MergeRootfile().

std::vector< int > theStyles

Definition at line 29 of file compareAlignments.cc.

Referenced by compareAlignments(), and MergeRootfile().