CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
KalmanAlignmentDataCollector.cc
Go to the documentation of this file.
1 
3 
4 #include "TGraph.h"
5 #include "TNtuple.h"
6 #include "TFile.h"
7 #include "TH1F.h"
8 
9 using namespace std;
10 
11 
13 
14 
16 
17 
19 
20 
22 
23 
25 {
26  //if ( !theDataCollector ) theDataCollector = new KalmanAlignmentDataCollector();
27  return theDataCollector;
28 }
29 
30 
32 {
33  //if ( !theDataCollector ) theDataCollector = new KalmanAlignmentDataCollector();
34  theDataCollector->config( config );
35 }
36 
37 
38 void KalmanAlignmentDataCollector::fillHistogram( string histo_name, float data )
39 {
40  //if ( !theDataCollector ) theDataCollector = new KalmanAlignmentDataCollector();
41  theDataCollector->fillTH1F( histo_name, data );
42 }
43 
44 
45 void KalmanAlignmentDataCollector::fillHistogram( string histo_name, int histo_number, float data )
46 {
47  //if ( !theDataCollector ) theDataCollector = new KalmanAlignmentDataCollector();
48  theDataCollector->fillTH1F( histo_name, histo_number, data );
49 }
50 
51 
52 void KalmanAlignmentDataCollector::fillGraph( string graph_name, float x_data, float y_data )
53 {
54  //if ( !theDataCollector ) theDataCollector = new KalmanAlignmentDataCollector();
55  theDataCollector->fillTGraph( graph_name, x_data, y_data );
56 }
57 
58 
59 void KalmanAlignmentDataCollector::fillGraph( string graph_name, int graph_number, float x_data, float y_data )
60 {
61  //if ( !theDataCollector ) theDataCollector = new KalmanAlignmentDataCollector();
62  theDataCollector->fillTGraph( graph_name, graph_number, x_data, y_data );
63 }
64 
65 
66 void KalmanAlignmentDataCollector::fillNtuple( std::string ntuple_name, float data )
67 {
68  //if ( !theDataCollector ) theDataCollector = new KalmanAlignmentDataCollector();
69  theDataCollector->fillTNtuple( ntuple_name, data );
70 }
71 
72 
74 {
75  //if ( !theDataCollector ) theDataCollector = new KalmanAlignmentDataCollector();
77 }
78 
79 
81 {
82  //if ( !theDataCollector ) theDataCollector = new KalmanAlignmentDataCollector();
83  theDataCollector->writeToTFile( file_name, mode );
84 }
85 
86 
88 {
90 }
91 
92 
94 {
96 }
97 
98 
99 void KalmanAlignmentDataCollector::fillTH1F( string histo_name, float data )
100 {
101  if ( theHistoData.find( histo_name ) == theHistoData.end() )
102  {
103  theHistoData[histo_name] = vector< float > ( 1, data );
104  }
105  else
106  {
107  theHistoData[histo_name].push_back( data );
108  }
109 
110  return;
111 }
112 
113 
114 void KalmanAlignmentDataCollector::fillTH1F( string histo_name, int histo_number, float data )
115 {
116  string full_histo_name = histo_name + toString( histo_number );
117  fillTH1F( full_histo_name, data );
118 
119  return;
120 }
121 
122 
123 void KalmanAlignmentDataCollector::fillTGraph( string graph_name, float x_data, float y_data )
124 {
125  if ( theXGraphData.find( graph_name ) == theXGraphData.end() )
126  {
127  theXGraphData[graph_name] = vector< float > ( 1, x_data );
128  theYGraphData[graph_name] = vector< float > ( 1, y_data );
129  }
130  else
131  {
132  theXGraphData[graph_name].push_back( x_data );
133  theYGraphData[graph_name].push_back( y_data );
134  }
135 
136  return;
137 }
138 
139 
140 void KalmanAlignmentDataCollector::fillTGraph( string graph_name, int graph_number, float x_data, float y_data )
141 {
142  string full_graph_name = graph_name + toString( graph_number );
143  fillTGraph( full_graph_name, x_data, y_data );
144 
145  return;
146 }
147 
148 
149 void KalmanAlignmentDataCollector::fillTNtuple( std::string ntuple_name, float data )
150 {
151  if ( theNtupleData.find( ntuple_name ) == theNtupleData.end() )
152  {
153  theNtupleData[ntuple_name] = vector< float > ( 1, data );
154  }
155  else
156  {
157  theNtupleData[ntuple_name].push_back( data );
158  }
159 }
160 
161 
163 {
164  string fileName = theConfiguration.getUntrackedParameter< string >( "FileName", "KalmanAlignmentData.root" );
165  string fileMode = theConfiguration.getUntrackedParameter< string >( "Mode", "RECREATE" );
166  writeToTFile( fileName, fileMode );
167 }
168 
169 
170 void KalmanAlignmentDataCollector::writeToTFile( string file_name, string mode )
171 {
172  int nBins = theConfiguration.getUntrackedParameter< int >( "NBins", 200 );
173  double xMin = theConfiguration.getUntrackedParameter< double >( "XMin", -10. );
174  double xMax = theConfiguration.getUntrackedParameter< double >( "XMax", 10. );
175 
176  TFile* file = new TFile( file_name.c_str(), mode.c_str() );
177 
178  if ( !theHistoData.empty() )
179  {
180  map< string, vector< float > >::iterator itH = theHistoData.begin();
181  vector< float >::iterator itV;
182  TH1F* tempHisto;
183 
184  while ( itH != theHistoData.end() )
185  {
186  tempHisto = new TH1F( itH->first.c_str(), itH->first.c_str(), nBins, xMin, xMax );
187 
188  itV = itH->second.begin();
189  while ( itV != itH->second.end() ) {
190  tempHisto->Fill( *itV );
191  itV++;
192  }
193 
194  tempHisto->Write();
195  delete tempHisto;
196 
197  ++itH;
198  }
199  }
200 
201  if ( !theXGraphData.empty() )
202  {
203  map< string, vector< float > >::iterator itXG = theXGraphData.begin();
204  map< string, vector< float > >::iterator itYG = theYGraphData.begin();
205 
206  float* xData;
207  float* yData;
208 
209  TGraph* tempGraph;
210 
211  while ( itXG != theXGraphData.end() )
212  {
213  int nData = itXG->second.size();
214 
215  xData = new float[nData];
216  yData = new float[nData];
217 
218  for ( int iData = 0; iData < nData; iData++ )
219  {
220  xData[iData] = itXG->second[iData];
221  yData[iData] = itYG->second[iData];
222  }
223 
224  tempGraph = new TGraph( nData, xData, yData );
225  tempGraph->SetName( itXG->first.c_str() );
226  tempGraph->SetTitle( itXG->first.c_str() );
227 
228  tempGraph->Write();
229  delete tempGraph;
230 
231  delete[] xData;
232  delete[] yData;
233 
234  ++itXG;
235  ++itYG;
236  }
237  }
238 
239 
240  if ( !theNtupleData.empty() )
241  {
242  map< string, vector< float > >::iterator itN = theNtupleData.begin();
243 
244  TNtuple* ntuple;
245 
246  while ( itN != theNtupleData.end() )
247  {
248  ntuple = new TNtuple( itN->first.c_str(), itN->first.c_str(), itN->first.c_str() );
249 
250  vector< float >::iterator itD = itN->second.begin(), itDEnd = itN->second.end();
251  while ( itD != itDEnd )
252  {
253  ntuple->Fill( *itD );
254  ++itD;
255  }
256 
257  ntuple->Write();
258  delete ntuple;
259 
260  ++itN;
261  }
262  }
263 
264 
265  file->Write();
266  file->Close();
267  delete file;
268 
269  return;
270 }
271 
272 
274 {
275  theHistoData.clear();
276  theXGraphData.clear();
277  theYGraphData.clear();
278 }
279 
280 
282 {
283  char temp[10];
284  snprintf( temp, sizeof(temp), "%u", i );
285 
286  return string( temp );
287 }
std::map< std::string, std::vector< float > > theNtupleData
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
static void fillHistogram(std::string histo_name, float data)
void fillTGraph(std::string graph_name, float x_data, float y_data)
void fillTH1F(std::string histo_name, float data)
std::map< std::string, std::vector< float > > theXGraphData
std::map< std::string, std::vector< float > > theHistoData
A simple class that allows fast and easy histograming and the production of graphs.
static KalmanAlignmentDataCollector * theDataCollector
static void fillNtuple(std::string ntuple_name, float data)
void fillTNtuple(std::string ntuple_name, float data)
static KalmanAlignmentDataCollector * get(void)
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
void config(const edm::ParameterSet &config)
static void configure(const edm::ParameterSet &config)
static void fillGraph(std::string graph_name, float x_data, float y_data)
std::map< std::string, std::vector< float > > theYGraphData