ROOTPWA
plotPhase.C
Go to the documentation of this file.
1 
2 //
3 // Copyright 2009 Sebastian Neubert
4 //
5 // This file is part of rootpwa
6 //
7 // rootpwa is free software: you can redistribute it and/or modify
8 // it under the terms of the GNU General Public License as published by
9 // the Free Software Foundation, either version 3 of the License, or
10 // (at your option) any later version.
11 //
12 // rootpwa is distributed in the hope that it will be useful,
13 // but WITHOUT ANY WARRANTY; without even the implied warranty of
14 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 // GNU General Public License for more details.
16 //
17 // You should have received a copy of the GNU General Public License
18 // along with rootpwa. If not, see <http://www.gnu.org/licenses/>.
19 //
21 //-----------------------------------------------------------
22 // File and Version Information:
23 // $Id$
24 //
25 // Description:
26 // Draws phase angle of (wave A - wave B) from tree.
27 //
28 //
29 // Author List:
30 // Sebastian Neubert TUM (original author)
31 //
32 //
33 //-----------------------------------------------------------
34 
35 
36 #include <iostream>
37 #include <sstream>
38 
39 #include "TGraphErrors.h"
40 #include "TAxis.h"
41 #include "TPad.h"
42 #include "TLegend.h"
43 #include "TList.h"
44 
45 #include "reportingUtils.hpp"
46 #include "plotPhase.h"
47 
48 
49 using namespace std;
50 using namespace rpwa;
51 
52 
53 // signature with wave names
54 TMultiGraph*
55 plotPhase(const unsigned int nmbTrees, // number of fitResult trees
56  TTree** trees, // array of fitResult trees
57  const int waveIndexA, // index of first wave
58  const int waveIndexB, // index of second wave
59  const bool saveEps, // if set, EPS file with name wave ID is created
60  const int* graphColors, // array of colors for graph line and marker
61  const bool drawLegend, // if set legend is drawn
62  const string& graphTitle, // name and title of graph (default is wave IDs)
63  const char* drawOption, // draw option for graph
64  const string& selectExpr, // TTree::Draw() selection expression
65  const string& branchName) // fitResult branch name
66 {
67  for (unsigned int i = 0; i < nmbTrees; ++i)
68  if (!trees[i]) {
69  printErr << "null pointer to tree[" << i << "]. aborting." << endl;
70  return 0;
71  }
72  string waveNameA, waveNameB;
73  {
74  // get wave names (assume same wave set in all bins)
75  fitResult* massBin = new fitResult();
76  trees[0]->SetBranchAddress(branchName.c_str(), &massBin);
77  trees[0]->GetEntry(0);
78  waveNameA = massBin->waveName(waveIndexA).Data();
79  waveNameB = massBin->waveName(waveIndexB).Data();
80  }
81  printInfo << "plotting phase between waves '" << waveNameA << "' [" << waveIndexA << "] "
82  << "and '" << waveNameB << "' [" << waveIndexB << "]" << endl;
83  if (selectExpr != "")
84  cout << " using selection criterion '" << selectExpr << "'" << endl;
85 
86  // create multiGraph
87  TMultiGraph* graph = new TMultiGraph();
88  {
89  if (graphTitle != "") {
90  graph->SetName (graphTitle.c_str());
91  graph->SetTitle(graphTitle.c_str());
92  } else {
93  stringstream name;
94  name << "phase_" << waveNameA << "_" << waveNameB;
95  stringstream title;
96  title << "#Delta#varphi(" << waveNameA << " [" << waveIndexA << "], "
97  << waveNameB << " [" << waveIndexB << "])";
98  graph->SetName(name.str().c_str());
99  graph->SetTitle(title.str().c_str());
100  }
101  }
102 
103  // fill multiGraph
104  for (unsigned int i = 0; i < nmbTrees; ++i) {
105 
106  // get wave index for this tree (assume same wave set in all bins)
107  fitResult* massBin = new fitResult();
108  trees[i]->SetBranchAddress(branchName.c_str(), &massBin);
109  trees[i]->GetEntry(0);
110  const int waveIndexAThisTree = massBin->waveIndex(waveNameA);
111  const int waveIndexBThisTree = massBin->waveIndex(waveNameB);
112  if (waveIndexAThisTree < 0) {
113  printInfo << "cannot find wave '" << waveNameA << "' in tree '" << trees[i]->GetTitle() << "'. "
114  << "skipping." << endl;
115  continue;
116  }
117  if (waveIndexBThisTree < 0) {
118  printInfo << "cannot find wave '" << waveNameB << "' in tree '" << trees[i]->GetTitle() << "'. "
119  << "skipping." << endl;
120  continue;
121  }
122 
123  // build and run TTree::Draw() expression
124  stringstream drawExpr;
125  drawExpr << branchName << ".phase(" << waveIndexAThisTree << "," << waveIndexBThisTree << "):"
126  << branchName << ".phaseErr(" << waveIndexAThisTree << "," << waveIndexBThisTree << "):"
127  << branchName << ".massBinCenter() >> h" << waveIndexAThisTree << "_"
128  << waveIndexBThisTree << "_" << i;
129  cout << " running TTree::Draw() expression '" << drawExpr.str() << "' "
130  << "on tree '" << trees[i]->GetName() << "', '" << trees[i]->GetTitle() << "'" << endl;
131  trees[i]->Draw(drawExpr.str().c_str(), selectExpr.c_str(), "goff");
132 
133  // extract data from TTree::Draw() result and build graph
134  const int nmbBins = trees[i]->GetSelectedRows();
135  vector<double> x(nmbBins), xErr(nmbBins);
136  for (int j = 0; j < nmbBins; ++j) {
137  x [j] = trees[i]->GetV3()[j] * 0.001; // convert mass to GeV
138  xErr[j] = 0;
139  }
140  TGraphErrors* g = new TGraphErrors(nmbBins,
141  &(*(x.begin())), // mass
142  trees[i]->GetV1(), // phase
143  &(*(xErr.begin())), // mass error
144  trees[i]->GetV2()); // phase error
145 
146  // beautify graph
147  stringstream graphName;
148  graphName << graph->GetName() << "_" << i;
149  g->SetName (graphName.str().c_str());
150  g->SetTitle(graphName.str().c_str());
151  g->SetMarkerStyle(21);
152  g->SetMarkerSize(0.5);
153  if (graphColors) {
154  g->SetMarkerColor(graphColors[i]);
155  g->SetLineColor (graphColors[i]);
156  }
157  graph->Add(g);
158  }
159 
160  graph->SetMinimum(-240);
161  graph->SetMaximum( 240);
162  // draw graph
163  graph->Draw(drawOption);
164  graph->GetXaxis()->SetTitle("Mass [GeV]");
165  graph->GetYaxis()->SetTitle("Phase Angle [deg]");
166  gPad->Update();
167 
168  // add legend
169  if (drawLegend && (nmbTrees > 1)) {
170  TLegend* legend = new TLegend(0.65, 0.80, 0.99, 0.99);
171  legend->SetFillColor(10);
172  legend->SetBorderSize(1);
173  legend->SetMargin(0.2);
174  for (unsigned int i = 0; i < nmbTrees; ++i) {
175  TGraph* g = static_cast<TGraph*>(graph->GetListOfGraphs()->At(i));
176  legend->AddEntry(g, trees[i]->GetTitle(), "LPE");
177  }
178  legend->Draw();
179  }
180 
181  // create EPS file
182  if (saveEps)
183  gPad->SaveAs(((string)graph->GetName() + ".eps").c_str());
184 
185  return graph;
186 }