GET library
RCoBo_ChannelMeanRMS.C

Example macro file that reads a reduced-CoBo file to build the average and the RMS for each channel (only one is displayed).

RCoBo_ChannelMeanRMS.png

This example uses the functions from the GETRoot library. It can be executed in a ROOT session with:

* [ROOT] .L RCoBo_ChannelMeanRMS.C
* [ROOT] PlotRCoBo_ChannelMeanRMS ( "../../data/RCoBo/CoBo_baseline.graw", 1, 17 );
*
//======================================================================
// ROOT macros for RCoBo system
// - read a test file (from reduced CoBo)
// - compute mean sample for all channels
// - compute RMS sample for all channels
// - display result for 1 channel
//
//======================================================================
/*
These macros require to load the GET libraries
gROOT->Macro ( (getenv("GET_BASE")+string("/root/GETRootLibs.C")).c_str() );
To execute the macro
.L RCoBo_ChannelMeanRMS.C
PlotRCoBo_ChannelMeanRMS ( "../../data/RCoBo/CoBo_baseline.graw", 1, 17 );
*/
//----------------------------------------------------------------------
// This function reads all events from a R-CoBo run file
// (a baseline run).
// For all channels, it creates the average signal.
// The the file is read a second time to compute the RMS of the
// signal.
//
// The average and the RMS are computed using functions from the
// GETRoot aditionnal library.
//
// The result is plotted for one channel.
void PlotRCoBo_ChannelMeanRMS ( const string & data_file, u_short aget_id, u_short chan_id )
{
GSetVerboseLevel(3);
getProcessedInfoCount = 5;
// create the GET system (reduced CoBo)
// - 512 samples / channel
// - 0.01 us sampling time
u_int n = 512;
double dt = 0.01;
GETRCoBo rcobo(n,dt);
// get the reference for the channel to plot
GETSample & sample_out = rcobo[0][0][aget_id][chan_id].OutSample();
//------------------------------------------------------------
// Compute the average for all channels
cout << endl
<< "Computing averages" << endl;
// call GET/ROOT function to build average samples
// (all events, no output sample processing)
u_int nr = GETRCoBoMean ( rcobo, data_file, 0, false );
cout << endl
<< "Number of events read: " << nr << endl
<< endl;
if (nr > 0)
{
// memorize the sample for the channel to plot
// (the data in rcobo object will be overwritten when computing RMS...)
GETSample sample_mean ( sample_out );
//------------------------------------------------------------
// Compute the RMS for all channels
cout << "Computing RMS" << endl
<< endl;
// call GET/ROOT function to build RMS samples, with respect
// to current content of the channels (mean values)
GETRCoBoStdDeviation ( rcobo, data_file, 0, false );
//------------------------------------------------------------
// Plot the results for one channel
// compute the mean sample +/- the RMS
GETSample sample_inf = sample_mean - sample_out;
GETSample sample_sup = sample_mean + sample_out;
// the GETSample is a RRealSampleFFT object (see GFFT documentation)
TGraph * gr_mean = sample_mean.CreateFunctionGraph ( );
gr_mean->SetTitle ( "Average + RMS channel output" );
gr_mean->SetLineColor ( kRed+1 );
gr_mean->SetLineWidth ( 2 );
gr_mean->SetLineStyle ( 1 );
gr_mean->GetXaxis()->SetRangeUser ( 0., (n-1)*dt );
TGraph * gr_inf = sample_inf.CreateFunctionGraph ( );
gr_inf->SetLineColor ( kGreen+1 );
gr_inf->SetLineWidth ( 1 );
gr_inf->SetLineStyle ( 1 );
TGraph * gr_sup = sample_sup.CreateFunctionGraph ( );
gr_sup->SetLineColor ( kGreen+1 );
gr_sup->SetLineWidth ( 1 );
gr_sup->SetLineStyle ( 1 );
TCanvas * canvas = new TCanvas ( "PlotRMS", "PLOTRMS", 800, 400 );
gr_mean->Draw ( "AL" );
gr_inf->Draw ( "L" );
gr_sup->Draw ( "L" );
// create an EPS file
gSystem->mkdir ( "plots", kTRUE );
canvas->SaveAs ( "plots/RCoBo_ChannelMeanRMS.eps" );
}
}