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