Diana Software
ScanEvents.C

example on how to read the output root files produced by MRootFileWriter

Author
Marco Vignati

This macro is foreseen to read files produced by diana_tutorial.cfg and other tutorials. Run it directly on the output file:

root[0] .L pkg/modtutorial/ScanEvents.C
root[1] ScanEvents("data/Anal_000810_C_p001.root")

or on the .list file

root[0] .L pkg/modtutorial/ScanEvents.C
root[1] ScanEvents("data/Analyzed.list")

Optionally one can also compile the macro with ACLiC, the CINT interface to the C++ compiler.

root[0] .x lib/diana-root.C
root[1] .L pkg/modtutorial/ScanEvents.C++
root[2] ScanEvents("data/Analyzed.list")
// include files below are needed only if the macro is being compiled with ACLiC
#if !defined(__CINT__) || defined(__MAKECINT__)
#include "TCanvas.h"
#include "TGraph.h"
#include "TFile.h"
#include "TAxis.h"
// files below can be included only if diana-root.C has been executed,
// which sets the DIANA software include paths.
#include "QChain.hh"
#include "QHeader.hh"
#include "QPulse.hh"
#include "QPulseInfo.hh"
#include "QBaseType.hh"
#include "QBaselineData.hh"
#include "QBool.hh"
#endif
// Avoid putting Diana:: in front of objects defined in the Diana namespace
using namespace Diana;
using namespace std;
void ScanEvents(const char* filename)
{
// QChain is a TChain of TTrees named "qtree".
// Both .list and .root files can be added to a QChain
QChain* chain = new QChain();
if(!chain->Add(filename)){
cout<<"Error while reading file: "<<filename<<endl;
return;
}
// Set branches of raw event QObject%s and of analysis QObject%s if present.
// Remember to add a dot after the branch name.
// Quantities below are the raw event
QHeader* header = 0;
QPulse* pulse = 0;
QPulseInfo* pulseInfo = 0;
chain->SetBranchAddress("DAQ@Header.",&header);
chain->SetBranchAddress("DAQ@Pulse.",&pulse);
chain->SetBranchAddress("DAQ@PulseInfo.",&pulseInfo);
// Quantities below are computed in diana_tutorial.cfg
QDouble* maxBaseline = 0;
QBaselineData* baselineData = 0;
QCountPulsesData* countPulses = 0;
if(chain->FindBranch("PulseBasicParameters@MaxBaseline.")) {
chain->SetBranchAddress("PulseBasicParameters@MaxBaseline.",&maxBaseline);
chain->SetBranchAddress("BaselineModule@BaselineData.",&baselineData);
chain->SetBranchAddress("BCountPulses@CountPulsesData.",&countPulses);
}
// Quantities below are computed in diana_tutorial_averagepulse.cfg
QBool* badPulse = 0;
QBool* amplFilter = 0;
if(chain->FindBranch("BadPulse@Passed.")) {
chain->SetBranchAddress("BadPulse@Passed.",&badPulse);
chain->SetBranchAddress("TutorialAmplitudeFilter@Passed.",&amplFilter);
}
// Get number of entries in the chain
Long64_t nentries = chain->GetEntries();
if(nentries == 0) {
cout<<"Empty files!"<<endl;
return;
}
cout<<"Number of entries: "<<nentries<<endl;
// Get first entry to determine run number
chain->GetEntry(0);
Int_t run = header->GetRun();
// use the QGlobalDataManager to read the QRunData object contained in the output file
dm.SetOwner("ScanEvents");
dm.EnableCache(false);
dm.Get("DAQ",&rHandle,chain->GetFile()->GetName());
// if not present in files, try to load from DB.
if(!rHandle.IsValid()) {
dm.Get("DAQ",&rHandle,"DB");
}
// if not found print an error
if(!rHandle.IsValid()) {
cout<<rHandle.GetError()<<endl;
cout<<"RunData not valid, not converting from ADC to mV"<<endl;
}
// Create the pointer for the TGraph used to draw raw samples
TGraph* pulseGraph = 0;
// Create the canvas for the raw pulse
TCanvas* pulseCanvas = new TCanvas("PulseCanvas","Raw pulse");
// loop on TChain entries
for(Long64_t entry = 0; entry < nentries; entry++) {
// get QChain entry
chain->GetEntry(entry);
// Print on screen the raw event and Max-Baseline
cout<<"##############################"<<endl;
cout<<"Entry: "<<entry<<endl;
cout<<"******* Dump of DAQ@Header:"<<endl;
cout<<*header;
cout<<"******* Dump of DAQ@PulseInfo:"<<endl;
cout<<*pulseInfo;
cout<<"******* Dump of DAQ@Pulse"<<endl;
cout<<*pulse;
// show maxBaseline only if present
if(maxBaseline) {
cout<<"******* Dump of analysis quantities:"<<endl;
cout<<"Max-Baseline : IsValid = "<<maxBaseline->IsValid()<<" Value: "<<*maxBaseline<<" [mV]"<<endl;
cout<<"Baseline slope : IsValid = "<<baselineData->IsValid()<<" Value: "<<baselineData->fBaselineSlope<<" [mV/ADC sample]"<<endl;
cout<<"Number of pulses: IsValid = "<<countPulses->IsValid()<<" Value: "<<countPulses->GetNumberOfPulses()<<endl;
}
if(badPulse) {
cout<<"BadPulse@Passed : IsValid = "<<badPulse->IsValid()<<" Value: "<<*badPulse<<endl;
cout<<"AmplFilt@Passed : IsValid = "<<amplFilter->IsValid()<<" Value: "<<*amplFilter<<endl;
}
// Get ADC2mV from QChannelRunData object of this channel, which is contained in QRunData.
double ADC2mV = -1;
if(rHandle.IsValid()) {
const QRunData& runData = rHandle.Get();
const QChannelRunData& chanRunData = runData.GetChannelRunData(pulseInfo->GetChannelId());
ADC2mV = chanRunData.fADC2mV;
cout<<"ADC2mV = "<<ADC2mV<<endl;
}
// Draw the raw pulse.
// Copy the raw samples, which are stored in a QVectorI (int elements),
// into a QVector (double elements), because QVector has more functionalities than QVectorI.
QVector samplesmV;
samplesmV = pulse->GetSamplesADC();
if(ADC2mV > 0) samplesmV *= ADC2mV;
// Draw the pulse via QVector::GetGraph(). The TGraph object
// returned is owned by the caller so we have to delete
// the TGraph of the previous iteration every time.
if(pulseGraph) delete pulseGraph;
pulseGraph = samplesmV.GetGraph();
pulseGraph->GetXaxis()->SetTitle("sample index");
if(ADC2mV > 0) pulseGraph->GetYaxis()->SetTitle("mV");
else pulseGraph->GetYaxis()->SetTitle("ADC counts");
pulseGraph->Draw("AL");
// Refresh canvas content after drawing
pulseCanvas->Modified();
pulseCanvas->Update();
// wait for a user hit before switching to next event
cout<<"****************************"<<endl;
cout<<"Press any key to go to next event, \'q\' to interrupt"<<endl;
char p = getchar();
if(p == 'q') return;
}
}
QRunDataHandle rHandle(753)
QChannelRunData chanRunData
QGlobalDataManager dm
void ScanEvents(const char *filename)
Definition: ScanEvents.C:52
base types wrapped into a QObject. Currently implemented types are QInt QDouble and QFloat....
Definition: QBaseType.hh:17
baseline data
double fBaselineSlope
bool wrapped into a QObject
Definition: QBool.hh:17
TChain used in diana.
Definition: QChain.hh:23
Int_t Add(const char *filename, Long64_t nentries=kBigNumber)
Definition: QChain.cc:50
basic channel and run based info. Used in the QRunData object.
double fADC2mV
conversion: mV = ADC * fADC2mV
number of pulses and time interval beetwen peaks in the same acquired window
const int & GetNumberOfPulses() const
Object to manage I/O (DB, file, or memory) of diana global quantities.
void SetOwner(const std::string &owner)
set the module that is accessing this object
void EnableCache(const bool enable)
enable / disable caching of objects
QError Get(const std::string &owner, GlobalHandle< Q > *gh, const std::string &inSource, bool printError=true) const
Get an object using a global handle.
Raw event: basic information like run number and time.
Definition: QHeader.hh:16
int GetRun() const
destructor
Definition: QHeader.hh:22
bool IsValid() const
check wheter object is valid
Definition: QObject.hh:114
Raw event: bolometer channel, trigger positions and types.
Definition: QPulseInfo.hh:18
Raw event: sampled waveform.
Definition: QPulse.hh:22
global handle for QRunData
Basic run based info.
Definition: QRunData.hh:20
const QChannelRunData & GetChannelRunData(const int channel) const
get channel based run data quantities
Definition: QRunData.cc:339
Interface for vectors in Diana analysis.
Definition: QVector.hh:30
TGraph * GetGraph(double samplingFrequency=1., double scale=1.) const
Get a root TGraph (owned by the caller)
Definition: QVector.cc:761
the Diana namespace is needed because sometimes we use Qt libraries, that use same class names of our...