CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
compareAlignments.cc
Go to the documentation of this file.
1 #include <string.h>
2 #include <cstring>
3 #include "TChain.h"
4 #include "TFile.h"
5 #include "TH1.h"
6 #include "TTree.h"
7 #include "TKey.h"
8 #include "Riostream.h"
9 #include <vector>
10 #include <sstream>
11 #include "TCanvas.h"
12 #include "TLegend.h"
13 #include "TROOT.h"
14 #include "TPaveStats.h"
15 #include "TObjArray.h"
16 #include "TObjString.h"
17 #include "TStyle.h"
18 #include "TEnv.h"
19 
20 
21 TList *FileList;
22 TList *LabelList;
23 TFile *Target;
24 std::vector< std::string > lowestlevels;
25 std::vector< int > theColors;
26 std::vector< int > theStyles;
27 
28 void MergeRootfile( TDirectory *target, TList *sourcelist, TList *labellist );
29 void nicePad(Int_t logx,Int_t logy);
30 void SetMinMaxRange(TObjArray *hists);
31 
32 void ColourStatsBoxes(TObjArray *hists);
33 
34 void compareAlignments(TString namesandlabels="readFromFile")
35 {
36  cout << "Comparing using: >"<<namesandlabels<<"<"<<endl;
37 
38  gStyle->SetOptStat(111110);
39  gStyle->SetTitleFillColor(10);
40  gStyle->SetTitleBorderSize(0);
41 
42  Target = TFile::Open( "result.root", "RECREATE" );
43  FileList = new TList();
44  LabelList = new TList();
45 
46  int formatCounter = 1;
47  //TObjArray* stringarray = namesandlabels.Tokenize(",");
48  TObjArray *nameandlabelpairs = namesandlabels.Tokenize(",");
49  for (Int_t i = 0; i < nameandlabelpairs->GetEntries(); ++i) {
50  TObjArray *aFileLegPair = TString(nameandlabelpairs->At(i)->GetName()).Tokenize("=");
51 
52  if(aFileLegPair->GetEntries() == 2) {
53  TFile* currentFile = TFile::Open(aFileLegPair->At(0)->GetName());
54  if( currentFile != NULL && !currentFile->IsZombie() ){
55  FileList->Add( currentFile ); // 2
56  if(TString(aFileLegPair->At(1)->GetName()).Contains("|")){
57  TObjArray* formatedLegendEntry = TString(aFileLegPair->At(1)->GetName()).Tokenize("|");
58  LabelList->Add( formatedLegendEntry->At(0) );
59  if(formatedLegendEntry->GetEntries() > 1){
60  theColors.push_back(atoi(formatedLegendEntry->At(1)->GetName()));
61 
62  if(formatedLegendEntry->GetEntries() > 2)
63  theStyles.push_back(atoi(formatedLegendEntry->At(2)->GetName()));
64  else
65  theStyles.push_back( formatCounter );
66  }else{
67  std::cout <<"if you give a \"|\" in the legend name you will need to at least give a int for the color"<<std::endl;
68  }
69  formatCounter++;
70  }else{
71  LabelList->Add( aFileLegPair->At(1) );
72  theColors.push_back(formatCounter);
73  theStyles.push_back(formatCounter);
74  formatCounter++;
75  }
76  }else{
77  std::cout << "Could not open: "<<aFileLegPair->At(0)->GetName()<<std::endl;
78  }
79  }
80  else {
81  std::cout << "Please give file name and legend entry in the following form:\n"
82  << " filename1=legendentry1,filename2=legendentry2[|color[|style]]"<<std::endl;
83 
84  }
85 
86  }
87 
88  // ************************************************************
89  // List of Files
90  //FileList->Add( TFile::Open("../test/AlignmentValidation_Elliptical.root") );
91  //FileList->Add( TFile::Open("../test/AlignmentValidation_10pb.root") ); // 2
92  //FileList->Add( TFile::Open("../test/AlignmentValidation_custom.root") ); // 2
93  // ************************************************************
94 
95  // put here the lowest level up to which you want to combine the
96  // histogramms
97  lowestlevels.push_back("TPBLadder");
98  lowestlevels.push_back("TPEPanel");
99  lowestlevels.push_back("TIBHalfShell");
100  lowestlevels.push_back("TIDRing");
101  lowestlevels.push_back("TOBRod");
102  lowestlevels.push_back("TECSide");
103 // lowestlevels.push_back("Det");
104 
105 
106 
108 
109 }
110 
111 void MergeRootfile( TDirectory *target, TList *sourcelist, TList *labellist ) {
112 
113  if( sourcelist->GetSize() == 0){
114  std::cout<< "Cowardly refuse to merge empty SourceList! " <<std::endl;
115  return;
116  }
117 
118  TString path( (char*)strstr( target->GetPath(), ":" ) );
119  path.Remove( 0, 2 );
120 
121  TFile *first_source = (TFile*)sourcelist->First();
122  TObjString *first_label = (TObjString*)labellist->First();
123 
124  first_source->cd( path );
125  TDirectory *current_sourcedir = gDirectory;
126  //gain time, do not add the objects in the list in memory
127  Bool_t status = TH1::AddDirectoryStatus();
128  TH1::AddDirectory(kFALSE);
129 
130  // loop over all keys in this directory
131 
132  TIter nextkey( current_sourcedir->GetListOfKeys() );
133  TKey *key, *oldkey=0;
134  while ( (key = (TKey*)nextkey())) {
135 
136  //keep only the highest cycle number for each key
137  if (oldkey && !strcmp(oldkey->GetName(),key->GetName())) continue;
138 
139  // read object from first source file
140  first_source->cd( path );
141  TObject *obj = key->ReadObj();
142 
143  if ( obj->IsA()->InheritsFrom( TH1::Class() ) ) {
144  // descendant of TH1 -> merge it
145  TCanvas c(obj->GetName(),obj->GetName(),500,500);
146  c.SetFillColor(10);
147 
148  bool is2d = false;
149  if(strstr(obj->ClassName() ,"TH2") != NULL )
150  is2d = true;
151  TH1 *h1 = static_cast<TH1*>(obj);
152 
153  int q = 1;
154  TObjArray *histarray = new TObjArray;
155 
156  h1->SetLineStyle(theStyles.at(q-1));
157  h1->SetLineWidth(2);
158 
159  h1->SetLineColor(theColors.at(q-1));
160  h1->GetYaxis()->SetTitleOffset(1.5);
161  if(strstr(h1->GetName(),"summary") != NULL )
162  h1->Draw("x0e1*H");
163  else if(is2d)
164  h1->Draw();
165  else
166  h1->Draw();
167 
168  TLegend leg(0.2,0.85,0.775,0.93);
169  leg.AddEntry(h1,first_label->String().Data(),"L");
170  leg.SetBorderSize(0);
171  leg.SetFillColor(10);
172  // loop over all source files and add the content of the
173  // correspondant histogram to the one pointed to by "h1"
174  TFile *nextsource = (TFile*)sourcelist->After( first_source );
175  TObjString *nextlabel = (TObjString*)labellist->After( labellist->First() );
176 
177  histarray->Add(h1);
178  while ( nextsource ) {
179 
180  // make sure we are at the correct directory level by cd'ing to path
181 
182  nextsource->cd( path );
183  TKey *key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(h1->GetName());
184  if (key2) {
185  ++q;
186  TH1 *h2 = (TH1*)key2->ReadObj();
187 
188  if(!is2d){
189  h2->SetLineStyle(theStyles.at(q-1));
190  h2->SetLineWidth(2);
191  }
192 
193  h2->SetLineColor(theColors.at(q-1));
194  std::stringstream newname;
195  newname << h2->GetName() << q;
196 
197  h2->SetName(newname.str().c_str());
198  if(strstr(newname.str().c_str(),"summary") != NULL )
199  h2->DrawClone("x0*He1sames");
200  else if(is2d)
201  h2->DrawClone("sames");
202  else
203  h2->DrawClone("sames");
204  leg.AddEntry(c.FindObject(h2->GetName()),nextlabel->String().Data(),"L");
205  histarray->Add(c.FindObject(h2->GetName()));
206  delete h2;
207  } else {
208  std::cerr << "Histogram "<< key2->GetTitle() << " is not present in file " << nextsource->GetName() << std::endl;
209  }
210 
211  nextsource = (TFile*)sourcelist->After( nextsource );
212  nextlabel = (TObjString*)labellist->After(nextlabel);
213  }
214  nicePad(0,0);
215  leg.Draw();
216  c.Update();
217  if(strstr(h1->GetName(),"summary") == NULL )
218  SetMinMaxRange(histarray);
219  ColourStatsBoxes(histarray);
220  target->cd();
221  c.Write();
222  histarray->Delete();
223 
224 
225 
226 
227  } else if ( obj->IsA()->InheritsFrom( TDirectory::Class() ) ) {
228  // it's a subdirectory
229 
230  std::string dirname = obj->GetName();
231  for( std::vector< std::string >::const_iterator lowlevelit = lowestlevels.begin(),
232  lowlevelitend = lowestlevels.end(); lowlevelit != lowlevelitend; ++lowlevelit)
233  if( dirname.find(*lowlevelit) != std::string::npos )
234  return;
235 
236  // create a new subdir of same name and title in the target file
237  target->cd();
238  TDirectory *newdir = target->mkdir( obj->GetName(), obj->GetTitle() );
239 
240  // newdir is now the starting point of another round of merging
241  // newdir still knows its depth within the target file via
242  // GetPath(), so we can still figure out where we are in the recursion
243  MergeRootfile( newdir, sourcelist, labellist );
244 
245 
246  } else {
247 
248  // object is of no type that we know or can handle
249  cout << "Unknown object type, name: "
250  << obj->GetName() << " title: " << obj->GetTitle() << endl;
251  }
252 
253 
254  } // while ( ( TKey *key = (TKey*)nextkey() ) )
255 
256  // save modifications to target file
257  target->SaveSelf(kTRUE);
258  TH1::AddDirectory(status);
259 }
260 
261 
262 
263 
264 void nicePad(Int_t logx,Int_t logy)
265 {
266  gPad->SetBottomMargin(0.10);
267  gPad->SetRightMargin(0.1);
268  gPad->SetLeftMargin(0.15);
269  gPad->SetTopMargin(0.15);
270  gPad->SetTickx(1);
271  gPad->SetTicky(1);
272  if(logy==1)
273  {
274  gPad->SetLogy();
275  }
276  else
277  {
278  gPad->SetLogy(0);
279  }
280  if(logx==1)
281  {
282  gPad->SetLogx();
283  }
284  else
285  {
286  gPad->SetLogx(0);
287  }
288 }
289 
290 
291 void ColourStatsBoxes(TObjArray *hists)
292 {
293 
294  Double_t fStatsX1 = 0.85, fStatsX2 = 1., fStatsY1 = 0.85, fStatsY2 = 1.;
295  // colours stats boxes like hists' line colors and moves the next to each other
296  if (!hists) return;
297  Double_t x1 = fStatsX1, x2 = fStatsX2, y1 = fStatsY1, y2 = fStatsY2;
298  for (Int_t iH = 0; iH < hists->GetEntries(); ++iH) {
299  TH1 *h = static_cast<TH1*>(hists->At(iH));
300  if (!h) continue;
301  TObject *statObj = h->GetListOfFunctions()->FindObject("stats");
302  if (statObj && statObj->InheritsFrom(TPaveStats::Class())) {
303  TPaveStats *stats = static_cast<TPaveStats*>(statObj);
304  stats->SetLineColor(static_cast<TH1*>(hists->At(iH))->GetLineColor());
305  stats->SetTextColor(static_cast<TH1*>(hists->At(iH))->GetLineColor());
306  stats->SetFillColor(10);
307  stats->SetX1NDC(x1);
308  stats->SetX2NDC(x2);
309  stats->SetY1NDC(y1);
310  stats->SetY2NDC(y2);
311  y2 = y1 - 0.005; // shift down 2
312  y1 = y2 - (fStatsY2 - fStatsY1); // shift down 1
313  if (y1 < 0.) {
314  y1 = fStatsY1; y2 = fStatsY2; // restart y-positions
315  x2 = x1 - 0.005; // shift left 2
316  x1 = x2 - (fStatsX2 - fStatsX1); // shift left 1
317  if (x1 < 0.) { // give up, start again:
318  x1 = fStatsX1, x2 = fStatsX2, y1 = fStatsY1, y2 = fStatsY2;
319  }
320  }
321  //} else if (gStyle->GetOptStat() != 0) { // failure in case changed in list via TExec....
322  //this->Warning("ColourStatsBoxes", "No stats found for %s", hists->At(iH)->GetName());
323  }
324  }
325 }
326 
327 
328 void SetMinMaxRange(TObjArray *hists)
329 {
330  Double_t min = 100000;
331  Double_t max = -100000;
332  for (Int_t iH = 0; iH < hists->GetEntries(); ++iH) {
333  TH1 *h = static_cast<TH1*>(hists->At(iH));
334  if (!h) continue;
335  for(int i = 1; i <= h->GetNbinsX(); ++i) {
336  if(h->GetBinContent(i) + h->GetBinError(i) > max ) max = h->GetBinContent(i) + h->GetBinError(i);
337  if(h->GetBinContent(i) - h->GetBinError(i) < min ) min = h->GetBinContent(i) - h->GetBinError(i);
338  }
339  }
340 
341  TH1 *h_first = static_cast<TH1*>(hists->At(0));
342  h_first->SetMaximum(max+max*0.1);
343  if(min = 0.) {
344  min = -1111;
345  h_first->SetMinimum(min);
346  } else {
347  h_first->SetMinimum(min-min*0.1);
348  }
349 }
void compareAlignments(TString namesandlabels="readFromFile")
std::vector< int > theStyles
int i
Definition: DBlmapReader.cc:9
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
void MergeRootfile(TDirectory *target, TList *sourcelist, TList *labellist)
std::vector< int > theColors
#define NULL
Definition: scimark2.h:8
void ColourStatsBoxes(TObjArray *hists)
TList * LabelList
T min(T a, T b)
Definition: MathUtil.h:58
TFile * Target
string key
FastSim: produces sample of signal events, overlayed with premixed minbias events.
void SetMinMaxRange(TObjArray *hists)
std::vector< std::string > lowestlevels
tuple cout
Definition: gather_cfg.py:145
tuple status
Definition: ntuplemaker.py:245
void nicePad(Int_t logx, Int_t logy)
TList * FileList