Diana Software
ScanEvents.C
Go to the documentation of this file.
1 
28 // include files below are needed only if the macro is being compiled with ACLiC
29 #if !defined(__CINT__) || defined(__MAKECINT__)
30 #include "TCanvas.h"
31 #include "TGraph.h"
32 #include "TFile.h"
33 #include "TAxis.h"
34 // files below can be included only if diana-root.C has been executed,
35 // which sets the DIANA software include paths.
36 #include "QChain.hh"
37 #include "QHeader.hh"
38 #include "QPulse.hh"
39 #include "QPulseInfo.hh"
40 #include "QBaseType.hh"
41 #include "QRunDataHandle.hh"
42 #include "QGlobalDataManager.hh"
43 #include "QBaselineData.hh"
44 #include "QCountPulsesData.hh"
45 #include "QBool.hh"
46 #endif
47 
48 // Avoid putting Diana:: in front of objects defined in the Diana namespace
49 using namespace Diana;
50 using namespace std;
51 
52 void ScanEvents(const char* filename)
53 {
54  // QChain is a TChain of TTrees named "qtree".
55  // Both .list and .root files can be added to a QChain
56  QChain* chain = new QChain();
57  if(!chain->Add(filename)){
58  cout<<"Error while reading file: "<<filename<<endl;
59  return;
60  }
61  // Set branches of raw event QObject%s and of analysis QObject%s if present.
62  // Remember to add a dot after the branch name.
63  // Quantities below are the raw event
64  QHeader* header = 0;
65  QPulse* pulse = 0;
66  QPulseInfo* pulseInfo = 0;
67  chain->SetBranchAddress("DAQ@Header.",&header);
68  chain->SetBranchAddress("DAQ@Pulse.",&pulse);
69  chain->SetBranchAddress("DAQ@PulseInfo.",&pulseInfo);
70  // Quantities below are computed in diana_tutorial.cfg
71  QDouble* maxBaseline = 0;
72  QBaselineData* baselineData = 0;
73  QCountPulsesData* countPulses = 0;
74  if(chain->FindBranch("PulseBasicParameters@MaxBaseline.")) {
75  chain->SetBranchAddress("PulseBasicParameters@MaxBaseline.",&maxBaseline);
76  chain->SetBranchAddress("BaselineModule@BaselineData.",&baselineData);
77  chain->SetBranchAddress("BCountPulses@CountPulsesData.",&countPulses);
78  }
79  // Quantities below are computed in diana_tutorial_averagepulse.cfg
80  QBool* badPulse = 0;
81  QBool* amplFilter = 0;
82  if(chain->FindBranch("BadPulse@Passed.")) {
83  chain->SetBranchAddress("BadPulse@Passed.",&badPulse);
84  chain->SetBranchAddress("TutorialAmplitudeFilter@Passed.",&amplFilter);
85  }
86 
87  // Get number of entries in the chain
88  Long64_t nentries = chain->GetEntries();
89  if(nentries == 0) {
90  cout<<"Empty files!"<<endl;
91  return;
92  }
93  cout<<"Number of entries: "<<nentries<<endl;
94 
95  // Get first entry to determine run number
96  chain->GetEntry(0);
97  Int_t run = header->GetRun();
98  // use the QGlobalDataManager to read the QRunData object contained in the output file
100  dm.SetOwner("ScanEvents");
101  dm.EnableCache(false);
102  QRunDataHandle rHandle(run);
103  dm.Get("DAQ",&rHandle,chain->GetFile()->GetName());
104  // if not present in files, try to load from DB.
105  if(!rHandle.IsValid()) {
106  dm.Get("DAQ",&rHandle,"DB");
107  }
108  // if not found print an error
109  if(!rHandle.IsValid()) {
110  cout<<rHandle.GetError()<<endl;
111  cout<<"RunData not valid, not converting from ADC to mV"<<endl;
112  }
113 
114  // Create the pointer for the TGraph used to draw raw samples
115  TGraph* pulseGraph = 0;
116  // Create the canvas for the raw pulse
117  TCanvas* pulseCanvas = new TCanvas("PulseCanvas","Raw pulse");
118 
119  // loop on TChain entries
120  for(Long64_t entry = 0; entry < nentries; entry++) {
121  // get QChain entry
122  chain->GetEntry(entry);
123 
124  // Print on screen the raw event and Max-Baseline
125  cout<<"##############################"<<endl;
126  cout<<"Entry: "<<entry<<endl;
127  cout<<"******* Dump of DAQ@Header:"<<endl;
128  cout<<*header;
129  cout<<"******* Dump of DAQ@PulseInfo:"<<endl;
130  cout<<*pulseInfo;
131  cout<<"******* Dump of DAQ@Pulse"<<endl;
132  cout<<*pulse;
133  // show maxBaseline only if present
134  if(maxBaseline) {
135  cout<<"******* Dump of analysis quantities:"<<endl;
136  cout<<"Max-Baseline : IsValid = "<<maxBaseline->IsValid()<<" Value: "<<*maxBaseline<<" [mV]"<<endl;
137  cout<<"Baseline slope : IsValid = "<<baselineData->IsValid()<<" Value: "<<baselineData->fBaselineSlope<<" [mV/ADC sample]"<<endl;
138  cout<<"Number of pulses: IsValid = "<<countPulses->IsValid()<<" Value: "<<countPulses->GetNumberOfPulses()<<endl;
139  }
140  if(badPulse) {
141  cout<<"BadPulse@Passed : IsValid = "<<badPulse->IsValid()<<" Value: "<<*badPulse<<endl;
142  cout<<"AmplFilt@Passed : IsValid = "<<amplFilter->IsValid()<<" Value: "<<*amplFilter<<endl;
143  }
144 
145  // Get ADC2mV from QChannelRunData object of this channel, which is contained in QRunData.
146  double ADC2mV = -1;
147  if(rHandle.IsValid()) {
148  const QRunData& runData = rHandle.Get();
149  const QChannelRunData& chanRunData = runData.GetChannelRunData(pulseInfo->GetChannelId());
150  ADC2mV = chanRunData.fADC2mV;
151  cout<<"ADC2mV = "<<ADC2mV<<endl;
152  }
153  // Draw the raw pulse.
154  // Copy the raw samples, which are stored in a QVectorI (int elements),
155  // into a QVector (double elements), because QVector has more functionalities than QVectorI.
156  QVector samplesmV;
157  samplesmV = pulse->GetSamplesADC();
158  if(ADC2mV > 0) samplesmV *= ADC2mV;
159  // Draw the pulse via QVector::GetGraph(). The TGraph object
160  // returned is owned by the caller so we have to delete
161  // the TGraph of the previous iteration every time.
162  if(pulseGraph) delete pulseGraph;
163  pulseGraph = samplesmV.GetGraph();
164  pulseGraph->GetXaxis()->SetTitle("sample index");
165  if(ADC2mV > 0) pulseGraph->GetYaxis()->SetTitle("mV");
166  else pulseGraph->GetYaxis()->SetTitle("ADC counts");
167  pulseGraph->Draw("AL");
168  // Refresh canvas content after drawing
169  pulseCanvas->Modified();
170  pulseCanvas->Update();
171 
172  // wait for a user hit before switching to next event
173  cout<<"****************************"<<endl;
174  cout<<"Press any key to go to next event, \'q\' to interrupt"<<endl;
175  char p = getchar();
176  if(p == 'q') return;
177  }
178 
179 }
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
the Diana namespace is needed because sometimes we use Qt libraries, that use same class names of our...