]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/EMCALLEDda.cxx
Corrected bug, the histogram vector was not cleared properly
[u/mrichter/AliRoot.git] / EMCAL / EMCALLEDda.cxx
1 /*
2   EMCAL DA for online calibration: for LED studies
3   
4   Contact: silvermy@ornl.gov
5   Run Type: PHYSICS or STANDALONE
6   DA Type: MON 
7   Number of events needed: continously accumulating for all runs, rate ~0.1-1 Hz
8   Input Files: argument list
9   Output Files: RESULT_FILE=EMCALLED.root, to be exported to the DAQ FXS
10   fileId:  FILE_ID=EMCALLED    
11   Trigger types used: CALIBRATION_EVENT (temporarily also PHYSICS_EVENT to start with)
12   [When we have real data files later, we should only use CALIBRATION_EVENT]
13 */
14 /*
15   This process reads RAW data from the files provided as command line arguments
16   and save results (class itself) in a file (named from RESULT_FILE define - 
17   see below).
18 */
19
20 #define RESULT_FILE  "EMCALLED.root"
21 #define FILE_ID "EMCALLED"
22 #define AliDebugLevel() -1
23 #define FILE_PEDClassName "emcCalibPedestal"
24 #define FILE_SIGClassName "emcCalibSignal"
25 const int kNRCU = 4;
26 /* LOCAL_DEBUG is used to bypass daq* calls that do not work locally */
27 //#define LOCAL_DEBUG 1 // comment out to run normally                                                            
28 extern "C" {
29 #include <daqDA.h>
30 }
31 #include "event.h" /* in $DATE_COMMON_DEFS/; includes definition of event types */
32 #include "monitor.h" /* in $DATE_MONITOR_DIR/; monitor* interfaces */
33
34 #include "stdio.h"
35 #include "stdlib.h"
36
37 // ROOT includes
38 #include <TFile.h> 
39 #include <TROOT.h> 
40 #include <TPluginManager.h> 
41 #include <TSystem.h> 
42
43 //
44 //AliRoot includes
45 //
46 #include "AliRawReader.h"
47 #include "AliRawReaderDate.h"
48 #include "AliRawEventHeaderBase.h"
49 #include "AliCaloRawStream.h"
50 #include "AliCaloAltroMapping.h"
51 #include "AliLog.h"
52
53 //
54 // EMC calibration-helper algorithm includes
55 //
56 #include "AliCaloCalibPedestal.h"
57 #include "AliCaloCalibSignal.h"
58
59 /*
60   Main routine, EMC signal detector algorithm 
61   Arguments: list of DATE raw data files
62 */
63
64 int main(int argc, char **argv) {
65
66   AliLog::SetClassDebugLevel("AliCaloRawStream",-5);
67   AliLog::SetClassDebugLevel("AliRawReaderDate",-5);
68   AliLog::SetModuleDebugLevel("RAW",-5);
69
70   if (argc<2) {
71     printf("Wrong number of arguments\n");
72     return -1;  
73   }
74
75   /* magic line - for TStreamerInfo */
76   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
77                                         "*",
78                                         "TStreamerInfo",
79                                         "RIO",
80                                         "TStreamerInfo()"); 
81
82   int i, status;
83
84   /* log start of process */
85   printf("EMCAL DA started - %s\n",__FILE__);
86
87   /* declare monitoring program */
88   status=monitorDeclareMp( __FILE__ );
89   if (status!=0) {
90     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
91     return -1;
92   }
93   /* define wait event timeout - 1s max */
94   monitorSetNowait();
95   monitorSetNoWaitNetworkTimeout(1000);
96
97   /* Retrieve mapping files from DAQ DB */ 
98  const char* mapFiles[kNRCU] = {"RCU0A.data","RCU1A.data","RCU0C.data","RCU1C.data"};
99
100   for(Int_t iFile=0; iFile<kNRCU; iFile++) {
101     int failed = daqDA_DB_getFile(mapFiles[iFile], mapFiles[iFile]);
102     if(failed) { 
103       printf("Cannot retrieve file %d : %s from DAQ DB. Exit now..\n",
104              iFile, mapFiles[iFile]);
105 #ifdef LOCAL_DEBUG
106 #else
107       return -1;
108 #endif
109     }
110   }
111   
112   /* Open mapping files */
113   AliCaloAltroMapping *mapping[kNRCU];
114   TString path = "./";
115   path += "RCU";
116   TString path2;
117 TString side[] = {"A","C"};//+ and - pseudorapidity supermodules
118  for(Int_t j = 0; j < 2; j++){
119   for(Int_t i = 0; i < 2; i++) {
120     path2 = path;
121     path2 += i;
122     path2 += side[j]; 
123     path2 += ".data";
124     mapping[i] = new AliCaloAltroMapping(path2.Data());
125   }
126  }
127
128   /* set up our analysis classes */  
129   AliCaloCalibPedestal * calibPedestal = new AliCaloCalibPedestal(AliCaloCalibPedestal::kEmCal); 
130   calibPedestal->SetAltroMapping( mapping );
131   AliCaloCalibSignal * calibSignal = new AliCaloCalibSignal(AliCaloCalibSignal::kEmCal); 
132   calibSignal->SetAltroMapping( mapping );
133
134   AliRawReader *rawReader = NULL;
135   int nevents=0;
136
137   /* loop over RAW data files */
138   for ( i=1; i<argc; i++ ) {
139
140     /* define data source : this is argument i */
141     printf("Processing file %s\n", argv[i]);
142     status=monitorSetDataSource( argv[i] );
143     if (status!=0) {
144       printf("monitorSetDataSource() failed. Error=%s. Exiting ...\n", monitorDecodeError(status));
145       return -1;
146     }
147
148     /* read until EOF */
149     struct eventHeaderStruct *event;
150     eventTypeType eventT;
151
152     for ( ; ; ) { // infinite loop
153
154       /* check shutdown condition */
155       if (daqDA_checkShutdown()) {break;}
156
157       /* get next event (blocking call until timeout) */
158       status=monitorGetEventDynamic((void **)&event);
159       if (status==MON_ERR_EOF) {
160         printf ("End of File %d (%s) detected\n", i, argv[i]);
161         break; /* end of monitoring file has been reached */
162       }
163       if (status!=0) {
164         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
165         break;
166       }
167
168       /* retry if got no event */
169       if (event==NULL) {
170         continue;
171       }
172       eventT = event->eventType; /* just convenient shorthand */
173
174       /* skip start/end of run events */
175       if ( (eventT != physicsEvent) && (eventT != calibrationEvent) ) {
176         continue;
177       }
178
179       nevents++; // count how many acceptable events we have
180
181       //  Signal calibration
182       rawReader = new AliRawReaderDate((void*)event);
183       calibSignal->ProcessEvent(rawReader);
184       rawReader->Reset();
185       calibPedestal->ProcessEvent(rawReader);
186       delete rawReader;
187
188       /* free resources */
189       free(event);    
190
191     } //until EOF
192   } // loop over files
193
194   //
195   // write class to rootfile
196   //
197
198   printf ("%d physics/calibration events processed.\n",nevents);
199
200   TFile f(RESULT_FILE, "recreate");
201   if (!f.IsZombie()) { 
202     f.cd();
203     calibPedestal->Write(FILE_PEDClassName);
204     calibSignal->Write(FILE_SIGClassName);
205     f.Close();
206     printf("Objects saved to file \"%s\" as \"%s\" and \"%s\".\n", 
207            RESULT_FILE, FILE_PEDClassName, FILE_SIGClassName); 
208   } 
209   else {
210     printf("Could not save the object to file \"%s\".\n", 
211            RESULT_FILE);
212   }
213
214   //
215   // closing down; see if we can delete our analysis helper(s) also
216   //
217   delete calibPedestal;
218   delete calibSignal;
219   for(Int_t iFile=0; iFile<kNRCU; iFile++) {
220     if (mapping[iFile]) delete mapping[iFile];
221   }
222
223   /* store the result file on FES */
224 #ifdef LOCAL_DEBUG
225 #else
226   status = daqDA_FES_storeFile(RESULT_FILE, FILE_ID);
227 #endif
228
229   return status;
230 }