]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/PMDda.cxx
added protection
[u/mrichter/AliRoot.git] / PMD / PMDda.cxx
1 /*
2
3 DAcase2.c
4
5 This program connects to the DAQ data source passed as argument
6 and populates local "./result.txt" file with the ids of events received
7 during the run.
8
9 The program exits when being asked to shut down (daqDA_checkshutdown)
10 or End of Run event.
11
12 Messages on stdout are exported to DAQ log system.
13
14 contact: alice-datesupport@cern.ch
15
16 */
17 extern "C" {
18 #include <daqDA.h>
19 }
20
21 #include "event.h"
22 #include "monitor.h"
23 //#include "daqDA.h"
24
25 #include <Riostream.h>
26 #include <stdio.h>
27 #include <stdlib.h>
28
29 //AliRoot
30 #include "AliRawReaderDate.h"
31 #include "AliPMDCalibPedestal.h"
32 #include "AliPMDCalibGain.h"
33
34 //ROOT
35 #include "TFile.h"
36 #include "TH1F.h"
37 #include "TBenchmark.h"
38 #include "TTree.h"
39 #include "TROOT.h"
40 #include "TPluginManager.h"
41
42
43
44 /* Main routine
45       Arguments: 
46       1- monitoring data source
47 */
48 int main(int argc, char **argv) {
49   
50     /* magic line from Rene */
51     gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
52                                           "*",
53                                           "TStreamerInfo",
54                                           "RIO",
55                                           "TStreamerInfo()");
56
57     
58     AliPMDCalibPedestal calibped;
59     AliPMDCalibGain calibgain;
60
61   TTree *ped  = new TTree("ped","PMD Pedestal tree");
62   TTree *gain = new TTree("gain","PMD Gain tree");
63
64   TH1F::AddDirectory(0);
65   
66       
67   // decoding the events
68   
69   int status;
70
71   if (argc!=2) {
72     printf("Wrong number of arguments\n");
73     return -1;
74   }
75
76   /* open result file */
77   FILE *fp=NULL;
78   fp=fopen("./result.txt","a");
79   if (fp==NULL) {
80     printf("Failed to open file\n");
81     return -1;
82   }
83
84   /* define data source : this is argument 1 */  
85   status=monitorSetDataSource( argv[1] );
86   if (status!=0) {
87     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
88     return -1;
89   }
90
91   /* declare monitoring program */
92   status=monitorDeclareMp( __FILE__ );
93   if (status!=0) {
94     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
95     return -1;
96   }
97
98   /* define wait event timeout - 1s max */
99   monitorSetNowait();
100   monitorSetNoWaitNetworkTimeout(1000);
101   
102   /* log start of process */
103   printf("DA example case2 monitoring program started\n");  
104
105   /* init some counters */
106   int nevents_physics=0;
107   int nevents_total=0;
108
109   struct eventHeaderStruct *event;
110   eventTypeType eventT;
111   Int_t iev=0;
112
113   /* main loop (infinite) */
114   for(;;) {
115     
116     /* check shutdown condition */
117     if (daqDA_checkShutdown()) {break;}
118     
119     /* get next event (blocking call until timeout) */
120     status=monitorGetEventDynamic((void **)&event);
121     if (status==MON_ERR_EOF) {
122       printf ("End of File detected\n");
123       break; /* end of monitoring file has been reached */
124     }
125     
126     if (status!=0) {
127       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
128       break;
129     }
130     
131     /* retry if got no event */
132     if (event==NULL) {
133       continue;
134     }
135
136     iev++; 
137
138    /* use event - here, just write event id to result file */
139     eventT=event->eventType;
140     switch (event->eventType){
141       
142       /* START OF RUN */
143     case START_OF_RUN:
144       break;
145       /* END START OF RUN */
146       
147     /* END OF RUN */
148     case END_OF_RUN:
149       break;
150       
151     case PHYSICS_EVENT:
152       printf(" event number = %i \n",iev);
153       AliRawReader *rawReader = new AliRawReaderDate((void*)event);
154       calibped.ProcessEvent(rawReader);
155
156       calibgain.ProcessEvent(rawReader);
157
158       delete rawReader;
159       rawReader = 0x0;
160
161     }
162        
163     /* free resources */
164     free(event);
165     
166     /* exit when last event received, no need to wait for TERM signal */
167     if (eventT==END_OF_RUN) {
168       printf("EOR event detected\n");
169       calibped.Analyse(ped);
170       calibgain.Analyse(gain);
171
172       break;
173     }
174   }
175   
176   //write the Run level file   
177   TFile * fileRun = new TFile ("outPMDdaRun.root","RECREATE"); 
178   TBenchmark *bench = new TBenchmark();
179   bench->Start("PMD");
180   bench->Stop("PMD");
181   bench->Print("PMD");
182   fileRun->Close();
183
184   /* write report */
185   fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
186
187
188   TFile * pedRun = new TFile ("pmd_ped.root","RECREATE"); 
189   ped->Write();
190   pedRun->Close();
191
192   TFile * gainRun = new TFile ("pmd_calib.root","RECREATE"); 
193   gain->Write();
194   gainRun->Close();
195
196
197
198   /* close result file */
199   fclose(fp);
200
201
202   return status;
203 }