/* ======================================================================== *\ ! ! * ! * This file is part of MARS, the MAGIC Analysis and Reconstruction ! * Software. It is distributed to you in the hope that it can be a useful ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes. ! * It is distributed WITHOUT ANY WARRANTY. ! * ! * Permission to use, copy, modify and distribute this software and its ! * documentation for any purpose is hereby granted without fee, ! * provided that the above copyright notice appear in all copies and ! * that both that copyright notice and this permission notice appear ! * in supporting documentation. It is provided "as is" without express ! * or implied warranty. ! * ! ! ! Author(s): Thomas Bretz, 11/2003 ! ! Copyright: MAGIC Software Development, 2000-2004 ! ! \* ======================================================================== */ ///////////////////////////////////////////////////////////////////////////// // // MHVsTime // // Use this class if you want to display any rule vs time (or event number) // // Axis titles // =========== // // 1) If no other title is given the rule for the y-axis is used. // 2) If the MH3 has a non-default title (MH3::SetTitle called) // this title is set as the histogram title. It can be used to overwrite // the axis titles. For more information see TH1::SetTitle, eg. // SetTitle("MyGraph;;Counts"); // The title for the x-axis is ignored and set automatically (MAKE SURE // YOU HAVE TWO SEMICOLON!) // // eg. // MHVsTime hist("MHillas.fAlpha"); // MHVsTime hist("MPointintPos.GetAbsErr"); // MHVsTime hist("MPointintPos.GetAbsErr*kRad2Deg"); // // To set a maximum number of data-points (eg. to display the last 20min // only) call SetMaxPts(200) // // SetMaxPts(-1) disables this feature. // ///////////////////////////////////////////////////////////////////////////// #include "MHVsTime.h" #include // tolower #include #include #include #include #include #include "MLog.h" #include "MLogManip.h" #include "MTime.h" #include "MParList.h" #include "MDataChain.h" #include "MRawEvtHeader.h" ClassImp(MHVsTime); using namespace std; const TString MHVsTime::gsDefName = "MHVsTime"; const TString MHVsTime::gsDefTitle = "Container for a graph vs time/evtnumber"; // -------------------------------------------------------------------------- // // Default constructor. For more informations about a valid rule // see MDataChain. // MHVsTime::MHVsTime(const char *rule, const char *error) : fGraph(NULL), fData(NULL), fError(0), fScale(1), fMaxPts(-1), fNumEvents(1), fUseEventNumber(0) { fName = gsDefName; fTitle = gsDefTitle; if (!rule) return; fData = new MDataChain(rule); if (error) fError = new MDataChain(error); fGraph = error ? new TGraphErrors : new TGraph; fGraph->SetPoint(0, 0, 0); // Dummy point! fGraph->SetEditable(); // Used as flag: First point? yes/no fGraph->SetMarkerStyle(kFullDotMedium); } // -------------------------------------------------------------------------- // // Deletes the histogram // MHVsTime::~MHVsTime() { if (fGraph) delete fGraph; if (fData) delete fData; } // -------------------------------------------------------------------------- // // Return the data members used by the data chain to be used in // MTask::AddBranchToList // TString MHVsTime::GetDataMember() const { return fData ? fData->GetDataMember() : (TString)""; } // -------------------------------------------------------------------------- // // PreProcess the MDataChain. Create a new TGraph. Delete an old one if // already allocated. // Bool_t MHVsTime::SetupFill(const MParList *plist) { if (!fGraph || !fData) { *fLog << err << "ERROR - MHVsTime cannot be used with its default constructor!" << endl; return kFALSE; } if (!fData->PreProcess(plist)) return kFALSE; if (fError && !fError->PreProcess(plist)) return kFALSE; fGraph->Set(1); fGraph->SetPoint(0, 0, 0); // Dummy point! fGraph->SetEditable(); // Used as flag: First point? yes/no TString title(fData ? GetRule() : (TString)"Graph"); title += " vs "; title += fUseEventNumber ? "Event Number" : "Time"; fGraph->SetNameTitle(fName, fTitle==gsDefTitle?title:fTitle); fMean = 0; fN = 0; return kTRUE; } // -------------------------------------------------------------------------- // // Set the name of the TGraph and the MHVsTime container // void MHVsTime::SetName(const char *name) { fGraph->SetName(name); MParContainer::SetName(name); } // -------------------------------------------------------------------------- // // Set the title of the TGraph and the MHVsTime container // void MHVsTime::SetTitle(const char *title) { fGraph->SetTitle(title); MParContainer::SetTitle(title); } // -------------------------------------------------------------------------- // // Set the next data point. If the graph exceeds fMaxPts remove the first // Bool_t MHVsTime::Fill(const MParContainer *par, const Stat_t w) { Double_t t = 0; if (fUseEventNumber) { const MRawEvtHeader *h = dynamic_cast(par); t = h ? h->GetDAQEvtNumber() : fGraph->GetN(); } else { const MTime *tm = dynamic_cast(par); if (!tm) { *fLog << err << dbginf << "No MTime found..." << endl; return kFALSE; } // If the time is not valid skip this entry if (!*tm) return kTRUE; // Do not fill events with equal time if (*tm==fLast || *tm==MTime()) return kTRUE; fLast = *tm; t = tm->GetAxisTime(); } const Double_t v = fData->GetValue(); const Double_t e = fError ? fError->GetValue() : 0; //*fLog << all << "ADD " << v << " " << e << endl; fMean += v; fMeanErr += e; fN++; if (fN==fNumEvents) { if (fMaxPts>0 && fGraph->GetN()>fMaxPts || fGraph->IsEditable()) { fGraph->RemovePoint(0); fGraph->SetEditable(kFALSE); } fGraph->SetPoint(fGraph->GetN(), t, fMean/fN*fScale); if (fError) static_cast(fGraph)->SetPointError(fGraph->GetN()-1, 0, fMeanErr/fN*fScale); fMean = 0; fMeanErr = 0; fN = 0; } return kTRUE; } void MHVsTime::Paint(Option_t *opt) { // SetPoint deletes the histogram! if (fUseEventNumber) fGraph->GetHistogram()->SetXTitle("Event Number"); else { fGraph->GetHistogram()->SetXTitle("Time"); fGraph->GetHistogram()->GetXaxis()->SetLabelSize(0.033); fGraph->GetHistogram()->GetXaxis()->SetTimeFormat("%H:%M:%S %F1995-01-01 00:00:00 GMT"); fGraph->GetHistogram()->GetXaxis()->SetTimeDisplay(1); } if (TestBit(kIsLogy)) gPad->SetLogy(); // This is a workaround if the TGraph has only one point. // Otherwise MStatusDisplay::Update hangs. gPad->GetListOfPrimitives()->Remove(fGraph); gPad->GetListOfPrimitives()->Add(fGraph, fGraph->GetN()<2 ? "A" : opt); } // -------------------------------------------------------------------------- // // This displays the TGraph like you expect it to be (eg. time on the axis) // It should also make sure, that the displayed time always is UTC and // not local time... // void MHVsTime::Draw(Option_t *opt) { if (!fGraph) return; if (fGraph->GetN()==0) return; TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(fGraph); pad->SetBorderMode(0); TString str(opt); if (!str.Contains("A")) str += "A"; if (!str.Contains("P")) str += "P"; if (str.Contains("same", TString::kIgnoreCase)) { str.ReplaceAll("same", ""); str.ReplaceAll("A", ""); } AppendPad(str); fGraph->Draw(str); } // -------------------------------------------------------------------------- // // Used to rebuild a MHVsTime object of the same type (data members, // dimension, ...) // MParContainer *MHVsTime::New() const { MHVsTime *h=new MHVsTime(fData ? (const char*)GetRule() : NULL); h->SetScale(fScale); if (fUseEventNumber) h->SetUseEventNumber(); h->SetMaxPts(fMaxPts); return h; } TString MHVsTime::GetRule() const { return fData ? fData->GetRule() : (TString)""; } // -------------------------------------------------------------------------- // // Returns the total number of bins in a histogram (excluding under- and // overflow bins) // Int_t MHVsTime::GetNbins() const { return fGraph->GetN(); } /* TH1 *MHVsTime::GetHist() { return fGraph ? fGraph->GetHistogram() : 0; } const TH1 *MHVsTime::GetHist() const { return fGraph ? fGraph->GetHistogram() : 0; } TH1 *MHVsTime::GetHistByName(const TString name) { return GetHist(); } */