source: trunk/MagicSoft/Mars/mhvstime/MHVsTime.cc@ 7971

Last change on this file since 7971 was 7971, checked in by tbretz, 18 years ago
*** empty log message ***
File size: 9.8 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:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2004
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#include "MHVsTime.h"
54
55#include <ctype.h> // tolower
56#include <fstream>
57
58#include <TPad.h>
59#include <TStyle.h>
60#include <TCanvas.h>
61
62#include <TGraphErrors.h>
63
64#include "MLog.h"
65#include "MLogManip.h"
66
67#include "MTime.h"
68#include "MParList.h"
69#include "MDataChain.h"
70#include "MRawEvtHeader.h"
71
72ClassImp(MHVsTime);
73
74using namespace std;
75
76const TString MHVsTime::gsDefName = "MHVsTime";
77const TString MHVsTime::gsDefTitle = "Container for a graph vs time/evtnumber";
78
79// --------------------------------------------------------------------------
80//
81// Default constructor. For more informations about a valid rule
82// see MDataChain.
83//
84MHVsTime::MHVsTime(const char *rule, const char *error)
85 : fGraph(NULL), fData(NULL), fError(NULL), fScale(1), fMaxPts(-1),
86 fNumEvents(1), fUseEventNumber(0)
87{
88 fName = gsDefName;
89 fTitle = gsDefTitle;
90
91 if (!rule)
92 return;
93
94 fData = new MDataChain(rule);
95
96 if (error)
97 fError = new MDataChain(error);
98
99 fGraph = error ? new TGraphErrors : new TGraph;
100 fGraph->SetPoint(0, 0, 0); // Dummy point!
101 fGraph->SetEditable(); // Used as flag: First point? yes/no
102 fGraph->SetMarkerStyle(kFullDotMedium);
103}
104
105// --------------------------------------------------------------------------
106//
107// Deletes the histogram
108//
109MHVsTime::~MHVsTime()
110{
111 if (fGraph)
112 delete fGraph;
113
114 if (fData)
115 delete fData;
116
117 if (fError)
118 delete fError;
119}
120
121// --------------------------------------------------------------------------
122//
123// Call SetMinimum of fGraph
124//
125void MHVsTime::SetMinimum(Double_t min)
126{
127 if (fGraph)
128 fGraph->SetMinimum(min);
129}
130
131// --------------------------------------------------------------------------
132//
133// Return the data members used by the data chain to be used in
134// MTask::AddBranchToList
135//
136TString MHVsTime::GetDataMember() const
137{
138 return fData ? fData->GetDataMember() : (TString)"";
139}
140
141// --------------------------------------------------------------------------
142//
143// PreProcess the MDataChain. Create a new TGraph. Delete an old one if
144// already allocated.
145//
146Bool_t MHVsTime::SetupFill(const MParList *plist)
147{
148 if (!fGraph || !fData)
149 {
150 *fLog << err << "ERROR - MHVsTime cannot be used with its default constructor!" << endl;
151 return kFALSE;
152 }
153
154 if (!fData->PreProcess(plist))
155 return kFALSE;
156
157 if (fError && !fError->PreProcess(plist))
158 return kFALSE;
159
160 fGraph->Set(1);
161 fGraph->SetPoint(0, 0, 0); // Dummy point!
162 fGraph->SetEditable(); // Used as flag: First point? yes/no
163
164 TString title(fData ? GetRule() : (TString)"Graph");
165 title += " vs ";
166 title += fUseEventNumber ? "Event Number" : "Time";
167
168 fGraph->SetNameTitle(fName, fTitle==gsDefTitle?title:fTitle);
169
170 fMean = 0;
171 fN = 0;
172
173 return kTRUE;
174}
175
176// --------------------------------------------------------------------------
177//
178// Set the name of the TGraph and the MHVsTime container
179//
180void MHVsTime::SetName(const char *name)
181{
182 fGraph->SetName(name);
183 MH::SetName(name);
184}
185
186// --------------------------------------------------------------------------
187//
188// Set the title of the TGraph and the MHVsTime container
189//
190void MHVsTime::SetTitle(const char *title)
191{
192 fGraph->SetTitle(title);
193 MH::SetTitle(title);
194}
195
196// --------------------------------------------------------------------------
197//
198// Set the next data point. If the graph exceeds fMaxPts remove the first
199//
200Bool_t MHVsTime::Fill(const MParContainer *par, const Stat_t w)
201{
202 Double_t t = 0;
203 if (fUseEventNumber)
204 {
205 const MRawEvtHeader *h = dynamic_cast<const MRawEvtHeader*>(par);
206 t = h ? h->GetDAQEvtNumber() : fGraph->GetN();
207 }
208 else
209 {
210 const MTime *tm = dynamic_cast<const MTime*>(par);
211 if (!tm)
212 {
213 *fLog << err << dbginf << "No MTime found..." << endl;
214 return kFALSE;
215 }
216 // If the time is not valid skip this entry
217 if (!*tm)
218 return kTRUE;
219
220 // Do not fill events with equal time
221 if (*tm==fLast || *tm==MTime())
222 return kTRUE;
223
224 fLast = *tm;
225
226 t = tm->GetAxisTime();
227 }
228
229 const Double_t v = fData->GetValue();
230 const Double_t e = fError ? fError->GetValue() : 0;
231
232 //*fLog << all << "ADD " << v << " " << e << endl;
233
234 fMean += v;
235 fMeanErr += e;
236 fN++;
237
238 if (fN==fNumEvents)
239 {
240 if (fMaxPts>0 && fGraph->GetN()>fMaxPts || fGraph->IsEditable())
241 {
242 fGraph->RemovePoint(0);
243 fGraph->SetEditable(kFALSE);
244 }
245
246 fGraph->SetPoint(fGraph->GetN(), t, fMean/fN*fScale);
247
248 if (fError)
249 static_cast<TGraphErrors*>(fGraph)->SetPointError(fGraph->GetN()-1, 0, fMeanErr/fN*fScale);
250
251 fMean = 0;
252 fMeanErr = 0;
253 fN = 0;
254 }
255
256 return kTRUE;
257}
258
259void MHVsTime::Paint(Option_t *opt)
260{
261 /*
262 *fLog << all << fGraph << " " << gPad->GetName() << " ----> Paint " << opt << endl;
263 TListIter Next(gPad->GetListOfPrimitives()); TObject *obj;
264 while ((obj=Next())) *fLog << obj << " " << obj->GetName() << " " << obj->ClassName() << " " << Next.GetOption() << endl;
265 */
266 if (!fGraph)
267 return;
268
269// *fLog << all << fGraph->GetN() << " " << opt << endl;
270
271 if (fGraph->GetN()==0)
272 return;
273
274 TString str(opt);
275 if (!str.Contains("A"))
276 str += "A";
277 if (!str.Contains("P"))
278 str += "P";
279 if (str.Contains("same", TString::kIgnoreCase))
280 {
281 str.ReplaceAll("same", "");
282 str.ReplaceAll("A", "");
283 }
284
285 // This is not done automatically anymore since root 5.12/00
286 // and it is necessary to force a proper update of the axis.
287 TH1 *h = fGraph->GetHistogram();
288 if (h)
289 {
290 delete h;
291 fGraph->SetHistogram(0);
292 h = fGraph->GetHistogram();
293 }
294 if (h)
295 {
296 TAxis *axe = h->GetXaxis();
297 // SetPoint deletes the histogram!
298 if (fUseEventNumber)
299 axe->SetTitle("Event Number");
300 else
301 {
302 axe->SetTitle("Time");
303 axe->SetLabelSize(0.033);
304 axe->SetTimeFormat("%H:%M %F1995-01-01 00:00:00 GMT");
305 axe->SetTimeDisplay(1);
306 }
307 }
308
309 if (TestBit(kIsLogy))
310 gPad->SetLogy();
311
312 // This is a workaround if the TGraph has only one point.
313 // Otherwise MStatusDisplay::Update hangs.
314// gPad->GetListOfPrimitives()->Remove(fGraph);
315// fGraph->Draw(fGraph->GetN()<2 ? "A" : str.Data());
316 //gPad->GetListOfPrimitives()->Add(fGraph, fGraph->GetN()<2 ? "A" : opt);
317 // AppendPad(str);
318 // fGraph->Draw(str);
319
320}
321
322// --------------------------------------------------------------------------
323//
324// This displays the TGraph like you expect it to be (eg. time on the axis)
325// It should also make sure, that the displayed time always is UTC and
326// not local time...
327//
328void MHVsTime::Draw(Option_t *opt)
329{
330 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(fGraph);
331 pad->SetBorderMode(0);
332 AppendPad(opt);
333
334 fGraph->Draw("A");
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.