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