source: branches/MarsISDCBranchBasedOn17887/mhvstime/MHVsTime.cc@ 18066

Last change on this file since 18066 was 12880, checked in by tbretz, 13 years ago
Removed some old workarounds hopefully not necessary anymor with newer root version -- tested with 5.32.00 so far.
File size: 9.7 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Bretz, 11/2003 <mailto:thomas.bretz@epfl.ch>
19!
20! Copyright: MAGIC Software Development, 2000-2012
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MHVsTime
28//
29// Use this class if you want to display any rule vs time (or event number)
30//
31// Axis titles
32// ===========
33//
34// 1) If no other title is given the rule for the y-axis is used.
35// 2) If the MH3 has a non-default title (MH3::SetTitle called)
36// this title is set as the histogram title. It can be used to overwrite
37// the axis titles. For more information see TH1::SetTitle, eg.
38// SetTitle("MyGraph;;Counts");
39// The title for the x-axis is ignored and set automatically (MAKE SURE
40// YOU HAVE TWO SEMICOLON!)
41//
42// eg.
43// MHVsTime hist("MHillas.fAlpha");
44// MHVsTime hist("MPointintPos.GetAbsErr");
45// MHVsTime hist("MPointintPos.GetAbsErr*kRad2Deg");
46//
47// To set a maximum number of data-points (eg. to display the last 20min
48// only) call SetMaxPts(200)
49//
50// SetMaxPts(-1) disables this feature.
51//
52//
53// Class Version 2:
54// ----------------
55// + MData *fData; // Object from which the data is filled
56// - MDataChain *fData; // Object from which the data is filled
57// + MData *fError; // Object from which the error is filled
58// - MDataChain *fError; // Object from which the error is filled
59//
60// Class Version 3:
61// ----------------
62// + Double_t fMinimum; // User defined minimum
63// + Double_t fMaximum; // User defined maximum
64//
65/////////////////////////////////////////////////////////////////////////////
66#include "MHVsTime.h"
67
68#include <ctype.h> // tolower
69#include <fstream>
70
71#include <TPad.h>
72#include <TStyle.h>
73#include <TCanvas.h>
74
75#include <TH1.h>
76#include <TGraphErrors.h>
77
78#include "MLog.h"
79#include "MLogManip.h"
80
81#include "MTime.h"
82#include "MParList.h"
83#include "MDataPhrase.h"
84#include "MRawEvtHeader.h"
85
86ClassImp(MHVsTime);
87
88using namespace std;
89
90const TString MHVsTime::gsDefName = "MHVsTime";
91const TString MHVsTime::gsDefTitle = "Container for a graph vs time/evtnumber";
92
93// --------------------------------------------------------------------------
94//
95// Default constructor. For more informations about a valid rule
96// see MDataPhrase.
97//
98MHVsTime::MHVsTime(const char *rule, const char *error)
99 : fGraph(NULL), fData(NULL), fError(NULL), fScale(1), fMaxPts(-1),
100 fNumEvents(1), fMinimum(-1111), fMaximum(-1111), fUseEventNumber(0)
101{
102 fName = gsDefName;
103 fTitle = gsDefTitle;
104
105 if (!rule)
106 return;
107
108 fData = new MDataPhrase(rule);
109
110 if (error)
111 fError = new MDataPhrase(error);
112
113 fGraph = error ? new TGraphErrors : new TGraph;
114 fGraph->SetMarkerStyle(kFullDotMedium);
115}
116
117// --------------------------------------------------------------------------
118//
119// Deletes the histogram
120//
121MHVsTime::~MHVsTime()
122{
123 if (fGraph)
124 delete fGraph;
125
126 if (fData)
127 delete fData;
128
129 if (fError)
130 delete fError;
131}
132
133// --------------------------------------------------------------------------
134//
135// Return the data members used by the data chain to be used in
136// MTask::AddBranchToList
137//
138TString MHVsTime::GetDataMember() const
139{
140 return fData ? fData->GetDataMember() : (TString)"";
141}
142
143// --------------------------------------------------------------------------
144//
145// PreProcess the MDataPhrase. Create a new TGraph. Delete an old one if
146// already allocated.
147//
148Bool_t MHVsTime::SetupFill(const MParList *plist)
149{
150 if (!fGraph || !fData)
151 {
152 *fLog << err << "ERROR - MHVsTime cannot be used with its default constructor!" << endl;
153 return kFALSE;
154 }
155
156 if (!fData->PreProcess(plist))
157 return kFALSE;
158
159 if (fError && !fError->PreProcess(plist))
160 return kFALSE;
161
162 TString title(fData ? GetRule() : (TString)"Graph");
163 title += " vs ";
164 title += fUseEventNumber ? "Event Number" : "Time";
165
166 fGraph->SetNameTitle(fName, fTitle==gsDefTitle?title:fTitle);
167
168 fMean = 0;
169 fN = 0;
170
171 fMin = FLT_MAX;
172 fMax = -FLT_MAX;
173
174 return kTRUE;
175}
176
177// --------------------------------------------------------------------------
178//
179// Set the name of the TGraph and the MHVsTime container
180//
181void MHVsTime::SetName(const char *name)
182{
183 fGraph->SetName(name);
184 MH::SetName(name);
185}
186
187// --------------------------------------------------------------------------
188//
189// Set the title of the TGraph and the MHVsTime container
190//
191void MHVsTime::SetTitle(const char *title)
192{
193 fGraph->SetTitle(title);
194 MH::SetTitle(title);
195}
196
197// --------------------------------------------------------------------------
198//
199// Set the next data point. If the graph exceeds fMaxPts remove the first
200//
201Int_t MHVsTime::Fill(const MParContainer *par, const Stat_t w)
202{
203 Double_t t = 0;
204 if (fUseEventNumber)
205 {
206 const MRawEvtHeader *h = dynamic_cast<const MRawEvtHeader*>(par);
207 t = h ? h->GetDAQEvtNumber() : fGraph->GetN();
208 }
209 else
210 {
211 const MTime *tm = dynamic_cast<const MTime*>(par);
212 if (!tm)
213 {
214 *fLog << err << dbginf << "No MTime found..." << endl;
215 return kERROR;
216 }
217 // If the time is not valid skip this entry
218 if (!*tm)
219 return kTRUE;
220
221 // Do not fill events with equal time
222 if (*tm==fLast || *tm==MTime())
223 return kTRUE;
224
225 fLast = *tm;
226
227 t = tm->GetAxisTime();
228 }
229
230 const Double_t v = fData->GetValue();
231 const Double_t e = fError ? fError->GetValue() : 0;
232
233 //*fLog << all << "ADD " << v << " " << e << endl;
234
235 fMean += v;
236 fMeanErr += e;
237 fN++;
238
239 if (fN==fNumEvents)
240 {
241 /*
242 if ((fMaxPts>0 && fGraph->GetN()>fMaxPts) || fGraph->IsEditable())
243 {
244 fGraph->RemovePoint(0);
245 fGraph->SetEditable(kFALSE);
246 }*/
247
248 const Double_t val = fMean/fN*fScale;
249
250 fGraph->SetPoint(fGraph->GetN(), t, val);
251
252 if (fError)
253 static_cast<TGraphErrors*>(fGraph)->SetPointError(fGraph->GetN()-1, 0, fMeanErr/fN*fScale);
254
255 fMin = TMath::Min(fMin, val);
256 fMax = TMath::Max(fMax, val);
257
258 fMean = 0;
259 fMeanErr = 0;
260 fN = 0;
261 }
262
263 return kTRUE;
264}
265
266// --------------------------------------------------------------------------
267//
268// Set Minimum and Maximum;
269Bool_t MHVsTime::Finalize()
270{
271 const Double_t add = (fMax-fMin)*0.15;
272
273 if (fMinimum==-1111)
274 fGraph->SetMinimum(fMin-add);
275 if (fMaximum==-1111)
276 fGraph->SetMaximum(fMax+add);
277
278 return kTRUE;
279}
280
281// --------------------------------------------------------------------------
282//
283// This displays the TGraph like you expect it to be (eg. time on the axis)
284// It should also make sure, that the displayed time always is UTC and
285// not local time...
286//
287void MHVsTime::Draw(Option_t *opt)
288{
289 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(fGraph);
290 pad->SetBorderMode(0);
291 pad->SetGrid();
292
293 AppendPad("");
294
295 // This is not done automatically anymore since root 5.12/00
296 // and it is necessary to force a proper update of the axis.
297 TH1 *h = fGraph->GetHistogram();
298
299 if (h)
300 {
301 TAxis *axe = h->GetXaxis();
302 // SetPoint deletes the histogram!
303 if (fUseEventNumber)
304 axe->SetTitle("Event Number");
305 else
306 {
307 axe->SetTitle("Time");
308 axe->SetLabelSize(0.033);
309 axe->SetTimeFormat("%H:%M %F1995-01-01 00:00:00 GMT");
310 axe->SetTimeDisplay(1);
311 }
312 }
313
314 if (TestBit(kIsLogy))
315 gPad->SetLogy();
316
317 // If this is set to early the plot remains empty in root 5.12/00
318 if (fMinimum!=-1111)
319 fGraph->SetMinimum(fMinimum);
320 if (fMaximum!=-1111)
321 fGraph->SetMaximum(fMaximum);
322
323 TString str(opt);
324 if (!str.Contains("A"))
325 str += "A";
326 if (!str.Contains("P"))
327 str += "P";
328 if (str.Contains("same", TString::kIgnoreCase))
329 {
330 str.ReplaceAll("same", "");
331 str.ReplaceAll("A", "");
332 }
333
334 fGraph->Draw(str.Data());
335}
336
337// --------------------------------------------------------------------------
338//
339// Used to rebuild a MHVsTime object of the same type (data members,
340// dimension, ...)
341//
342MParContainer *MHVsTime::New() const
343{
344 MHVsTime *h=new MHVsTime(fData ? (const char*)GetRule() : NULL);
345 h->SetScale(fScale);
346 if (fUseEventNumber)
347 h->SetUseEventNumber();
348 h->SetMaxPts(fMaxPts);
349 return h;
350}
351
352TString MHVsTime::GetRule() const
353{
354 return fData ? fData->GetRule() : (TString)"";
355}
356
357// --------------------------------------------------------------------------
358//
359// Returns the total number of bins in a histogram (excluding under- and
360// overflow bins)
361//
362Int_t MHVsTime::GetNbins() const
363{
364 return fGraph->GetN();
365}
366/*
367TH1 *MHVsTime::GetHist()
368{
369 return fGraph ? fGraph->GetHistogram() : 0;
370}
371
372const TH1 *MHVsTime::GetHist() const
373{
374 return fGraph ? fGraph->GetHistogram() : 0;
375}
376
377TH1 *MHVsTime::GetHistByName(const TString name)
378{
379 return GetHist();
380}
381*/
Note: See TracBrowser for help on using the repository browser.