Logo coherent WaveBurst  
Library Reference Guide
Logo
DrawGWaveArray.C
Go to the documentation of this file.
1 // This example shows how to produce time,freq,time-freq plots of
2 // a signal stored in a wavearray structure
3 
4 {
5  //#define SAVE_PLOT // uncomment to save plot to file
6 
7  #define PLOT_TYPE GWAT_TIME // draw time shifted signal (black)
8  //#define PLOT_TYPE GWAT_FFT // draw on the same plot the original signal (red)
9  //#define PLOT_TYPE GWAT_TF // get pointer to watplot object
10 
11  #define SIG_FREQ 10 // signal frequency
12  #define TSHIFT 0.05 // time shift (sec)
13 
14 
15  gwavearray<double> gw(16384); // instantiation with size 16384 samples
16  gw.rate(16384); // set sample rate = 16384 Hz
17  gw.start(0); // set start time = 0
18 
19  // fill with a sinusoid @ 100 Hz
20  double dt=1/gw.rate();
21  for(int i=0;i<gw.size();i++) gw[i]=sin(2*PI*SIG_FREQ*dt*i);
22 
23  // draw methods
24 
25  //gw.Draw(GWAT_TIME); // draw signal in time domain
26  //gw.Draw(GWAT_FFT); // draw signal in frequency domain
27  //gw.Draw(GWAT_TF); // draw signal in time-frequency domain
28 
29 
30  // apply a time shift of 0.005 sec
31  //gw.TimeShift(TSHIFT);
32  //gw.Draw(GWAT_TIME,"SAME",kRed);
33 
34  // create a new time shifted signal @ SIG_FREQ Hz
36  u.resize(16384);
37  u.rate(16384);
38  u.start(0);
39  for(int i=0;i<u.size();i++) u[i]=sin(2*PI*SIG_FREQ*(dt*i+TSHIFT));
40 
41  // create a new gwavearray using the wavearray u
42  gwavearray<double> gu(&u);
43  //gwavearray<double> gu = u;
44 
45  gu.Draw(PLOT_TYPE); // draw time shifted signal (black)
46  gu.Draw(&gw,PLOT_TYPE,"SAME",kRed); // draw on the same plot the original signal (red)
47 
48 #ifdef SAVE_PLOT
49  if(PLOT_TYPE==GWAT_TF) {
50  CWB::STFT* stft = gu.GetSTFT();
51  stft->Print("gwavearray_plot.png");
52  }
53 #endif
54 
55  watplot* plot = gu.GetWATPLOT(); // get pointer to watplot object
56 
57  if((PLOT_TYPE==GWAT_TIME) || (PLOT_TYPE==GWAT_FFT)) {
58 
59  // set x,y axis titles and plot title
60  char gtitle[256];
61  sprintf(gtitle,"Original(Red) - Time Shifted(Black)");
62  if(PLOT_TYPE==GWAT_TIME) plot->gtitle(gtitle,"time(sec)","amplitude");
63  if(PLOT_TYPE==GWAT_FFT) plot->gtitle(gtitle,"frequency(hz)","amplitude");
64 
65  // save plot to file
66 #ifdef SAVE_PLOT
67  TString gfile="gwavearray_plot.png";
68  (*plot) >> gfile;
69 
70  cout << "created plot file name : " << gfile << endl;
71  exit(0);
72 #endif
73 
74  }
75 
76  return plot->canvas; // used by THtml
77 }
void gtitle(TString title="", TString xtitle="", TString ytitle="")
Definition: watplot.cc:1256
virtual void rate(double r)
Definition: wavearray.hh:141
void Print(TString pname)
Definition: STFT.cc:315
TString("c")
plot gtitle(gtitle,"frequency (Hz)","strain/#sqrt{Hz}")
virtual void start(double s)
Definition: wavearray.hh:137
i drho i
#define PI
Definition: watfun.hh:32
virtual size_t size() const
Definition: wavearray.hh:145
TCanvas * canvas
Definition: watplot.hh:192
gWSeries< double > gw(w)
x plot
double dt
#define TSHIFT
#define SIG_FREQ
#define PLOT_TYPE
Definition: gwat.hh:30
Definition: gwat.hh:31
sprintf(tfres,"(1/%g)x(%g) (sec)x(Hz)", 2 *df, df)
TString gfile
virtual void resize(unsigned int)
Definition: wavearray.cc:463
CWB::STFT * stft
Definition: ChirpMass.C:121
exit(0)