ROOTANA
Loading...
Searching...
No Matches
analyzer_example.cxx
Go to the documentation of this file.
1// Simple example of how to define a class that derives from
2// TRootanaEventLoop, thereby abstracting a way a lot of the
3// default stuff associated with analyzer.cxx
4
5#include <stdio.h>
6#include <iostream>
7
9#include "TH1D.h"
10#include "TV1190Data.hxx"
11#include "TV792Data.hxx"
12
13class MyTestLoop: public TRootanaEventLoop {
14
15 int nnn ;
16
17public:
18
20 nnn = 0;
21 };
22
23 virtual ~MyTestLoop() {};
24
25 void BeginRun(int transition,int run,int time){
26 std::cout << "Custom: begin run " << run << std::endl;
27 }
28
29 void EndRun(int transition,int run,int time){
30 std::cout << "Custom end run " << run <<std::endl;
31 }
32
34
35
36 void *ptr;
37 int size = dataContainer.GetMidasData().LocateBank(NULL, "FR10", &ptr);
38 if (ptr){
39 nnn++;
40 if(nnn%100 == 0){
41 std::cout << "Current run : " << GetCurrentRunNumber() << std::endl;
42 std::cout << "Trying to handle this event... " << size << " " << nnn << std::endl;
43 }
44
46 char sname[256];
47 sprintf(sname, "size%d", 0);
49 if (!hsize){
50 printf("Create [%s]\n", sname);
51 hsize = new TH1D(sname, sname, 600, 0, 6000);
52 }
53 hsize->Fill(size);
54 }
55
56
57 TV1190Data *v1190 = dataContainer.GetEventData<TV1190Data>("TDC0");
58 if(v1190){
59
60 std::cout << "TDC measurements for V1190" << std::endl;
61 std::vector<TDCMeasurement>& measurements = v1190->GetMeasurements();
62 for(unsigned int i = 0; i < measurements.size(); i++){
64
65 std::cout << "Measurement: " << tdcmeas.GetMeasurement() << " for tdc/chan " <<
66 tdcmeas.GetTDCNumber() << "/"<< tdcmeas.GetChannel() << std::endl;
67
68 }
69
70 }
71
72 TV792Data *v792 = dataContainer.GetEventData<TV792Data>("ADC0");
73 if(v792 ){
74
75 v792->Print();
76
77 }
78
79 return true;
80 }
81
82
83 // Describe some other command line argument
84 void Usage(void){
85 std::cout << "\t-D option: this is a fun new option " << std::endl;
86 }
87
88 // Define some other command line argument
89 bool CheckOption(std::string option){
90 const char* arg = option.c_str();
91 if (strncmp(arg,"-D",2)==0){
92 std::cout << arg+2 << std::endl;
93 std::cout << "I'm happy with this flag!" << std::endl;
94 return true;
95 }
96
97 return false;
98 }
99
100};
101
102
103
104
105
106
107int main(int argc, char *argv[])
108{
109
111 return MyTestLoop::Get().ExecuteLoop(argc, argv);
112
113}
114
int main(int argc, char *argv[])
bool ProcessMidasEvent(TDataContainer &dataContainer)
void Usage(void)
void EndRun(int transition, int run, int time)
virtual ~MyTestLoop()
void BeginRun(int transition, int run, int time)
bool CheckOption(std::string option)
uint32_t GetMeasurement() const
Get the TDC measurement.
virtual TObject * Get(const char *namecycle)
virtual Bool_t cd(const char *path=0)
static TRootanaEventLoop & Get(void)
int GetCurrentRunNumber() const
Current Run Number.
int ExecuteLoop(int argc, char *argv[])
Method to actually process the Midas information, either as file or online.
static void CreateSingleton()
TDirectory * fOnlineHistDir
TDirectory for online histograms.
std::vector< TDCMeasurement > & GetMeasurements()
Get the Vector of TDC Measurements.
Class for storing data from CAEN V792 module.
Definition TV792Data.hxx:53
void Print()
Print the bank contents in a structured way.
Definition TV792Data.cxx:32