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

Last change on this file since 4833 was 4828, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 7.3 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// eg.
32// MHVsTime hist("MHillas.fAlpha");
33// MHVsTime hist("MPointintPos.GetAbsErr");
34// MHVsTime hist("MPointintPos.GetAbsErr*kRad2Deg");
35//
36// To set a maximum number of data-points (eg. to display the last 20min
37// only) call SetMaxPts(200)
38//
39// SetMaxPts(-1) disables this feature.
40//
41/////////////////////////////////////////////////////////////////////////////
42#include "MHVsTime.h"
43
44#include <ctype.h> // tolower
45#include <fstream>
46
47#include <TPad.h>
48#include <TStyle.h>
49#include <TCanvas.h>
50
51#include <TGraph.h>
52
53#include "MLog.h"
54#include "MLogManip.h"
55
56#include "MTime.h"
57#include "MParList.h"
58#include "MDataChain.h"
59#include "MRawEvtHeader.h"
60
61ClassImp(MHVsTime);
62
63using namespace std;
64
65static const TString gsDefName = "MHVsTime";
66static const TString gsDefTitle = "Container for a graph vs time/evtnumber";
67
68// --------------------------------------------------------------------------
69//
70// Default constructor. For more informations about a valid rule
71// see MDataChain.
72//
73MHVsTime::MHVsTime(const char *rule)
74 : fGraph(NULL), fData(NULL), fScale(1), fMaxPts(-1), fUseEventNumber(0)
75{
76 fName = gsDefName;
77 fTitle = gsDefTitle;
78
79 if (!rule)
80 return;
81
82 fGraph = new TGraph;
83 fData = new MDataChain(rule);
84
85 fGraph->SetMarkerStyle(kFullDotMedium);
86}
87
88// --------------------------------------------------------------------------
89//
90// Deletes the histogram
91//
92MHVsTime::~MHVsTime()
93{
94 if (fGraph)
95 delete fGraph;
96
97 if (fData)
98 delete fData;
99}
100
101// --------------------------------------------------------------------------
102//
103// Return the data members used by the data chain to be used in
104// MTask::AddBranchToList
105//
106TString MHVsTime::GetDataMember() const
107{
108 return fData ? fData->GetDataMember() : (TString)"";
109}
110
111// --------------------------------------------------------------------------
112//
113// PreProcess the MDataChain. Create a new TGraph. Delete an old one if
114// already allocated.
115//
116Bool_t MHVsTime::SetupFill(const MParList *plist)
117{
118 // reset histogram (necessary if the same eventloop is run more than once)
119 //fGraph->Reset();
120
121 if (fData && !fData->PreProcess(plist))
122 return kFALSE;
123
124 if (fGraph)
125 {
126 delete fGraph;
127 fGraph = new TGraph;
128 }
129
130 TString title(fData ? GetRule() : (TString)"Histogram");
131 title += " vs ";
132 title += fUseEventNumber ? "Event Number" : "Time";
133
134 fGraph->SetNameTitle(fName, title);
135
136 fMean = 0;
137 fN = 0;
138
139 return kTRUE;
140}
141
142// --------------------------------------------------------------------------
143//
144// Set the name of the TGraph and the MHVsTime container
145//
146void MHVsTime::SetName(const char *name)
147{
148 fGraph->SetName(name);
149 MParContainer::SetName(name);
150}
151
152// --------------------------------------------------------------------------
153//
154// Set the title of the TGraph and the MHVsTime container
155//
156void MHVsTime::SetTitle(const char *title)
157{
158 fGraph->SetTitle(title);
159 MParContainer::SetTitle(title);
160}
161
162// --------------------------------------------------------------------------
163//
164// Set the next data point. If the graph exceeds fMaxPts remove the first
165//
166Bool_t MHVsTime::Fill(const MParContainer *par, const Stat_t w)
167{
168 Double_t t = 0;
169 if (fUseEventNumber)
170 {
171 const MRawEvtHeader *h = dynamic_cast<const MRawEvtHeader*>(par);
172 t = h ? h->GetDAQEvtNumber() : fGraph->GetN();
173 }
174 else
175 {
176 const MTime *tm = dynamic_cast<const MTime*>(par);
177 if (!tm)
178 {
179 *fLog << err << dbginf << "No MTime found..." << endl;
180 return kFALSE;
181 }
182 // If the time is not valid skip this entry
183 if (!*tm)
184 return kTRUE;
185 t = tm->GetAxisTime();
186 }
187
188 const Double_t v = fData->GetValue()*fScale;
189
190 if (fMaxPts>0 && fGraph->GetN()>fMaxPts)
191 fGraph->RemovePoint(0);
192
193 fMean += v;
194 fN++;
195
196 if (fN==fNumEvents)
197 {
198 fGraph->SetPoint(fGraph->GetN(), t, fMean/fN);
199 fMean = 0;
200 fN = 0;
201 }
202
203 return kTRUE;
204}
205
206// --------------------------------------------------------------------------
207//
208// This displays the TGraph like you expect it to be (eg. time on the axis)
209// It should also make sure, that the displayed time always is UTC and
210// not local time...
211//
212void MHVsTime::Draw(Option_t *opt)
213{
214 if (fGraph->GetN()==0)
215 return;
216
217 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(fGraph);
218 pad->SetBorderMode(0);
219
220 AppendPad("");
221
222 TString str(opt);
223
224 if (fUseEventNumber)
225 fGraph->GetHistogram()->SetXTitle("Event Number");
226 else
227 {
228 fGraph->GetHistogram()->SetXTitle("Time");
229 fGraph->GetHistogram()->GetXaxis()->SetTimeFormat("%H:%M:%S %F1995-01-01 00:00:00");
230 fGraph->GetHistogram()->GetXaxis()->SetTimeDisplay(1);
231 fGraph->GetHistogram()->GetXaxis()->SetLabelSize(0.033);
232 }
233 fGraph->GetHistogram()->SetYTitle(GetRule());
234
235 if (!str.Contains("A"))
236 str += "A";
237 if (!str.Contains("P"))
238 str += "P";
239
240 if (str.Contains("same", TString::kIgnoreCase))
241 {
242 str.ReplaceAll("same", "");
243 str.ReplaceAll("A", "");
244 }
245
246 fGraph->Draw(str);
247 if (fGraph->TestBit(kIsLogy))
248 pad->SetLogy();
249
250 pad->Modified();
251 pad->Update();
252}
253
254// --------------------------------------------------------------------------
255//
256// Used to rebuild a MHVsTime object of the same type (data members,
257// dimension, ...)
258//
259MParContainer *MHVsTime::New() const
260{
261 MHVsTime *h=new MHVsTime(fData ? (const char*)GetRule() : NULL);
262 h->SetScale(fScale);
263 if (fUseEventNumber)
264 h->SetUseEventNumber();
265 h->SetMaxPts(fMaxPts);
266 return h;
267}
268
269TString MHVsTime::GetRule() const
270{
271 return fData ? fData->GetRule() : (TString)"";
272}
273
274// --------------------------------------------------------------------------
275//
276// Returns the total number of bins in a histogram (excluding under- and
277// overflow bins)
278//
279Int_t MHVsTime::GetNbins() const
280{
281 return fGraph->GetN();
282}
283/*
284TH1 *MHVsTime::GetHist()
285{
286 return fGraph ? fGraph->GetHistogram() : 0;
287}
288
289const TH1 *MHVsTime::GetHist() const
290{
291 return fGraph ? fGraph->GetHistogram() : 0;
292}
293
294TH1 *MHVsTime::GetHistByName(const TString name)
295{
296 return GetHist();
297}
298*/
Note: See TracBrowser for help on using the repository browser.