//*-- Modified : v0r51 07/11/00 by Damir Buskulic//*-- Modified : v0r50 28/09/00 by Damir Buskulic//*-- Modified : v0r49 21/06/00 by Damir Buskulic//*-- Modified : v0r48 30/05/00 by Damir Buskulic//*-- Modified : v0r47 14/04/00 by Damir Buskulic//*-- Modified : v0r45 11/11/99 by Damir Buskulic//*-- Modified : dev version 10/12/98 by Damir Buskulic//*-- Author : Damir Buskulic 10/9/98//////////////////////////////////////////////////////////////////////////// //// VSPlot //// //// class describing a plot (representation) of a series //// ////////////////////////////////////////////////////////////////////////////
#include "VManagerFrameL.h"
#include "VStyle.h"
#include "VSPlot.h"
#include "VGPSTime.h"
#include "TPad.h"
#include "TH1.h"
#include "TPaveStats.h"
#include "TLatex.h"
#include "TAxis.h"
const Int_t kNoStats = BIT(9);
const Int_t kAxisRange = BIT(11);
ClassImp(VSPlot)
//______________________________________________________________________________VSPlot::VSPlot() : TH1D()
{}
//______________________________________________________________________________VSPlot::VSPlot(Text_t* name, Text_t* title, Int_t nbinsx, Axis_t xlow, Axis_t xup, VManager* mng)
:TH1D(name,title,nbinsx,xlow,xup)
{
mManager = mng;
BuildVSP();
}
//______________________________________________________________________________VSPlot::VSPlot(Text_t* name, Text_t* title, Int_t nbinsx, Axis_t xlow, Axis_t xup)
:TH1D(name,title,nbinsx,xlow,xup)
{
mManager = 0;
BuildVSP();
}
//______________________________________________________________________________VSPlot::VSPlot(Text_t* name, Text_t* title, Int_t nbinsx, Axis_t* xbins)
:TH1D(name,title,nbinsx,xbins)
{
mManager = 0;
BuildVSP();
}
//______________________________________________________________________________VSPlot::VSPlot(Text_t* name, Text_t* title, Text_t* asciifilename)
:TH1D(name,title,1,0,1)
{
// Special constructor. Takes data from an ascii file and includes this data// in a VSPlot. The ascii file should have the format// Time1 value1// Time2 value2// etc...double time,value;
FILE* asciifile;
int i;
int arraysize=1000;
TArrayD timearr(arraysize);
TArrayD valuearr(arraysize);
// Builds standard parametersmManager = 0;
BuildVSP();
// The lowest time and step are determined from the ascii fileDouble_t loweststep = 100000000000.;
Double_t oldstarttime = 0.;
Double_t timax = 0.;
Double_t timin = 100000000000.;
// opens the ascii file and reads data from it
asciifile = fopen(asciifilename,"r");
if (!asciifile) return;
int nentries = 0;
int okreadfile;
okreadfile = fscanf(asciifile,"%lf %lf %*[^n]",&time,&value);
while(okreadfile>0) {
timearr[nentries] = time;
valuearr[nentries] = value;
// Sets the min and max time
if (time>timax) timax = time;
if (time<timin) timin = time;
// See if the time step between values is the lowest one. There must be at least// some order in the frames to make it work.double currentstep = time - oldstarttime;
if (currentstep < loweststep) loweststep = currentstep;
oldstarttime = time;
nentries++;
if (nentries>=arraysize) {
arraysize *= 2; // double the array size
timearr.Set(arraysize);
valuearr.Set(arraysize); // This is non destructive
}
okreadfile = fscanf(asciifile,"%lf %lf %*[^n]",&time,&value);
}
// Sets the start time to be the min timemGPSStart.SetTimeD(timin);
// Finds the number of pointsint nsteps = (int)((timax-timin)/loweststep)+1;
// Resets the plot part and sets the limits and number of points
Reset();
SetBins(nsteps,0,timax-timin+loweststep);
// Sets the values
for (i = 0; i<nentries;i++) {
int istep = (int)((timearr[i] - timin + loweststep/2)/loweststep);
SetBinContent(istep+1, valuearr[i]);
}
if (fclose(asciifile)) Warning("VSPlot","Error in closing asciifile");
}
//______________________________________________________________________________VSPlot::VSPlot(VSPlot &vsplot)
{
((VSPlot&)vsplot).Copy(*this);
vsplot.AdoptManager(0);
}
//______________________________________________________________________________voidVSPlot::Copy(TObject& newvsp)
{
TH1D::Copy((VSPlot&)newvsp);
((VSPlot&)newvsp).mGPSStart = mGPSStart;
((VSPlot&)newvsp).mStep = mStep;
((VSPlot&)newvsp).mLeapS = mLeapS;
((VSPlot&)newvsp).mLocalTime = mLocalTime;
((VSPlot&)newvsp).fStartTimeDisplay = fStartTimeDisplay;
((VSPlot&)newvsp).mStartTimeXOffset = mStartTimeXOffset;
((VSPlot&)newvsp).mStartTimeYOffset = mStartTimeYOffset;
((VSPlot&)newvsp).mUnitX = mUnitX.Data();
((VSPlot&)newvsp).mUnitY = mUnitY.Data();
}
//______________________________________________________________________________voidVSPlot::BuildVSP()
{
// Builds the VSPlot and deletes VSPlots having the same name and maintained// by the same manager
if (mManager) {
// VSPlot *vsold = (VSPlot*) mManager->GetListOfPlots()->FindObject(GetName());// if (vsold!=0) {// if (gVEGAWarnings) {// Warning("BuildVSP","Replacing with deletion existing series plot : %s",GetName());// }// delete vsold;// }
SetBit(kCanDelete);
mManager->GetListOfPlots()->Add(this);
}
// Sets the options
if (gVStyle) {
fStartTimeDisplay = gVStyle->GetStartTimeDisplay();
mStartTimeXOffset = gVStyle->GetStartTimeXOffset();
mStartTimeYOffset = gVStyle->GetStartTimeYOffset();
fXaxis.SetTimeDisplay(gVStyle->GetTimeOnXAxis());
fXaxis.SetTimeFormat(VGPSTime::ExpandFormat(gVStyle->GetTimeFormat()));
} else {
fStartTimeDisplay = 1;
mStartTimeXOffset = 0.02;
mStartTimeYOffset = 0.02;
}
mLocalTime = 0;
mStep = 0;
mLeapS = 0;
SetGPSStart(0);
}
//______________________________________________________________________________VSPlot::~VSPlot()
{
// VSPlot destructor. Also removes the entry concerning this plot// in the original serieAbandonManager();
}
//______________________________________________________________________________voidVSPlot::AbandonManager()
{
// Removes this plot from the maintenance of the VEGA manager
if (mManager) {
mManager->RemovePlot(GetName());
mManager=0;
}
}
//______________________________________________________________________________voidVSPlot::AdoptManager(VManager* mgr)
{
// Adopt the manager as the maintainer of this plot
if (mManager) AbandonManager();
mManager=mgr;
if (mManager) mManager->GetListOfPlots()->Add(this);
}
//______________________________________________________________________________Stat_tVSPlot::GetPlotRMS()
{
// Return the Root Mean Square of the values of this seriesStat_t x, rms2, stats[10];
GetStats(stats);
if (stats[0] == 0) return 0;
Double_t N = GetXaxis()->GetLast() - GetXaxis()->GetFirst() +1;
x = stats[0]/N;
rms2 = TMath::Abs(stats[1]/N - x*x);
return TMath::Sqrt(rms2);
}
//______________________________________________________________________________Stat_tVSPlot::GetAverageLevel()
{
// Return the average level of the values of this seriesStat_t stats[10];
GetStats(stats);
if (stats[0] == 0) return 0;
Double_t N = GetXaxis()->GetLast() - GetXaxis()->GetFirst() +1;
return stats[0]/N;
}
//______________________________________________________________________________voidVSPlot::GetStats(Stat_t *stats)
{
// fill the array stats from the contents of this plot/series // The array stats must be correctly dimensionned in the calling program. // stats[0] = sumw // stats[1] = sumw2 // stats[2] = sumwx // stats[3] = sumwx2 // stats[4] = sumwy if 2-d or 3-d // stats[5] = sumwy2 // stats[6] = sumwz if 3-d // stats[7] = sumwz2 // // If no axis-subrange is specified (via TAxis::SetRange), the array stats // is simply a copy of the statistics quantities computed at filling time. // If a sub-range is specified, the function recomputes these quantities // from the bin contents in the current axis range. // Loop on bins (possibly including underflows/overflows)Int_t bin, binx;
Stat_t w;
Double_t x;
if (fDimension == 1) {
if (fTsumw == 0 || fXaxis.TestBit(kAxisRange)) {
for (bin=0;bin<4;bin++) stats[bin] = 0;
for (binx=fXaxis.GetFirst();binx<=fXaxis.GetLast();binx++) {
x = fXaxis.GetBinCenter(binx);
w = GetBinContent(binx);
stats[0] += w;
stats[1] += w*w;
stats[2] += w*x;
stats[3] += w*x*x;
}
} else {
stats[0] = fTsumw;
stats[1] = fTsumw2;
stats[2] = fTsumwx;
stats[3] = fTsumwx2;
}
}
}
// ______________________________________________________________________________Int_tVSPlot::ShiftAppend(Int_t nval, Double_t* values)
{
// Shift the contents of the VSPlot.// "nval" is the length of the "values" input array, nval values are removed// from the beginning of the VSPlot and the new nval values are appended// to the end. If nval is bigger than the number of points in the plot,// nothing is done.// The start GPS time is adjusted accordingly.// Returns 1 in case of success, 0 otherwiseInt_t i;
Int_t nsamples = GetNbinsX();
if (nval>nsamples) return 0;
// Shift the part of the plot that is not deleted
for (i=nval+1;i<=nsamples;i++) {
// fArray has nsamples+2 elements and fArray[0] is the number of underflows// in an histogram (the base class of VSPlot) -> to be changed if change// the base class
fArray[i-nval] = fArray[i];
}
// Add the new valuesInt_t nshift = nsamples-nval+1;
for (i=0;i<nval;i++) {
fArray[i+nshift] = values[i];
}
// reset the statistics
fTsumw = 0;
// adjust the GPS time.// printf("mGPSstart = %d, %dn",mGPSStart.GetSec(), mGPSStart.GetNsec());mGPSStart += (double)nval*mStep;
return 1;
}
//______________________________________________________________________________TH1* VSPlot::DrawCopy(Option_t *option)
{
TString opt = option;
opt.ToLower();
if (gPad && !opt.Contains("same")) gPad->Clear();
VSPlot *newtvs = (VSPlot*)Clone();
newtvs->SetDirectory(0);
newtvs->SetBit(kCanDelete);
newtvs->AppendPad(opt.Data());
return newtvs;
}
// ______________________________________________________________________________Int_tVSPlot::DistancetoPrimitive(Int_t px, Int_t py)
{
// -*-*-*-*-*-*-*-*-*Compute distance from point px,py to a line*-*-*-*-*-*// ===========================================// Compute the closest distance of approach from point px,py to elements// of a VSPlot. Special treatment due to x axis being time axis.// in the case of "sameti" option, have to shift the plots so they// coincide in time.// The distance is computed in pixels units.Int_t dist;
if (!gPad) return 9999;
TString opt = GetDrawOption();
Double_t oldstart = GetXaxis()->GetXmin();
Double_t oldend = GetXaxis()->GetXmax();
if (opt.Contains("sameti")) {
// find the base plot
TObject* obj;
TListIter next(gPad->GetListOfPrimitives());
while ((obj = next()))
if (!strcmp(obj->IsA()->GetName(),"VSPlot")) break;
if (obj) {
VSPlot* baseplot = (VSPlot*)obj;
Double_t startx = GetGPSStart().GetTimeD() - baseplot->GetGPSStart().GetTimeD();
Double_t endx = startx + GetXaxis()->GetXmax() - GetXaxis()->GetXmin();
GetXaxis()->SetLimits(startx,endx);
}
}
dist = TH1::DistancetoPrimitive(px,py);
GetXaxis()->SetLimits(oldstart,oldend);
return dist;
}
// ______________________________________________________________________________voidVSPlot::Paint(Option_t *option)
{
//*-*-*-*-*-*-*-*-*Control routine to paint the VSPlots*-*-*-*-*-*-*//*-* ===============================================//// This function is automatically called by TCanvas::Update.// (see TH1::Draw for the list of options)
TString t;
Double_t tosave=0;
Double_t timeoffsave=0;
TString opt=option;
opt.ToLower();
if (!fPainter) fPainter = TVirtualHistPainter::HistPainter(this);
if (gPad->GetLogx()) {
// bmsave = gPad->GetBottomMargin();// gPad->SetBottomMargin(bmsave+0.02);
tosave = GetTitleOffset();
SetTitleOffset(tosave+0.35);
}
#if ROOT_VERSION_CODE <= ROOT_VERSION(3,5,2)
Double_t utctimeT;
// If axis painting with time display, uses the frame start time as time offset
if (GetXaxis()->GetTimeDisplay()) {
timeoffsave = gVStyle->GetTimeOffset();
// Transform gps time into utc time
utctimeT = (mGPSStart.GetSec() - mLeapS + 19 + 315964800);
if (strstr(gVStyle->GetStartTimeFormat(),"local")) {
utctimeT += mLocalTime;
}
gVStyle->SetTimeOffset(utctimeT);
}
#endif
if (fPainter) {
if (gPad->GetListOfPrimitives() && opt.Contains("sameti")) {
// superimpose, use true time
opt.ReplaceAll("sameti","same");
// find the base plot
TObject* obj;
TListIter next(gPad->GetListOfPrimitives());
while ((obj = next()))
if (!strcmp(obj->IsA()->GetName(),"VSPlot")) break;
if (obj) {
VSPlot* baseplot = (VSPlot*)obj;
Double_t oldstartx, oldendx;
Double_t startx = GetGPSStart().GetTimeD() - baseplot->GetGPSStart().GetTimeD();
Double_t endx = startx + GetXaxis()->GetXmax() - GetXaxis()->GetXmin();
oldstartx = GetXaxis()->GetXmin();
oldendx = GetXaxis()->GetXmax();
GetXaxis()->SetLimits(startx,endx);
fPainter->Paint(opt.Data());
GetXaxis()->SetLimits(oldstartx,oldendx);
}
} else if (TestBit(kNoStats)) {
fPainter->Paint(option);
} else {
SetBit(kNoStats);
fPainter->Paint(option);
TString opt = option;
opt.ToLower();
if (!opt.Contains("same")) {
// Painting the statistics box
TF1 *fit = (TF1*)fFunctions->First();
PaintStatVS(gStyle->GetOptStat(),fit);
// Painting the start time below axis
if (fStartTimeDisplay) {
TLatex* textstime = new TLatex();
textstime->SetNDC(kTRUE);
textstime->SetTextColor(fXaxis.GetLabelColor());
textstime->SetTextFont(fXaxis.GetLabelFont());
textstime->SetTextSize(fXaxis.GetLabelSize());
textstime->SetX(mStartTimeXOffset);
textstime->SetY(mStartTimeYOffset);
textstime->SetTextAlign(12);
if (strstr(gVStyle->GetStartTimeFormat(),"local") && mLocalTime>0 && mLocalTime<85300) {
t = "local (this site) T0: ";
} else {
t = "GPS T0: ";
}
t = t + VGPSTime::ExpressTime(mGPSStart.GetTimeD(), gVStyle->GetStartTimeFormat(), mLeapS,mLocalTime);
textstime->SetTitle(t.Data());
textstime->Paint();
delete textstime;
}
}
ResetBit(kNoStats);
}
}
if (gPad->GetLogx()) {
// gPad->SetBottomMargin(bmsave);
SetTitleOffset(tosave);
}
#if ROOT_VERSION_CODE <= ROOT_VERSION(3,5,2)
// Restore ref. time in gVStyle
if (GetXaxis()->GetTimeDisplay()) gVStyle->SetTimeOffset(timeoffsave);
#endif
}
// ______________________________________________________________________________voidVSPlot::PaintStatVS(Int_t dostat, TF1 *fit)
{
// *-*-*-*-*-*-*-*-*-*-*-*Draw the statistics box*-*-*-*-*-*-*-*-*-*-*-*-*-*-*// *-* =======================// The type of information printed in the VSPlot statistics box// can be selected via gVStyle->SetOptStat(mode).// The parameter mode can be = rasn (default = 1111)// n = 1; name of histogram is printed// s = 1; start time of the vector// a = 1; average value printed// r = 1; rms printed// Example: gStyle->SetOptStat(11);// print only name of vsplot and start time.// // The type of information about fit parameters printed in the histogram// statistics box can be selected via the parameter mode.// The parameter mode can be = pcev (default = 0111)// v = 1; print name/values of parameters// e = 1; print errors (if e=1, v must be 1)// c = 1; print Chisquare/Number of degress of freedom// p = 1; print Probability// Example: gStyle->SetOptFit(1011);// print fit probability, parameter names/values and errors.//
static char t[64];
Int_t dofit;
TPaveStats *stats = 0;
TIter next(GetListOfFunctions());
TObject *obj;
while ((obj = next())) {
if (obj->InheritsFrom(TPaveStats::Class())) {
stats = (TPaveStats*)obj;
break;
}
}
if (stats) {
dofit = stats->GetOptFit();
dostat = stats->GetOptStat();
} else {
dofit = gStyle->GetOptFit();
}
if (!dofit) fit = 0;
if (dofit == 1) dofit = 111;
if (dostat == 1) dostat = 1111;
Int_t print_name = dostat%10;
Int_t print_start = (dostat/10)%10;
Int_t print_average = (dostat/100)%10;
Int_t print_rms = (dostat/1000)%10;
Int_t nlines = print_name + print_start + print_average + print_rms;
Int_t print_fval = dofit%10;
Int_t print_ferrors = (dofit/10)%10;
Int_t print_fchi2 = (dofit/100)%10;
Int_t print_fprob = (dofit/1000)%10;
Int_t nlinesf = print_fval + print_fchi2 + print_fprob;
if (fit) nlinesf += fit->GetNpar();
// *-*- Pavetext with statisticsBool_t done = kFALSE;
if (!dostat && !fit) {
if (stats) { delete stats; GetListOfFunctions()->Remove(stats); }
return;
}
Double_t statw = 1.4*gStyle->GetStatW();
if (fit) statw = 1.8*gStyle->GetStatW();
Double_t stath = 0.25*(nlines+nlinesf)*gStyle->GetStatH();
if (stats) {
stats->Clear();
done = kTRUE;
} else {
stats = new TPaveStats(
gStyle->GetStatX()-statw,
gStyle->GetStatY()-stath,
gStyle->GetStatX(),
gStyle->GetStatY(),"brNDC");
stats->SetParent(GetListOfFunctions());
stats->SetOptFit(dofit);
stats->SetOptStat(dostat);
stats->SetFillColor(gStyle->GetStatColor());
stats->SetFillStyle(gStyle->GetStatStyle());
stats->SetBorderSize(gStyle->GetStatBorderSize());
stats->SetTextFont(gStyle->GetStatFont());
stats->SetName("stats");
stats->SetTextColor(1);
stats->SetTextAlign(12);
stats->SetBit(kCanDelete);
}
if (print_name)
stats->AddText( GetName() );
if (print_start) {
sprintf(t,"Start = %s",VGPSTime::ExpressTime(mGPSStart.GetTimeD(),gVStyle->GetTimeFormat(), mLeapS,mLocalTime));
stats->AddText(t);
}
char textstats[50];
if (print_average) {
sprintf(textstats,"Avg Lvl = %s%s","%",stats->GetStatFormat());
sprintf(t,textstats,GetAverageLevel());
// sprintf(t," = %6.4g",);
stats->AddText(t);
}
if (print_rms) {
sprintf(textstats,"RMS = %s%s","%",stats->GetStatFormat());
sprintf(t,textstats,GetAverageLevel());
sprintf(t,"RMS = %6.4g",GetPlotRMS());
stats->AddText(t);
}
// Draw Fit parameters
if (fit) {
Int_t ndf = fit->GetNDF();
sprintf(textstats,"#chi^{2} / ndf = %s%s / %d","%",stats->GetFitFormat(),ndf);
sprintf(t,textstats,(Float_t)fit->GetChisquare());
if (print_fchi2) stats->AddText(t);
if (print_fprob) {
sprintf(textstats,"Prob = %s%s","%",stats->GetFitFormat());
sprintf(t,textstats,(Float_t)TMath::Prob(fit->GetChisquare(),ndf));
stats->AddText(t);
}
if (print_fval || print_ferrors) {
for (Int_t ipar=0;ipar<fit->GetNpar();ipar++) {
if (print_ferrors) {
sprintf(textstats,"%-8s = %s%s #pm %s%s ",fit->GetParName(ipar),"%",stats->GetFitFormat(),"%",stats->GetFitFormat());
sprintf(t,textstats,(Float_t)fit->GetParameter(ipar)
,(Float_t)fit->GetParError(ipar));
} else {
sprintf(textstats,"%-8s = %s%s ",fit->GetParName(ipar),"%",stats->GetFitFormat());
sprintf(t,textstats,(Float_t)fit->GetParameter(ipar));
}
t[63] = 0;
stats->AddText(t);
}
}
}
if (!done) {GetListOfFunctions()->Add(stats);}
stats->Paint();
}
//______________________________________________________________________________char *VSPlot::GetObjectInfo(Int_t px, Int_t py) const
{
// Redefines TObject::GetObjectInfo.// Displays the vsplot info (time, bin number, contents)// corresponding to cursor position px,py//
if (!gPad) return (char*)"";
static char info[64];
Double_t x = gPad->PadtoX(gPad->AbsPixeltoX(px));
Double_t y = gPad->PadtoY(gPad->AbsPixeltoY(py));
Double_t x1 = gPad->PadtoX(gPad->AbsPixeltoX(px+1));
const char *drawOption = GetDrawOption();
Double_t xmin, xmax, uxmin,uxmax;
Double_t ymin, ymax, uymin,uymax;
Int_t binx,biny,binmin,binx1;
Double_t diffstart;
// If logarithmic scale, probably frequency axis, return standard GetInfo
if (gPad->GetLogx()) return TH1::GetObjectInfo(px,py);
// Get the start point for this plot if it was drawn with the "sameti" option
diffstart = 0;
TString opt = GetDrawOption();
if (opt.Contains("sameti")) {
// First find the base plot
TObject* obj;
TListIter next(gPad->GetListOfPrimitives());
while ((obj = next()))
if (!strcmp(obj->IsA()->GetName(),"VSPlot")) break;
if (obj) {
VSPlot* baseplot = (VSPlot*)obj;
VGPSTime gpsstartobj = mGPSStart;
diffstart = gpsstartobj.GetTimeD() - baseplot->GetGPSStart().GetTimeD();
}
}
if (gPad->IsVertical()) {
binx = fXaxis.FindFixBin(x);
binmin = fXaxis.GetFirst();
binx1 = fXaxis.FindFixBin(x1);
// special case if more than 1 bin in x per pixel
if (binx1-binx>1 && GetDimension() == 1) {
Double_t binval=GetBinContent(binx);
Int_t binnear=binx;
for (Int_t ibin=binx+1; ibin<binx1; ibin++) {
Double_t binvaltmp = GetBinContent(ibin);
if (TMath::Abs(y-binvaltmp) < TMath::Abs(y-binval)) {
binval=binvaltmp;
binnear=ibin;
}
}
binx = binnear;
}
} else {
x1 = gPad->PadtoY(gPad->AbsPixeltoY(py+1));
binx = fXaxis.FindFixBin(y);
binmin = fXaxis.GetFirst();
binx1 = fXaxis.FindFixBin(x1);
// special case if more than 1 bin in x per pixel
if (binx1-binx>1 && GetDimension() == 1) {
Double_t binval=GetBinContent(binx);
Int_t binnear=binx;
for (Int_t ibin=binx+1; ibin<binx1; ibin++) {
Double_t binvaltmp = GetBinContent(ibin);
if (TMath::Abs(x-binvaltmp) < TMath::Abs(x-binval)) {
binval=binvaltmp;
binnear=ibin;
}
}
binx = binnear;
}
}
VGPSTime gpsstartobj = mGPSStart;
Double_t gpsstart = gpsstartobj.GetTimeD();
sprintf(info,"GPS=%.4f, D=(%s), binx=%d, value=%g",x-diffstart+gpsstart,
VGPSTime::ExpressTime(floor(gpsstart+x-diffstart),"local",mLeapS,mLocalTime),binx,GetBinContent(binx));
return info;
}
//______________________________________________________________________________voidVSPlot::SetGPSStart(Double_t gpsstart)
{
// Set GPS time of the start of the X axistime_t utctime;
mGPSStart = gpsstart;
#if ROOT_VERSION_CODE > ROOT_VERSION(3,5,2)
if (gpsstart>0) {
// Transform gps time into utc time
utctime = (mGPSStart.GetSec() - mLeapS + 19 + 315964800);
if (strstr(gVStyle->GetStartTimeFormat(),"local") && mLocalTime>0 && mLocalTime<85300) {
utctime += mLocalTime;
}
struct tm* tg;
tg = gmtime(&utctime);
tg->tm_isdst = -1;
GetXaxis()->SetTimeOffset(mktime(tg));
}
#endif
}
//______________________________________________________________________________voidVSPlot::SetGPSStart(UInt_t gpsstartS, UInt_t gpsstartNS)
{
// Set GPS time of the start of the X axistime_t utctime;
mGPSStart.SetSec(gpsstartS); mGPSStart.SetNsec(gpsstartNS);
#if ROOT_VERSION_CODE > ROOT_VERSION(3,5,2)
if (gpsstartS>0) {
// Transform gps time into utc time
utctime = (mGPSStart.GetSec() - mLeapS + 19 + 315964800);
if (strstr(gVStyle->GetStartTimeFormat(),"local") && mLocalTime>0 && mLocalTime<85300) {
utctime += mLocalTime;
}
struct tm* tg;
tg = gmtime(&utctime);
tg->tm_isdst = -1;
GetXaxis()->SetTimeOffset(mktime(tg));
}
#endif
}
- ROOT page - VEGA page - Class index - Top of the page This page has been automatically generated. If you have any comments or suggestions
about the page layout send a mail to
, or
contact
with any questions or problems regarding ROOT or VEGA.