Two histograms added (expert mode):
[u/mrichter/AliRoot.git] / ZDC / ZDCPEDESTALda.cxx
CommitLineData
f3ca8792 1/*
2
039a08c3 3This program reads the DAQ data files passed as argument using the monitoring library.
4
5It computes the average event size and populates local "./result.txt" file with the
6result.
7
8The program reports about its processing progress.
f3ca8792 9
10Messages on stdout are exported to DAQ log system.
11
b3b2b9dd 12DA for ZDC standalone pedestal runs
ea628de7 13
b3b2b9dd 14Contact: Chiara.Oppedisano@to.infn.it
442e1b18 15Link:
b3b2b9dd 16Run Type: STANDALONE_PEDESTAL_RUN
442e1b18 17DA Type: LDC
b3b2b9dd 18Number of events needed: no constraint (tipically ~10^3)
442e1b18 19Input Files:
20Output Files: ZDCPedestal.dat, ZDCChMapping.dat
b3b2b9dd 21Trigger Types Used: Standalone Trigger
f3ca8792 22
23*/
218f916a 24#define PEDDATA_FILE "ZDCPedestal.dat"
25#define MAPDATA_FILE "ZDCChMapping.dat"
f3ca8792 26
27#include <stdio.h>
039a08c3 28#include <stdlib.h>
f3ca8792 29#include <Riostream.h>
30
31// DATE
f3ca8792 32#include <event.h>
33#include <monitor.h>
039a08c3 34#include <daqDA.h>
f3ca8792 35
36//ROOT
37#include <TRandom.h>
38#include <TH1F.h>
39#include <TH2F.h>
40#include <TProfile.h>
41#include <TF1.h>
42#include <TFile.h>
442e1b18 43#include <TFitter.h>
f3ca8792 44
45//AliRoot
46#include <AliRawReaderDate.h>
442e1b18 47#include <AliRawEventHeaderBase.h>
f3ca8792 48#include <AliZDCRawStream.h>
49
50
51/* Main routine
039a08c3 52 Arguments: list of DATE raw data files
f3ca8792 53*/
54int main(int argc, char **argv) {
442e1b18 55
56 TFitter *minuitFit = new TFitter(4);
57 TVirtualFitter::SetFitter(minuitFit);
f3ca8792 58
039a08c3 59 int status = 0;
60
61 /* log start of process */
63604324 62 printf("\n ZDC PEDESTAL program started\n");
039a08c3 63
64 /* check that we got some arguments = list of files */
65 if (argc<2) {
66 printf("Wrong number of arguments\n");
67 return -1;
68 }
69
f3ca8792 70 // --- Histograms for ADC pedestals
b3b2b9dd 71 // [22 signal channels + 2 reference PTMs] x 2 gain chains
f3ca8792 72 //
039a08c3 73 TH1F::AddDirectory(0);
b3b2b9dd 74 int const kNChannels = 24;
75 TH1F *hPedhg[kNChannels], *hPedOutOfTimehg[kNChannels];
76 TH2F *hPedCorrhg[kNChannels];
77 TH1F *hPedlg[kNChannels], *hPedOutOfTimelg[kNChannels];
78 TH2F *hPedCorrlg[kNChannels];
f3ca8792 79 //
b3b2b9dd 80 char namhist1hg[50], namhist2hg[50], namhist3hg[50];
81 char namhist1lg[50], namhist2lg[50], namhist3lg[50];
82 for(Int_t j=0; j<kNChannels; j++){
83 if(j<5){ // ZN1
84 sprintf(namhist1hg,"PedZN1hg_%d",j);
85 sprintf(namhist2hg,"PedZN1hgOutOfTime_%d",j);
86 sprintf(namhist3hg,"PedCorrZN1hg_%d",j);
87 //
88 sprintf(namhist1lg,"PedZN1lg_%d",j);
89 sprintf(namhist2lg,"PedZN1lgOutOfTime_%d",j);
90 sprintf(namhist3lg,"PedCorrZN1lg_%d",j);
91 }
92 else if(j>=5 && j<10){ // ZP1
93 sprintf(namhist1hg,"PedZP1hg_%d",j-5);
94 sprintf(namhist2hg,"PedZP1hgOutOfTime_%d",j-5);
95 sprintf(namhist3hg,"PedCorrZP1hg_%d",j-5);
96 //
97 sprintf(namhist1lg,"PedZP1lg_%d",j-5);
98 sprintf(namhist2lg,"PedZP1lgOutOfTime_%d",j-5);
99 sprintf(namhist3lg,"PedCorrZP1lg_%d",j-5);
f3ca8792 100 }
b3b2b9dd 101 else if(j>=10 && j<12){ // ZEM
102 sprintf(namhist1hg,"PedZEMhg_%d",j-10);
103 sprintf(namhist2hg,"PedZEMhgOutOfTime_%d",j-10);
104 sprintf(namhist3hg,"PedCorrZEMhg_%d",j-10);
105 //
106 sprintf(namhist1lg,"PedZEMlg_%d",j-10);
107 sprintf(namhist2lg,"PedZEMlgOutOfTime_%d",j-10);
108 sprintf(namhist3lg,"PedCorrZEMlg_%d",j-10);
f3ca8792 109 }
b3b2b9dd 110 else if(j>=12 && j<17){ // ZN2
111 sprintf(namhist1hg,"PedZN2hg_%d",j-12);
112 sprintf(namhist2hg,"PedZN2hgOutOfTime_%d",j-12);
113 sprintf(namhist3hg,"PedCorrZN2hg_%d",j-12);
114 //
115 sprintf(namhist1lg,"PedZN2lg_%d",j-12);
116 sprintf(namhist2lg,"PedZN2lgOutOfTime_%d",j-12);
117 sprintf(namhist3lg,"PedCorrZN2lg_%d",j-12);
f3ca8792 118 }
b3b2b9dd 119 else if(j>=17 && j<22){ // ZP2
120 sprintf(namhist1hg,"PedZP2hg_%d",j-17);
121 sprintf(namhist2hg,"PedZP2hgOutOfTime_%d",j-17);
122 sprintf(namhist3hg,"PedCorrZP2hg_%d",j-17);
123 //
124 sprintf(namhist1lg,"PedZP2lg_%d",j-17);
125 sprintf(namhist2lg,"PedZP2lgOutOfTime_%d",j-17);
126 sprintf(namhist3lg,"PedCorrZP2lg_%d",j-17);
f3ca8792 127 }
b3b2b9dd 128 else if(j>=22 && j<24){ //Reference PMs
129 sprintf(namhist1hg,"PedRefhg_%d",j-22);
130 sprintf(namhist2hg,"PedRefhgOutOfTime_%d",j-22);
131 sprintf(namhist3hg,"PedCorrRefhg_%d",j-22);
132 //
133 sprintf(namhist1lg,"PedReflg_%d",j-22);
134 sprintf(namhist2lg,"PedReflgOutOfTime_%d",j-22);
135 sprintf(namhist3lg,"PedCorrReflg_%d",j-22);
f3ca8792 136 }
b3b2b9dd 137 // --- High gain chain histos
442e1b18 138 hPedhg[j] = new TH1F(namhist1hg, namhist1hg, 200,0., 200.);
139 hPedOutOfTimehg[j] = new TH1F(namhist2hg, namhist2hg, 200,0., 200.);
b3b2b9dd 140 hPedCorrhg[j] = new TH2F(namhist3hg,namhist3hg,100,0.,200.,100,0.,200.);
141 // --- Low gain chain histos
142 hPedlg[j] = new TH1F(namhist1lg, namhist1lg, 100,0., 600.);
143 hPedOutOfTimelg[j] = new TH1F(namhist2lg, namhist2lg, 100,0., 600.);
144 hPedCorrlg[j] = new TH2F(namhist3lg,namhist3lg,100,0.,600.,100,0.,600.);
f3ca8792 145 }
146
f3ca8792 147
148 /* open result file */
149 FILE *fp=NULL;
150 fp=fopen("./result.txt","a");
151 if (fp==NULL) {
152 printf("Failed to open file\n");
153 return -1;
154 }
155
442e1b18 156 FILE *mapFile4Shuttle;
f3ca8792 157
039a08c3 158 /* report progress */
159 daqDA_progressReport(10);
f3ca8792 160
f3ca8792 161
f3ca8792 162 /* init some counters */
163 int nevents_physics=0;
164 int nevents_total=0;
165
039a08c3 166 /* read the data files */
167 int n;
442e1b18 168 for(n=1;n<argc;n++){
039a08c3 169
170 status=monitorSetDataSource( argv[n] );
f3ca8792 171 if (status!=0) {
039a08c3 172 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
173 return -1;
f3ca8792 174 }
175
039a08c3 176 /* report progress */
177 /* in this example, indexed on the number of files */
178 daqDA_progressReport(10+80*n/argc);
179
180 /* read the file */
181 for(;;) {
182 struct eventHeaderStruct *event;
183 eventTypeType eventT;
184
185 /* get next event */
186 status=monitorGetEventDynamic((void **)&event);
442e1b18 187 if(status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
188 if(status!=0) {
039a08c3 189 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
190 return -1;
191 }
192
193 /* retry if got no event */
442e1b18 194 if(event==NULL) {
039a08c3 195 break;
196 }
442e1b18 197
198 // Initalize raw-data reading and decoding
199 AliRawReader *reader = new AliRawReaderDate((void*)event);
200 reader->Select("ZDC");
201 // --- Reading event header
202 //UInt_t evtype = reader->GetType();
203 //printf("\n\t ZDCPEDESTALda -> ev. type %d\n",evtype);
204 //printf("\t ZDCPEDESTALda -> run # %d\n",reader->GetRunNumber());
205 //
206 AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
207
f3ca8792 208
039a08c3 209 /* use event - here, just write event id to result file */
210 eventT=event->eventType;
442e1b18 211
212 Int_t ich=0, adcMod[48], adcCh[48], sigCode[48], det[48], sec[48];
213 if(eventT==START_OF_DATA){
214
215 if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
216 else{
217 while(rawStreamZDC->Next()){
218 if(rawStreamZDC->IsChMapping()){
219 adcMod[ich] = rawStreamZDC->GetADCModFromMap(ich);
220 adcCh[ich] = rawStreamZDC->GetADCChFromMap(ich);
221 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
222 det[ich] = rawStreamZDC->GetDetectorFromMap(ich);
223 sec[ich] = rawStreamZDC->GetTowerFromMap(ich);
224 ich++;
225 }
226 }
227 }
228 // --------------------------------------------------------
229 // --- Writing ascii data file for the Shuttle preprocessor
218f916a 230 mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
442e1b18 231 for(Int_t i=0; i<ich; i++){
232 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",i,
233 adcMod[i],adcCh[i],sigCode[i],det[i],sec[i]);
234 //
235 //printf("ZDCPEDESTALDA.cxx -> ch.%d mod %d, ch %d, code %d det %d, sec %d\n",
236 // i,adcMod[i],adcCh[i],sigCode[i],det[i],sec[i]);
237 }
238 fclose(mapFile4Shuttle);
239 }
240
039a08c3 241 if(eventT==PHYSICS_EVENT){
242 fprintf(fp,"Run #%lu, event size: %lu, BC:%u, Orbit:%u, Period:%u\n",
243 (unsigned long)event->eventRunNb,
244 (unsigned long)event->eventSize,
245 EVENT_ID_GET_BUNCH_CROSSING(event->eventId),
246 EVENT_ID_GET_ORBIT(event->eventId),
247 EVENT_ID_GET_PERIOD(event->eventId));
442e1b18 248
249 // --- Reading data header
250 reader->ReadHeader();
039a08c3 251 const AliRawDataHeader* header = reader->GetDataHeader();
442e1b18 252 if(header){
b3b2b9dd 253 UChar_t message = header->GetAttributes();
442e1b18 254 if(message & 0x20){ // PEDESTAL RUN
255 //printf("\t STANDALONE_PEDESTAL RUN raw data found\n");
256 }
257 else{
258 printf("\t NO STANDALONE_PEDESTAL RUN raw data found\n");
259 return -1;
260 }
261 }
262 else{
263 printf("\t ATTENTION! No Raw Data Header found!!!\n");
264 return -1;
265 }
266
267 if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
268 //
269 // ----- Setting ch. mapping -----
270 for(Int_t jk=0; jk<48; jk++){
271 rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
272 rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
273 rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
274 rawStreamZDC->SetMapDet(jk, det[jk]);
275 rawStreamZDC->SetMapTow(jk, sec[jk]);
276 }
277 //
278 Int_t iraw=0;
039a08c3 279 Int_t RawADChg[kNChannels], RawADCoothg[kNChannels];
280 Int_t RawADClg[kNChannels], RawADCootlg[kNChannels];
281 for(Int_t j=0; j<kNChannels; j++){
282 RawADChg[j]=0; RawADCoothg[j]=0;
283 RawADClg[j]=0; RawADCootlg[j]=0;
284 }
442e1b18 285 //
039a08c3 286 while(rawStreamZDC->Next()){
442e1b18 287 Int_t index=-1;
288 if(rawStreamZDC->IsADCDataWord()){
b3b2b9dd 289 if(rawStreamZDC->GetSector(1)!=5){ // Physics signals
c35ed519 290 if(rawStreamZDC->GetSector(0)==1) index = rawStreamZDC->GetSector(1); // *** ZNC
291 else if(rawStreamZDC->GetSector(0)==2) index = rawStreamZDC->GetSector(1)+5; // *** ZPC
b3b2b9dd 292 else if(rawStreamZDC->GetSector(0)==3) index = rawStreamZDC->GetSector(1)+9; // *** ZEM
c35ed519 293 else if(rawStreamZDC->GetSector(0)==4) index = rawStreamZDC->GetSector(1)+12; // *** ZNA
294 else if(rawStreamZDC->GetSector(0)==5) index = rawStreamZDC->GetSector(1)+17; // *** ZPA
f3ca8792 295 }
b3b2b9dd 296 else{ // Reference PMs
297 index = (rawStreamZDC->GetSector(0)-1)/3+22;
f3ca8792 298 }
b3b2b9dd 299 //
442e1b18 300 /*printf("\t iraw %d index %d det %d quad %d res %d ADC %d\n", iraw, index,
b3b2b9dd 301 rawStreamZDC->GetSector(0), rawStreamZDC->GetSector(1),
302 rawStreamZDC->GetADCGain(), rawStreamZDC->GetADCValue());
303 */
442e1b18 304 //
305 if(iraw<2*kNChannels){ // --- In-time pedestals (1st 48 raw data)
b3b2b9dd 306 if(rawStreamZDC->GetADCGain()==0){
307 hPedhg[index]->Fill(rawStreamZDC->GetADCValue());
308 RawADChg[index] = rawStreamZDC->GetADCValue();
442e1b18 309 //
310 //printf("\t filling histo hPedhg[%d]\n",index);
b3b2b9dd 311 }
312 else{
313 hPedlg[index]->Fill(rawStreamZDC->GetADCValue());
314 RawADClg[index] = rawStreamZDC->GetADCValue();
442e1b18 315 //
316 //printf("\t filling histo hPedlg[%d]\n",index);
b3b2b9dd 317 }
442e1b18 318 }
319 else{ // --- Out-of-time pedestals
320 if(rawStreamZDC->GetADCGain()==0){
321 hPedOutOfTimehg[index]->Fill(rawStreamZDC->GetADCValue());
322 RawADCoothg[index] = rawStreamZDC->GetADCValue();
323 //
324 //printf("\t filling histo hPedOutOfTimehg[%d]\n",index);
039a08c3 325 }
442e1b18 326 else{
327 hPedOutOfTimelg[index]->Fill(rawStreamZDC->GetADCValue());
328 RawADCootlg[index] = rawStreamZDC->GetADCValue();
329 //
330 //printf("\t filling histo hPedOutOfTimelg[%d]\n",index);
039a08c3 331 }
442e1b18 332 }
333 iraw++;
334 }//IsADCDataWord()
335 //
336 if(iraw == 4*kNChannels){ // Last ADC channel -> Filling correlation histos
337 for(Int_t k=0; k<kNChannels; k++){
039a08c3 338 hPedCorrhg[k]->Fill(RawADCoothg[k], RawADChg[k]);
339 hPedCorrlg[k]->Fill(RawADCootlg[k], RawADClg[k]);
442e1b18 340 }
341 }
039a08c3 342 }
343 //
344 nevents_physics++;
345 //
346 delete reader;
347 delete rawStreamZDC;
348
349 }//(if PHYSICS_EVENT)
350 nevents_total++;
351
352 /* free resources */
353 free(event);
f3ca8792 354
f3ca8792 355 }
039a08c3 356 }
f3ca8792 357
358 /* Analysis of the histograms */
359 //
360 FILE *fileShuttle;
218f916a 361 fileShuttle = fopen(PEDDATA_FILE,"w");
f3ca8792 362 //
b3b2b9dd 363 Float_t MeanPed[2*kNChannels], MeanPedWidth[2*kNChannels],
63604324 364 MeanPedOOT[2*kNChannels], MeanPedWidthOOT[2*kNChannels];
b3b2b9dd 365 // --- In-time pedestals
366 TF1 *ADCfunchg[kNChannels];
367 for(Int_t i=0; i<kNChannels; i++){
368 hPedhg[i]->Fit("gaus","Q");
369 ADCfunchg[i] = hPedhg[i]->GetFunction("gaus");
442e1b18 370 MeanPed[i] = (Double_t) ADCfunchg[i]->GetParameter(1);
371 MeanPedWidth[i] = (Double_t) ADCfunchg[i]->GetParameter(2);
218f916a 372 fprintf(fileShuttle,"\t%f\t%f\n",MeanPed[i],MeanPedWidth[i]);
ea628de7 373 //printf("\t MeanPed[%d] = %f\n",i, MeanPed[i]);
b3b2b9dd 374 }
375 TF1 *ADCfunclg[kNChannels];
376 for(Int_t i=0; i<kNChannels; i++){
377 hPedlg[i]->Fit("gaus","Q");
378 ADCfunclg[i] = hPedlg[i]->GetFunction("gaus");
442e1b18 379 MeanPed[i+kNChannels] = (Double_t) ADCfunclg[i]->GetParameter(1);
380 MeanPedWidth[i+kNChannels] = (Double_t) ADCfunclg[i]->GetParameter(2);
218f916a 381 fprintf(fileShuttle,"\t%f\t%f\n",MeanPed[i+kNChannels],MeanPedWidth[i+kNChannels]);
ea628de7 382 //printf("\t MeanPed[%d] = %f\n",i+kNChannels, MeanPed[i+kNChannels]);
f3ca8792 383 }
384 // --- Out-of-time pedestals
b3b2b9dd 385 TF1 *ADCootfunchg[kNChannels];
386 for(Int_t i=0; i<kNChannels; i++){
387 hPedOutOfTimehg[i]->Fit("gaus","Q");
388 ADCootfunchg[i] = hPedOutOfTimehg[i]->GetFunction("gaus");
442e1b18 389 MeanPedOOT[i] = (Double_t) ADCootfunchg[i]->GetParameter(1);
390 MeanPedWidthOOT[i] = (Double_t) ADCootfunchg[i]->GetParameter(2);
218f916a 391 fprintf(fileShuttle,"\t%f\t%f\n",MeanPedOOT[i],MeanPedWidthOOT[i]);
ea628de7 392 //printf("\t MeanPedOOT[%d] = %f\n",i, MeanPedOOT[i]);
b3b2b9dd 393 }
394 TF1 *ADCootfunclg[kNChannels];
395 for(Int_t i=0; i<kNChannels; i++){
396 hPedOutOfTimelg[i]->Fit("gaus","Q");
397 ADCootfunclg[i] = hPedOutOfTimelg[i]->GetFunction("gaus");
442e1b18 398 MeanPedOOT[i+kNChannels] = (Double_t) ADCootfunclg[i]->GetParameter(1);
399 MeanPedWidthOOT[i+kNChannels] = (Double_t) ADCootfunclg[i]->GetParameter(2);
218f916a 400 fprintf(fileShuttle,"\t%f\t%f\n",MeanPedOOT[i+kNChannels],MeanPedWidthOOT[i+kNChannels]);
ea628de7 401 //printf("\t MeanPedOOT[%d] = %f\n",i+kNChannels, MeanPedOOT[i+kNChannels]);
f3ca8792 402 }
63604324 403
404 // ***************************************************
405 // Unless we have a narrow correlation to fit we
406 // don't fit and store in-time vs. out-of-time
407 // histograms -> mean pedstal subtracted!!!!!!
408 // ***************************************************
b3b2b9dd 409 // --- Correlations
63604324 410/*
411 Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels];
b3b2b9dd 412 TProfile *hPedCorrProfhg[kNChannels], *hPedCorrProflg[kNChannels];
413 TF1 *ffunchg[kNChannels], *ffunclg[kNChannels];
f3ca8792 414 char namhist4[50];
b3b2b9dd 415 for(int i=0;i<kNChannels;i++) {
416 sprintf(namhist4,"ADCHRvsOOT%d_Prof",i);
417 hPedCorrProfhg[i] = hPedCorrhg[i]->ProfileX(namhist4,-1,-1,"S");
418 hPedCorrProfhg[i]->SetName(namhist4);
419 hPedCorrProfhg[i]->Fit("pol1","Q");
420 ffunchg[i] = hPedCorrProfhg[i]->GetFunction("pol1");
442e1b18 421 CorrCoeff0[i] = (Double_t) ffunchg[i]->GetParameter(0);
422 CorrCoeff1[i] = (Double_t) ffunchg[i]->GetParameter(1);
c35ed519 423 fprintf(fileShuttle,"\t%d\t%f\t%f\n",i,CorrCoeff0[i],CorrCoeff1[i]);
ea628de7 424 //printf("\t CorrCoeff0[%d] = %f, CorrCoeff1[%d] = %f\n",i, CorrCoeff0[i], i, CorrCoeff1[i]);
b3b2b9dd 425 }
426 for(int i=0;i<kNChannels;i++) {
427 sprintf(namhist4,"ADCLRvsOOT%d_Prof",i);
428 hPedCorrProflg[i] = hPedCorrlg[i]->ProfileX(namhist4,-1,-1,"S");
429 hPedCorrProflg[i]->SetName(namhist4);
430 hPedCorrProflg[i]->Fit("pol1","Q");
431 ffunclg[i] = hPedCorrProflg[i]->GetFunction("pol1");
442e1b18 432 CorrCoeff0[i+kNChannels] = (Double_t) ffunclg[i]->GetParameter(0);
433 CorrCoeff1[i+kNChannels] = (Double_t) ffunclg[i]->GetParameter(1);
c35ed519 434 fprintf(fileShuttle,"\t%d\t%f\t%f\n",i+kNChannels,CorrCoeff0[i+kNChannels],CorrCoeff1[i+kNChannels]);
ea628de7 435 //printf("\t CorrCoeff0[%d] = %f, CorrCoeff1[%d] = %f\n",
436 // i+kNChannels, CorrCoeff0[i+kNChannels], i+kNChannels, CorrCoeff1[i+kNChannels]);
f3ca8792 437 }
63604324 438*/
f3ca8792 439 //
440 fclose(fileShuttle);
b3b2b9dd 441 //
039a08c3 442 for(Int_t j=0; j<kNChannels; j++){
443 delete hPedhg[j];
444 delete hPedOutOfTimehg[j];
445 delete hPedCorrhg[j];
446 delete hPedlg[j];
447 delete hPedOutOfTimelg[j];
448 delete hPedCorrlg[j];
449 }
f3ca8792 450
442e1b18 451 //delete minuitFit;
452 TVirtualFitter::SetFitter(0);
453
f3ca8792 454 /* write report */
455 fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
456
457 /* close result file */
458 fclose(fp);
039a08c3 459
460 /* report progress */
461 daqDA_progressReport(90);
462
442e1b18 463 /* store the result files on FES */
63604324 464 // [1] File with mapping
218f916a 465 status = daqDA_FES_storeFile(MAPDATA_FILE,MAPDATA_FILE);
442e1b18 466 if(status){
63604324 467 printf("Failed to export file %d to DAQ FES\n",status);
442e1b18 468 return -1;
469 }
63604324 470 // [2] File with pedestal data
218f916a 471 status = daqDA_FES_storeFile(PEDDATA_FILE,PEDDATA_FILE);
039a08c3 472 if(status){
63604324 473 printf("Failed to export file %d to DAQ FES\n",status);
474 return -1;
475 }
476
477 /* store the result files on DB */
218f916a 478 status = daqDA_DB_storeFile(PEDDATA_FILE,PEDDATA_FILE);
63604324 479 if(status){
480 printf("Failed to export file %d to DAQ DB\n",status);
039a08c3 481 return -1;
482 }
483
484 /* report progress */
485 daqDA_progressReport(100);
f3ca8792 486
f3ca8792 487 return status;
488}