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