//*-- Modified : v0r51 07/11/00 by D. Buskulic
//*-- Modified : v0r50 29/09/00 by D. Buskulic
//*-- Modified : v0r48 30/05/00 by D. Buskulic
//*-- Modified : v0r47 17/04/00 by D. Buskulic
//*-- Modified : v0r46 13/03/00 by D. Buskulic
//*-- Modified : v0.45 11/11/99 by D. Buskulic
//*-- Modified : v0.42 02/09/99 by D. Buskulic
//*-- Author :    Damir Buskulic   04/03/99

//////////////////////////////////////////////////////////////////////////
//                                                                      //
// VManagerFrameL                                                       //
//                                                                      //
// Frames manager class using the Framelib.                             //
//                                                                      //
//////////////////////////////////////////////////////////////////////////

#include <iostream.h>
#include "VManagerFrameL.h"
#include "VFrUtil.h"
#include "VStyle.h"
#include "TPad.h"
#include "TROOT.h"
#include "TSystem.h"

ClassImp(VManagerFrameL)

//______________________________________________________________________________
VManagerFrameL::VManagerFrameL() : VManager()
{
// ---- VManagerFrameL default constructor ----
//

   mPlots = new TList(this);
   mPlots2 = new TList(this);
   mHistos = new TList(this);

   gVM = this;
}

//______________________________________________________________________________
VManagerFrameL::~VManagerFrameL()
{
// ---- VManagerFrameL destructor ----
//

   mPlots->ForEach(VSPlot,AbandonManager)();
   delete mPlots;
   mPlots2->ForEach(VSPlot2,AbandonManager)();
   delete mPlots2;
   
   delete mHistos;
   
   if (gDebug)
      cerr <<"VManager dtor called " << endl;
}

//______________________________________________________________________________
 void VManagerFrameL::Draw(FrVect* vect, Double_t offset, Double_t dur, Option_t *option)
{
//         Create VSPlot and Draw it starting from an FrVect
//         =================================================
// The serie is plot from "offset" for a length "dur". The default values
// mean that the entire series is plotted.
// For the moment, options are the same as for histograms

   TString  opt;
   opt = option;
   opt.ToLower();
   
   UInt_t i,ix,iy,ibegin,iend,ndat;
   Double_t vsoffset,vsdur,vspbegin, vspend, bining;
   Double_t vspybegin,vspyend;
   
   VSPlot* vsp;
   VSPlot2* vsp2;

// Protect against bad values
   if (!vect) return;
   char* vectname = vect->name;

// Get info for number of samples, offset and length   
   ndat = vect->nx[0];
   if (!ndat) {
      Error("Draw","No data in this series, cannot draw");
      return ;
   }
   vsoffset = vect->startX[0]; // this is the initial offset of the vector
                 // with respect to begining of the frame
   vsdur = VGetLength(vect);
   if (vsdur==0) {
      vsdur = (float)ndat;
      if (gVEGAWarnings)
         Warning("Draw","Null length for serie. Using length=ndata");
   }

// Define begin and end of the plot. Correct it if necessary
   vspbegin = vsoffset;
   if (offset>-9e19) {vspbegin += offset;}
   if (dur<9e19) {vspend = vspbegin+dur;} else {vspend=vspbegin+vsdur;}
   if (vspbegin<vsoffset) {vspbegin = vsoffset;}
   if (vspend>vsoffset+vsdur) {vspend = vsoffset+vsdur;}
   
   bining = vsdur/ndat;

// Compute the numbers of first and last bin to be displayed
   ibegin = int((vspbegin-vsoffset)/bining);
   iend   = int((vspend-vsoffset)/bining);

// Correct begin and end values of the series for display
   vspbegin = vsoffset+bining*ibegin;
   vspend = vsoffset+bining*iend;

// Sets the offset time of the VSPlot
   Double_t vspoff = 0;
   if (vspbegin !=0) {
      if (TMath::Abs((vspend-vspbegin)/vspbegin) < 1e-4) {
         modf(vspbegin*1e-3,&vspoff);
         vspoff *= 1.e3;
         vspbegin = vspbegin - vspoff;
         vspend = vspend - vspoff;
      }
   }

// 
//  Treat the two cases : 1D (create a VSPlot) or 2D (create a VSPlot2)
//
   if (vect->nDim==1) {
// 1 dim case
      vsp = new VSPlot(vectname,vectname,iend-ibegin,0,vspend - vspbegin,this);
   // These settings should be set BEFORE the GPS start because setting
   // GPS start uses them.
//      vsp->SetOffset(vspoff);
      vsp->SetStep(vect->dx[0]);
      vsp->SetLocalTime(vect->localTime);
      vsp->SetLeapS(vect->ULeapS);
#if FR_VERS<6000
      vsp->SetGPSStart(vect->GTimeS + 10e-9*vect->GTimeN + vspbegin);
#else
      vsp->SetGPSStart(vect->GTime + vspbegin);
#endif
      if (vect->unitX) {
         if (strcmp(vect->unitX[0],"NONE"))vsp->GetXaxis()->SetTitle(vect->unitX[0]);
      }
      if (vect->unitY) {
         if (strcmp(vect->unitY,"NONE")) vsp->GetYaxis()->SetTitle(vect->unitY);
      }
   
      unsigned short typevect = vect->type;
      if (typevect == FR_VECT_C) {
         unsigned char* tdata= (unsigned char*)vect->data;
         for (i=0;i<iend-ibegin;i++) {vsp->SetBinContent(i+1,tdata[i+ibegin]);}
      } else if (typevect == FR_VECT_2U) {
         unsigned short* tdataUS= vect->dataUS;
         for (i=0;i<iend-ibegin;i++) {vsp->SetBinContent(i+1,tdataUS[i+ibegin]);}
      } else if (typevect == FR_VECT_2S) {
         short* tdataS= vect->dataS;
         for (i=0;i<iend-ibegin;i++) {vsp->SetBinContent(i+1,tdataS[i+ibegin]);}
      } else if (typevect == FR_VECT_4U) {
         unsigned int* tdataUI= vect->dataUI;
         for (i=0;i<iend-ibegin;i++) {vsp->SetBinContent(i+1,tdataUI[i+ibegin]);}
      } else if (typevect == FR_VECT_4S) {
         int* tdataI= vect->dataI;
         for (i=0;i<iend-ibegin;i++) {vsp->SetBinContent(i+1,tdataI[i+ibegin]);}
      } else if (typevect == FR_VECT_4R) {
         float* tdataF= vect->dataF;
         for (i=0;i<iend-ibegin;i++) {vsp->SetBinContent(i+1,tdataF[i+ibegin]);}
      } else if (typevect == FR_VECT_8R) {
         double* tdataD= vect->dataD;
         for (i=0;i<iend-ibegin;i++) {vsp->SetBinContent(i+1,tdataD[i+ibegin]);}
      } else {
         Warning("Draw","Type of data not valid for ploting, returning 0");
         delete vsp;
         return ;
      }
   
//     Do not put the plot into any directory to avoid name conflicts
      vsp->SetDirectory(0);
   
      if (!opt.Contains("goff")) vsp->Draw(option);
      
   } else if (vect->nDim==2) {
// 2 dim case
      vspybegin = vect->startX[1];
      if (vect->dx[1]) {
         vspyend = vspybegin + vect->nx[1]*vect->dx[1];
      } else {
         if (gVEGAWarnings)
            Warning("Draw","Null length in y for serie. Using length=ndata");
         vspyend = vspybegin + vect->nx[1];
      }
      vsp2 = new VSPlot2(vectname,vectname,iend-ibegin,0,vspend - vspbegin,vect->nx[1],vspybegin,vspyend,this);
   // These settings should be set BEFORE the GPS start because setting
   // GPS start uses them.
//      vsp2->SetOffset(vspoff);
      vsp2->SetStepX(vect->dx[0]);
      vsp2->SetStepY(vect->dx[1]);
      vsp2->SetLocalTime(vect->localTime);
      vsp2->SetLeapS(vect->ULeapS);
#if FR_VERS<6000
      vsp2->SetGPSStart(vect->GTimeS + 10e-9*vect->GTimeN + vspbegin);
#else
      vsp2->SetGPSStart(vect->GTime + vspbegin);
#endif
      if (vect->unitX) {
         if (vect->unitX[0]) {
            if (strcmp(vect->unitX[0],"NONE"))vsp2->GetXaxis()->SetTitle(vect->unitX[0]);
         }
         if (vect->unitX[0]) {
            if (strcmp(vect->unitX[1],"NONE"))vsp2->GetYaxis()->SetTitle(vect->unitX[1]);
         }
      }
      if (vect->unitY) {
         if (strcmp(vect->unitY,"NONE")) vsp2->GetZaxis()->SetTitle(vect->unitY);
      }
      
      unsigned short typevect = vect->type;
      if (typevect == FR_VECT_C) {
         unsigned char* tdata= (unsigned char*) vect->data;
         for (ix=0;ix<iend-ibegin;ix++) {
            for (iy=0;iy<vect->nx[1];iy++) {
// The data order of elements storing in the array corresponding to a
// histogram is the opposite of the C convention ! Have to exchange X and Y
               vsp2->SetBinContent(ix+1+(iy+1)*(vect->nx[0]+2),tdata[ix+ibegin+iy*vect->nx[0]]);}}
      } else if (typevect == FR_VECT_2U) {
         unsigned short* tdataUS= vect->dataUS;
         for (ix=0;ix<iend-ibegin;ix++) {
            for (iy=0;iy<vect->nx[1];iy++) {
               vsp2->SetBinContent(ix+1+(iy+1)*(vect->nx[0]+2),tdataUS[ix+ibegin+iy*vect->nx[0]]);}}
      } else if (typevect == FR_VECT_2S) {
         short* tdataS= vect->dataS;
         for (ix=0;ix<iend-ibegin;ix++) {
            for (iy=0;iy<vect->nx[1];iy++) {
               vsp2->SetBinContent(ix+1+(iy+1)*(vect->nx[0]+2),tdataS[ix+ibegin+iy*vect->nx[0]]);}}
      } else if (typevect == FR_VECT_4U) {
         unsigned int* tdataUI= vect->dataUI;
         for (ix=0;ix<iend-ibegin;ix++) {
            for (iy=0;iy<vect->nx[1];iy++) {
               vsp2->SetBinContent(ix+1+(iy+1)*(vect->nx[0]+2),tdataUI[ix+ibegin+iy*vect->nx[0]]);}}
      } else if (typevect == FR_VECT_4S) {
         int* tdataI= vect->dataI;
         for (ix=0;ix<iend-ibegin;ix++) {
            for (iy=0;iy<vect->nx[1];iy++) {
               vsp2->SetBinContent(ix+1+(iy+1)*(vect->nx[0]+2),tdataI[ix+ibegin+iy*vect->nx[0]]);}}
      } else if (typevect == FR_VECT_4R) {
         float* tdataF= vect->dataF;
         for (ix=0;ix<iend-ibegin;ix++) {
            for (iy=0;iy<vect->nx[1];iy++) {
               vsp2->SetBinContent(ix+1+(iy+1)*(vect->nx[0]+2),tdataF[ix+ibegin+iy*vect->nx[0]]);}}
      } else if (typevect == FR_VECT_8R) {
         double* tdataD= vect->dataD;
         for (ix=0;ix<iend-ibegin;ix++) {
            for (iy=0;iy<vect->nx[1];iy++) {
               vsp2->SetBinContent(ix+1+(iy+1)*(vect->nx[0]+2),tdataD[ix+ibegin+iy*vect->nx[0]]);}}
      } else {
         Warning("Draw","Type of data not valid for plotting, returning 0");
         delete vsp2;
         return ;
      }

//     Set the number of entries not null so that the plot is displayed correctly
      vsp2->SetEntries(1);
      
//     Do not put the plot into any directory to avoid name conflicts
      vsp2->SetDirectory(0);
   
      if (!opt.Contains("goff")) vsp2->Draw(option);
   } else {
      Warning("Draw","Cannot plot %d dimensions object",vect->nDim);
   }
   
//*-*- Create a default canvas if none exists
//    if (!gPad && !opt.Contains("goff")) {
//    	if (!gROOT->GetMakeDefCanvas()) return;
//    	(gROOT->GetMakeDefCanvas())();
//    }

// Set the default options for a "T" and "S" serie
   if (!opt.Contains("goff") && gVStyle->GetAutoLogPlot()!=0) {
      if (VGetTypeOfVect(vect)=="S" && vect->nDim==1) {
         gPad->SetLogx(1);
         gPad->SetLogy(1);
      } else {
         gPad->SetLogx(0);
         gPad->SetLogy(0);
      }
   }
}

//______________________________________________________________________________
 void VManagerFrameL::Draw(FrVect* vect, Option_t *option)
{
//         Create VSPlot and Draw it starting from an FrVect
//         =================================================
// The entire series is plotted.
// For the moment, options are the same as for histograms
   
   Draw(vect,-1e20,1e20,option);
}

//______________________________________________________________________________
 void VManagerFrameL::Draw(FrameH* frame, char *nameofvect, Double_t offset, Double_t dur, Option_t *option)
{
//     Create VSPlot and Draw it starting from a vector in a frame
//     ===========================================================
// The series (FrAdcData or FrProcData or FrSimData)
// is plotted from "offset" for a length "dur". The default values
// mean that the entire series is plotted.
// The string nameofvect indicates the type (adc, proc, sim) and
// the name of the series to be extracted. 
// See Draw(FrVect....) for more information.
// For the moment, options are the same as for histograms

   FrVect* drawvector;

// Protect against bad values
   if (!frame) return;
   
   drawvector = VGetVect(frame,nameofvect);
   if (drawvector) {
      Draw(drawvector,offset,dur,option);
   } else {
      Warning("DrawHist","No vector %s found in frame",nameofvect);
   }
}

//______________________________________________________________________________
 void VManagerFrameL::Draw(FrameH* frame, char *nameofvect, Option_t *option)
{
//     Create VSPlot and Draw it starting from a vector in a frame
//     ===========================================================
// The entire series (FrAdcData or FrProcData or FrSimData)
// is plotted.
// The string nameofvect indicates the type (adc, proc, sim) and
// the name of the serie to be extracted. 
// See VGetVect(FrameH....) for more information.
// For the moment, options are the same as for histograms
   
   Draw(frame,nameofvect,-1e20,1e20,option);
}

//______________________________________________________________________________
 void VManagerFrameL::DrawHist(FrVect* vect, Double_t offset, Double_t dur, const char *hwork, Option_t *option)
{
//   Create 1D Histogram of the data values of an FrVect vector and Draw it
//   ======================================================================
// The data of the series (FrVect vector) used to build the histogram goes from
// "offset" for a length "dur". The default values mean that the entire series
// is used.
// The name of an already existing histogram can be given, for example "histo1"
// If the name is preceeded by a "+", like "+histo1",
// the data is appended to this histogram.
// If this is not the case, the old histogram contents are erased before filling.
// If no name is given, the limits of the new histogram are the min and max
// of the data.
// Options are passed to the new histogram.
// New option "temp". If the option string contains the option "temp" and no
// name of an already existing histogram is given, the
// created histogram will have the name "htemp". The created histogram will be
// outside any directory.
// Otherwise, the name of the FrVect will be used for the histogram.
// It will replace any existing histogram in gDirectory that has the same name

   TString  opt;
   opt = option;
   opt.ToLower();

   char *hdefault = (char *)"htemp";
   char *hname, *namevect;
   char *hname0=0;
   Int_t i,j,hkeep;
   TH1 *oldh1 = 0;
   TH1 *oldhtemp = 0;
   Bool_t profile = kFALSE;
            
   Int_t ibegin,iend,ndat;
   Double_t vsoffset,vsdur,vspbegin, vspend, bining;
   Double_t minvs,maxvs;
   Int_t nbinsvs;
   
             char* tdata;
   unsigned short* tdataUS;
            short* tdataS;
   unsigned int  * tdataUI;
            int  * tdataI;
            float* tdataF;
           double* tdataD;

// Protect against bad values
   if (!vect) return;
   
// Get info for name, number of samples, offset and length of the vector
   i = strlen(vect->name);
   namevect = new char[i+1];
   strcpy(namevect,vect->name);  // Copy the name of the vector
   ndat = vect->nData;
   if (!ndat) {
      Error("DrawHist","No data in this series, cannot build the histogram");
      if (namevect) delete [] namevect;
      return ;
   }
   vsoffset = 0;
   vsdur = VGetLength(vect);
   if (vsdur==0) {
      vsdur = (float)ndat;
      if (gVEGAWarnings)
         Warning("DrawHist","Null length for serie. Using length=ndata");
   }

// Get the name of the histogram and a pointer to it
// Beware of a potential memory leak on the name of the histogram !
   i = strlen(hwork);
   hname0 = new char[i+1]; // Use a temp variable to avoid mem leak
   strncpy(hname0,hwork,i); hname0[i]=0;
   hname = hname0;
   if (strcmp(hname,"")) {
      hkeep  = 1;
      Bool_t hnameplus = kFALSE;
      while (*hname == ' ') hname++;
      if (*hname == '+') {
         hnameplus = kTRUE;
         hname++;
         while (*hname == ' ') hname++;
         j = strlen(hname)-1;
         while (j) {
            if (hname[j] != ' ') break;
            hname[j] = 0;
            j--;
         }
      }       
      oldh1 = (TH1*)gDirectory->Get(hname);
      if (oldh1 && !hnameplus) oldh1->Reset();
   } else if (opt.Contains("temp")) {
      hname  = hdefault;
      hkeep  = 0;
      if (gDirectory) {
         oldh1 = (TH1*)gDirectory->Get(hname);
         if (oldh1) { oldh1->Delete(); oldh1 = 0;}
      }
   } else {
      hname  = namevect;
      hkeep  = 1;
      if (gDirectory) {
         oldh1 = (TH1*)gDirectory->Get(hname);
         if (oldh1) { oldh1->Delete(); oldh1 = 0;}
      }
   }

// Delete the histo named htemp if there is one displayed in the active pad
// This is to avoid a memory leak since htemp is always SetDirectory(0)
// -> cannot access it in the global lists of objects
   if (gPad) {
      oldhtemp = (TH1*)gPad->GetPrimitive(hdefault);
      if (oldhtemp) { oldhtemp->Delete(); oldhtemp = 0;}
   }

// In case old histo exists, check dimensionality
   if (oldh1) {
      Int_t mustdelete = 0;
      if (oldh1->InheritsFrom("TProfile")) profile = kTRUE;
      if (opt.Contains("prof")) {
         if (!profile) mustdelete = 1;
      } else {
         if (oldh1->GetDimension() != 1) mustdelete = 1;
      }
      if (mustdelete) {
         if (gVEGAWarnings)
           Warning("DrawHist","Deleting old histogram with different dimensions");
         delete oldh1; oldh1 = 0;
      }
   }

// Define begin and end of the used portion of the serie.
// Correct it if necessary
   vspbegin = vsoffset;
   if (offset>-9e19) {vspbegin += offset;}
   if (dur<9e19) {vspend = vspbegin+dur;} else {vspend=vspbegin+vsdur;}
   if (vspbegin<vsoffset) {vspbegin = vsoffset;}
   if (vspend>vsoffset+vsdur) {vspend = vsoffset+vsdur;}
   
   bining = vsdur/ndat;

// Compute the numbers of first and last bin to be used
   ibegin = int((vspbegin-vsoffset)/bining);
   iend   = int((vspend-vsoffset)/bining);

//  Get min and max of the serie
   minvs = 1e20; maxvs = -1e20;

   unsigned short typevect = vect->type;
   if (typevect == FR_VECT_C) {
      tdata=vect->data;
      for (i=ibegin;i<iend;i++) {
         if (tdata[i]<minvs) minvs=tdata[i];
         if (tdata[i]>maxvs) maxvs=tdata[i];}
   } else if (typevect == FR_VECT_2U) {
      tdataUS=vect->dataUS;
      for (i=ibegin;i<iend;i++) {
         if (tdataUS[i]<minvs) minvs=tdataUS[i];
         if (tdataUS[i]>maxvs) maxvs=tdataUS[i];}
   } else if (typevect == FR_VECT_2S) {
      tdataS= vect->dataS;
      for (i=ibegin;i<iend;i++) {
         if (tdataS[i]<minvs) minvs=tdataS[i];
         if (tdataS[i]>maxvs) maxvs=tdataS[i];}
   } else if (typevect == FR_VECT_4U) {
      tdataUI= vect->dataUI;
      for (i=ibegin;i<iend;i++) {
         if (tdataUI[i]<minvs) minvs=tdataUI[i];
         if (tdataUI[i]>maxvs) maxvs=tdataUI[i];}
   } else if (typevect == FR_VECT_4S) {
      tdataI= vect->dataI;
      for (i=ibegin;i<iend;i++) {
         if (tdataI[i]<minvs) minvs=tdataI[i];
         if (tdataI[i]>maxvs) maxvs=tdataI[i];}
   } else if (typevect == FR_VECT_4R) {
      tdataF= vect->dataF;
      for (i=ibegin;i<iend;i++) {
         if (tdataF[i]<minvs) minvs=tdataF[i];
         if (tdataF[i]>maxvs) maxvs=tdataF[i];}
   } else if (typevect == FR_VECT_8R) {
      tdataD= vect->dataD;
      for (i=ibegin;i<iend;i++) {
         if (tdataD[i]<minvs) minvs=tdataD[i];
         if (tdataD[i]>maxvs) maxvs=tdataD[i];}
   } else {
      Warning("DrawHist","Type of data not valid for histograming");
      if (hname0) delete hname0;
      if (namevect) delete namevect;
      return;
   }

   if (minvs>9e19 || maxvs <-9e19) {minvs=0; maxvs=1;}

// If the option "same" is used, try to make a histogram with correct bounds
// i.e. those of the previously drawn histogram.
   if (!oldh1) {
      nbinsvs = 100;
      if (opt.Contains("same")) {
          TH1 *oldhtemp = (TH1*)gPad->GetPrimitive(hdefault);
          if (oldhtemp) {
             nbinsvs = oldhtemp->GetXaxis()->GetNbins();
             minvs   = oldhtemp->GetXaxis()->GetXmin();
             maxvs   = oldhtemp->GetXaxis()->GetXmax();
          } else {
             minvs   = gPad->GetUxmin();
             maxvs   = gPad->GetUxmax();
          }
      }
   }

// Build a new histogram if necessary
   TH1F *h1;
   if (oldh1) {
      h1 = (TH1F*)oldh1;
   } else {
      h1 = new TH1F(hname,hname,nbinsvs,minvs-0.05*(maxvs-minvs),maxvs+0.05*(maxvs-minvs));
      if (!hkeep) {
         h1->SetBit(kCanDelete);
         h1->SetDirectory(0);
      }
   }
   if (hname0) delete [] hname0;
   if (namevect) delete [] namevect;

// Fill the histogram with serie values
   if (typevect == FR_VECT_C) {
      tdata= vect->data;
      for (i=ibegin;i<iend;i++) {h1->Fill(tdata[i]);}
   } else if (typevect == FR_VECT_2U) {
      tdataUS= vect->dataUS;
      for (i=ibegin;i<iend;i++) {h1->Fill(tdataUS[i]);}
   } else if (typevect == FR_VECT_2S) {
      tdataS= vect->dataS;
      for (i=ibegin;i<iend;i++) {h1->Fill(tdataS[i]);}
   } else if (typevect == FR_VECT_4U) {
      tdataUI= vect->dataUI;
      for (i=ibegin;i<iend;i++) {h1->Fill(tdataUI[i]);}
   } else if (typevect == FR_VECT_4S) {
      tdataI= vect->dataI;
      for (i=ibegin;i<iend;i++) {h1->Fill(tdataI[i]);}
   } else if (typevect == FR_VECT_4R) {
      tdataF= vect->dataF;
      for (i=ibegin;i<iend;i++) {h1->Fill(tdataF[i]);}
   } else if (typevect == FR_VECT_8R) {
      tdataD= vect->dataD;
      for (i=ibegin;i<iend;i++) {h1->Fill(tdataD[i]);}
   } else {
      Warning("DrawHist","Type of data not valid for histograming");
      return;
   }

//*-*- Create a default canvas if none exists
//    if (!gPad && !gProofServ && !opt.Contains("goff")) {
//    	if (!gROOT->GetMakeDefCanvas()) return;
//    	(gROOT->GetMakeDefCanvas())();
//    }

// Draw it
   opt.ReplaceAll("temp","");  // remove the option temp if it exists
   if (!opt.Contains("goff")) h1->Draw(opt.Data());
}

//______________________________________________________________________________
 void VManagerFrameL::DrawHist(FrVect* vect, const char *hwork, Option_t *option)
{
//   Create 1D Histogram of the data values of an FrVect vector and Draw it
//   ======================================================================
// All the data of the series (FrVect vector) are used to build the histogram.
// The name of an already existing histogram can be given, for example "histo1"
// If the name is preceeded by a "+", like "+histo1",
// the data is appended to this histogram.
// If this is not the case, the old histogram contents are erased before filling.
// If no name is given, the limits of the new histogram are the min and max
// of the data.
// Options are passed to the new histogram.
   
   DrawHist(vect, -1e20,1e20,hwork,option);
}

//______________________________________________________________________________
 void VManagerFrameL::DrawHist(FrameH* frame, char* nameofvect, Double_t offset, Double_t dur, const char *hwork, Option_t *option)
{
// Creates 1D Histogram of the data values of a vector in a Frame and Draw it
// ==========================================================================
// The data of the series (FrAdcData or FrProcData or FrSimData)
// used to build the histogram goes from "offset" for a length "dur".
// The string nameofvect indicates the type and the name of the series
// to be extracted.
// The name of an already existing histogram can be given, for example "histo1"
// If the name is preceeded by a "+", like "+histo1",
// the data is appended to this histogram.
// If this is not the case, the old histogram contents are erased before filling.
// If no name is given, the limits of the new histogram are the min and max
// of the data.
// Options are passed to the new histogram.

   FrVect* drawvector = 0;

// Protect against bad values
   if (!frame) return;
   
   drawvector = VGetVect(frame,nameofvect);
   if (drawvector) {
      DrawHist(drawvector,offset,dur,hwork,option);
   } else {
      Warning("DrawHist","No vector %s found in frame",nameofvect);
   }
}

//______________________________________________________________________________
 void VManagerFrameL::DrawHist(FrameH* frame, char* nameofvect, const char *hwork, Option_t *option)
{
// Creates 1D Histogram of the data values of a vector in a Frame and Draw it
// ==========================================================================
// All the data of the series (FrAdcData or FrProcData or FrSimData)
// is used to build the histogram.
// The string nameofvect indicates the type (adc, proc, sim) and
// the name of the serie to be extracted.
// See VGetVect(FrameH...) for more information.
// The name of an already existing histogram can be given, for example "histo1"
// If the name is preceeded by a "+", like "+histo1",
// the data is appended to this histogram.
// If this is not the case, the old histogram contents are erased before filling.
// If no name is given, the limits of the new histogram are the min and max
// of the data.
// Options are passed to the new histogram.
   
   DrawHist(frame,nameofvect,-1e20,1e20,hwork,option);
}

//______________________________________________________________________________
 TH1F* VManagerFrameL::GetLastHisto()
{
// ---- Returns the last histogram calculated and displayed ----
//
   TObject* hobj;
   hobj = gDirectory->GetList()->Last();
   if (hobj->InheritsFrom("TH1")) return (TH1F*)hobj;
   return 0;
}

//______________________________________________________________________________
 void VManagerFrameL::ShiftPlot(char* name, FrVect* vect)
{
// Shift the contents of the VSPlot. The beginning of the plot
// is removed while the new vector is appended to the end
}

//______________________________________________________________________________
 void VManagerFrameL::RemovePlot(const char* name)
{
// Removes the plot named "name" from the list of plots handled
// by this manager

   TObject* vsp;
   vsp = mPlots->FindObject(name);
   mPlots->Remove(vsp);
}

//______________________________________________________________________________
 void VManagerFrameL::RemovePlot2(const char* name)
{
// Removes the plot named "name" from the list of plots handled
// by this manager

   TObject* vsp2;
   vsp2 = mPlots2->FindObject(name);
   mPlots2->Remove(vsp2);
}

//______________________________________________________________________________
 TAxis* VManagerFrameL::GetXaxis(const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) return vs->GetXaxis(); else return 0;
}

//______________________________________________________________________________
 TAxis* VManagerFrameL::GetYaxis(const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) return vs->GetYaxis(); else return 0;
}

//______________________________________________________________________________
 TAxis* VManagerFrameL::GetZaxis(const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) return vs->GetZaxis(); else return 0;
}

//______________________________________________________________________________
 void VManagerFrameL::SetBarOffset(Float_t baroff, const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) vs->SetBarOffset(baroff);
}

//______________________________________________________________________________
 void VManagerFrameL::SetBarWidth(Float_t barwidth, const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) vs->SetBarWidth(barwidth);
}

//______________________________________________________________________________
 void VManagerFrameL::SetFillColor(Color_t color, const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) vs->SetFillColor(color);
}

//______________________________________________________________________________
 void VManagerFrameL::SetFillStyle(Style_t styl, const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) vs->SetFillStyle(styl);
}

//______________________________________________________________________________
 void VManagerFrameL::SetLineColor(Color_t color, const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) vs->SetLineColor(color);
}

//______________________________________________________________________________
 void VManagerFrameL::SetLineStyle(Style_t styl, const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) vs->SetLineStyle(styl);
}

//______________________________________________________________________________
 void VManagerFrameL::SetLineWidth(Width_t width, const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) vs->SetLineWidth(width);
}

//______________________________________________________________________________
 void VManagerFrameL::SetLabelColor(Color_t color, Option_t* axis, const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) vs->SetLabelColor(color);
}

//______________________________________________________________________________
 void VManagerFrameL::SetLabelFont(Style_t font, Option_t* axis, const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) vs->SetLabelFont(font);
}

//______________________________________________________________________________
 void VManagerFrameL::SetLabelOffset(Float_t offset, Option_t* axis, const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) vs->SetLabelOffset(offset);
}

//______________________________________________________________________________
 void VManagerFrameL::SetLabelSize(Float_t size, Option_t* axis, const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) vs->SetLabelSize(size);
}

//______________________________________________________________________________
 void VManagerFrameL::SetMarkerColor(Color_t tcolor, const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) vs->SetMarkerColor(tcolor);
}

//______________________________________________________________________________
 void VManagerFrameL::SetMarkerSize(Size_t msize, const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) vs->SetMarkerSize(msize);
}

//______________________________________________________________________________
 void VManagerFrameL::SetMarkerStyle(Style_t mstyle, const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) vs->SetMarkerStyle(mstyle);
}

//______________________________________________________________________________
 void VManagerFrameL::SetNdivisions(Int_t n, Option_t* axis, const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) vs->SetNdivisions(n,axis);
}

//______________________________________________________________________________
 void VManagerFrameL::SetTickLength(Float_t length, Option_t* axis, const char* name)
{
   VSPlot *vs;
   if (!strlen(name)) vs = GetLastPlot(); else vs = GetPlot(name);
   if (vs) vs->SetTickLength(length, axis);
}



- 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.