]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ZDC/ZDCLASERda.cxx
Changed default values of cuts
[u/mrichter/AliRoot.git] / ZDC / ZDCLASERda.cxx
CommitLineData
ead118d8 1/*
2
3This program reads the DAQ data files passed as argument using the monitoring library.
4
ead118d8 5The program reports about its processing progress.
6
7Messages on stdout are exported to DAQ log system.
8
9DA for ZDC standalone pedestal runs
10
11Contact: Chiara.Oppedisano@to.infn.it
12Link:
13Run Type: STANDALONE_LASER_RUN
442e1b18 14DA Type: LDC
ead118d8 15Number of events needed: no constraint (tipically ~10^3)
19419b67 16Input Files: ZDCPedestal.dat
ead118d8 17Output Files: ZDCLaser.dat
18Trigger Types Used: Standalone Trigger
19
20*/
78beff0d 21#define PEDDATA_FILE "ZDCPedestal.dat"
218f916a 22#define MAPDATA_FILE "ZDCChMapping.dat"
65448375 23#define LASHISTO_FILE "ZDCLaserHisto.root"
218f916a 24#define LASDATA_FILE "ZDCLaserCalib.dat"
ead118d8 25
26#include <stdio.h>
27#include <stdlib.h>
28#include <Riostream.h>
29
30// DATE
31#include <event.h>
32#include <monitor.h>
33#include <daqDA.h>
34
35//ROOT
65448375 36#include <TROOT.h>
37#include <TPluginManager.h>
ead118d8 38#include <TH1F.h>
39#include <TF1.h>
40#include <TFile.h>
442e1b18 41#include <TFitter.h>
65448375 42#include "TMinuitMinimizer.h"
ead118d8 43
44//AliRoot
45#include <AliRawReaderDate.h>
442e1b18 46#include <AliRawEventHeaderBase.h>
ead118d8 47#include <AliZDCRawStream.h>
48
49
50/* Main routine
51 Arguments: list of DATE raw data files
52*/
53int main(int argc, char **argv) {
442e1b18 54
65448375 55 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
56 "*",
57 "TStreamerInfo",
58 "RIO",
59 "TStreamerInfo()");
60
61 TMinuitMinimizer m;
62 gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", "Minuit","TMinuitMinimizer",
63 "Minuit", "TMinuitMinimizer(const char *)");
64 TVirtualFitter::SetDefaultFitter("Minuit");
ead118d8 65
66 int status = 0;
9d0b658a 67 int const kNChannels = 24;
27afc0c8 68 int const kNScChannels = 32;
19419b67 69
70 Int_t ich=0;
71 Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
72 Int_t det[2*kNChannels], sec[2*kNChannels];
73 for(Int_t y=0; y<2*kNChannels; y++){
74 adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0;
75 }
ead118d8 76
77 /* log start of process */
9d0b658a 78 printf("\n ZDC LASER program started\n");
ead118d8 79
80 /* check that we got some arguments = list of files */
81 if (argc<2) {
82 printf("Wrong number of arguments\n");
83 return -1;
84 }
85
86 // --- Histograms for LASER runs
87 // 20 signal channels + 2 reference PTMs
88 //
89 TH1F::AddDirectory(0);
bd94dbdd 90 // --- Histos for reference PMTs (high gain chains)
19419b67 91 TH1F *hPMRefChg = new TH1F("hPMRefChg","hPMRefChg", 100,0.,1000.);
92 TH1F *hPMRefAhg = new TH1F("hPMRefAhg","hPMRefAhg", 100,0.,1000.);
9d0b658a 93 TH1F *hPMRefClg = new TH1F("hPMRefClg","hPMRefClg", 100,0.,4000.);
94 TH1F *hPMRefAlg = new TH1F("hPMRefAlg","hPMRefAlg", 100,0.,4000.);
78beff0d 95 //
9d0b658a 96 // --- Histos for detector PMTs
97 TH1F *hZNChg[5], *hZPChg[5], *hZNAhg[5], *hZPAhg[5], *hZEMhg[2];
98 TH1F *hZNClg[5], *hZPClg[5], *hZNAlg[5], *hZPAlg[5], *hZEMlg[2];
99 char hnamZNChg[20], hnamZPChg[20], hnamZNAhg[20], hnamZPAhg[20];
100 char hnamZNClg[20], hnamZPClg[20], hnamZNAlg[20], hnamZPAlg[20];
101 char hnamZEMhg[20], hnamZEMlg[20];
bd94dbdd 102 for(Int_t j=0; j<5; j++){
9d0b658a 103 sprintf(hnamZNChg,"ZNChg-tow%d",j);
104 sprintf(hnamZPChg,"ZPChg-tow%d",j);
105 sprintf(hnamZNAhg,"ZNAhg-tow%d",j);
106 sprintf(hnamZPAhg,"ZPAhg-tow%d",j);
107 //
19419b67 108 hZNChg[j] = new TH1F(hnamZNChg, hnamZNChg, 100, 0., 1000.);
109 hZPChg[j] = new TH1F(hnamZPChg, hnamZPChg, 100, 0., 1000.);
110 hZNAhg[j] = new TH1F(hnamZNAhg, hnamZNAhg, 100, 0., 1000.);
111 hZPAhg[j] = new TH1F(hnamZPAhg, hnamZPAhg, 100, 0., 1000.);
9d0b658a 112 //
113 sprintf(hnamZNClg,"ZNClg-tow%d",j);
114 sprintf(hnamZPClg,"ZPClg-tow%d",j);
115 sprintf(hnamZNAlg,"ZNAlg-tow%d",j);
116 sprintf(hnamZPAlg,"ZPAlg-tow%d",j);
bd94dbdd 117 //
9d0b658a 118 hZNClg[j] = new TH1F(hnamZNClg, hnamZNClg, 100, 0., 4000.);
119 hZPClg[j] = new TH1F(hnamZPClg, hnamZPClg, 100, 0., 4000.);
120 hZNAlg[j] = new TH1F(hnamZNAlg, hnamZNAlg, 100, 0., 4000.);
121 hZPAlg[j] = new TH1F(hnamZPAlg, hnamZPAlg, 100, 0., 4000.);
122 //
123 if(j<2){
124 sprintf(hnamZEMhg,"ZEM%dhg",j);
125 sprintf(hnamZEMlg,"ZEM%dlg",j);
126 //
19419b67 127 hZEMhg[j] = new TH1F(hnamZEMhg, hnamZEMhg, 100, 0., 1000.);
9d0b658a 128 hZEMlg[j] = new TH1F(hnamZEMlg, hnamZEMlg, 100, 0., 4000.);
129 }
bd94dbdd 130 }
ead118d8 131
132 /* open result file */
133 FILE *fp=NULL;
134 fp=fopen("./result.txt","a");
135 if (fp==NULL) {
136 printf("Failed to open file\n");
137 return -1;
138 }
9d0b658a 139 /* report progress */
140 daqDA_progressReport(10);
78beff0d 141
142 // *** To analyze LASER events you MUST have a pedestal data file!!!
9d0b658a 143 // *** -> check if a pedestal run has been analyzed
78beff0d 144 int read = 0;
9d0b658a 145 read = daqDA_DB_getFile(PEDDATA_FILE, PEDDATA_FILE);
78beff0d 146 if(read){
147 printf("\t ERROR!!! ZDCPedestal.dat file NOT FOUND in DAQ db!!!\n");
148 return -1;
149 }
150 else printf("\t ZDCPedestal.dat file retrieved from DAQ db\n");
151
152 FILE *filePed = fopen(PEDDATA_FILE,"r");
153 if (filePed==NULL) {
154 printf("\t ERROR!!! Can't open ZDCPedestal.dat file!!!\n");
155 return -1;
156 }
157
158 // 144 = 48 in-time + 48 out-of-time + 48 correlations
9d0b658a 159 Float_t readValues[2][6*kNChannels];
160 Float_t MeanPedhg[kNChannels], MeanPedlg[kNChannels];
161 Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels];
78beff0d 162 // ***************************************************
163 // Unless we have a narrow correlation to fit we
164 // don't fit and store in-time vs. out-of-time
165 // histograms -> mean pedstal subtracted!!!!!!
166 // ***************************************************
78beff0d 167 //
9d0b658a 168 for(int jj=0; jj<6*kNChannels; jj++){
78beff0d 169 for(int ii=0; ii<2; ii++){
170 fscanf(filePed,"%f",&readValues[ii][jj]);
171 }
9d0b658a 172 if(jj<kNChannels){
173 MeanPedhg[jj] = readValues[0][jj];
174 //printf("\t MeanPedhg[%d] = %1.1f\n",jj, MeanPedhg[jj]);
78beff0d 175 }
9d0b658a 176 else if(jj>=kNChannels && jj<2*kNChannels){
177 MeanPedlg[jj-kNChannels] = readValues[0][jj];
178 //printf("\t MeanPedlg[%d] = %1.1f\n",jj-kNChannels, MeanPedlg[jj-kNChannels]);
78beff0d 179 }
9d0b658a 180 else if(jj>4*kNChannels){
181 CorrCoeff0[jj-4*kNChannels] = readValues[0][jj];
182 CorrCoeff1[jj-4*kNChannels] = readValues[1][jj];;
78beff0d 183 }
78beff0d 184 }
9d0b658a 185
186 FILE *mapFile4Shuttle;
ead118d8 187
188 /* report progress */
9d0b658a 189 daqDA_progressReport(20);
ead118d8 190
191
192 /* init some counters */
193 int nevents_physics=0;
194 int nevents_total=0;
195
d3f256ca 196 struct eventHeaderStruct *event;
197 eventTypeType eventT;
198
ead118d8 199 /* read the data files */
200 int n;
9d0b658a 201 for(n=1;n<argc;n++) {
ead118d8 202
203 status=monitorSetDataSource( argv[n] );
204 if (status!=0) {
205 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
206 return -1;
207 }
208
209 /* report progress */
210 /* in this example, indexed on the number of files */
9d0b658a 211 daqDA_progressReport(20+70*n/argc);
ead118d8 212
213 /* read the file */
214 for(;;) {
ead118d8 215
216 /* get next event */
217 status=monitorGetEventDynamic((void **)&event);
218 if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
219 if (status!=0) {
220 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
221 return -1;
222 }
223
224 /* retry if got no event */
225 if (event==NULL) {
226 break;
227 }
228
442e1b18 229 // Initalize raw-data reading and decoding
230 AliRawReader *reader = new AliRawReaderDate((void*)event);
231 reader->Select("ZDC");
232 // --- Reading event header
78beff0d 233 //UInt_t evtype = reader->GetType();
234 //printf("\n\t ZDCLASERda -> ev. type %d\n",evtype);
235 //printf("\t ZDCLASERda -> run # %d\n",reader->GetRunNumber());
442e1b18 236 //
237 AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
238
239
240 /* use event - here, just write event id to result file */
241 eventT=event->eventType;
19419b67 242
27afc0c8 243 Int_t iScCh=0;
244 Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels];
245 Int_t scDet[kNScChannels], scSec[kNScChannels];
246 for(Int_t y=0; y<kNScChannels; y++){
247 scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0;
248 }
249 //
250 Int_t modNum=-1, modType=-1;
251
442e1b18 252 if(eventT==START_OF_DATA){
253
27afc0c8 254 rawStreamZDC->SetSODReading(kTRUE);
255
9d0b658a 256 // --------------------------------------------------------
257 // --- Writing ascii data file for the Shuttle preprocessor
258 mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
442e1b18 259 if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
260 else{
27afc0c8 261 while((rawStreamZDC->Next())){
262 if(rawStreamZDC->IsHeaderMapping()){ // mapping header
263 modNum = rawStreamZDC->GetADCModule();
264 modType = rawStreamZDC->GetModType();
265 }
266 if(rawStreamZDC->IsChMapping()){
267 if(modType==1){ // ADC mapping ----------------------
268 adcMod[ich] = rawStreamZDC->GetADCModFromMap(ich);
269 adcCh[ich] = rawStreamZDC->GetADCChFromMap(ich);
270 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
271 det[ich] = rawStreamZDC->GetDetectorFromMap(ich);
272 sec[ich] = rawStreamZDC->GetTowerFromMap(ich);
273 //
274 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
275 ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
276 //
277 //printf(" Mapping in DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
278 // ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
279 //
280 ich++;
281 }
282 else if(modType==2){ //VME scaler mapping --------------------
283 scMod[iScCh] = rawStreamZDC->GetScalerModFromMap(iScCh);
284 scCh[iScCh] = rawStreamZDC->GetScalerChFromMap(iScCh);
285 scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
286 scDet[iScCh] = rawStreamZDC->GetScDetectorFromMap(iScCh);
287 scSec[iScCh] = rawStreamZDC->GetScTowerFromMap(iScCh);
288 //
289 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
290 iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
291 //
292 //printf(" Mapping in DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
293 // iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
294 //
295 iScCh++;
296 }
442e1b18 297 }
298 }
299 }
442e1b18 300 fclose(mapFile4Shuttle);
9d0b658a 301 }// SOD event
ead118d8 302
d3f256ca 303 else if(eventT==PHYSICS_EVENT){
9d0b658a 304 // --- Reading data header
442e1b18 305 reader->ReadHeader();
ead118d8 306 const AliRawDataHeader* header = reader->GetDataHeader();
307 if(header) {
308 UChar_t message = header->GetAttributes();
a5747d8f 309 if((message & 0x30) == 0x30){ // DEDICATED LASER RUN
442e1b18 310 //printf("\t STANDALONE_LASER_RUN raw data found\n");
ead118d8 311 }
312 else{
9d0b658a 313 printf("ZDCLASERda.cxx -> NO STANDALONE_LASER_RUN raw data found\n");
ead118d8 314 return -1;
315 }
316 }
442e1b18 317 else{
318 printf("\t ATTENTION! No Raw Data Header found!!!\n");
319 return -1;
320 }
321
7f4bde92 322 rawStreamZDC->SetSODReading(kTRUE);
323
ead118d8 324 if (!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
442e1b18 325 //
326 // ----- Setting ch. mapping -----
9d0b658a 327 for(Int_t jk=0; jk<2*kNChannels; jk++){
442e1b18 328 rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
329 rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
330 rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
331 rawStreamZDC->SetMapDet(jk, det[jk]);
332 rawStreamZDC->SetMapTow(jk, sec[jk]);
333 }
ead118d8 334 //
335 while(rawStreamZDC->Next()){
336 Int_t index=-1;
bd94dbdd 337 Int_t detector = rawStreamZDC->GetSector(0);
9d0b658a 338 Int_t sector = rawStreamZDC->GetSector(1);
bd94dbdd 339
340 if(rawStreamZDC->IsADCDataWord() && !(rawStreamZDC->IsUnderflow())
341 && !(rawStreamZDC->IsOverflow()) && detector!=-1){
342
9d0b658a 343 if(sector!=5){ // Physics signals
19419b67 344 if(detector==1) index = sector; // *** ZNC
9d0b658a 345 else if(detector==2) index = sector+5; // *** ZPC
346 else if(detector==3) index = sector+9; // *** ZEM
347 else if(detector==4) index = sector+12;// *** ZNA
348 else if(detector==5) index = sector+17;// *** ZPA
bd94dbdd 349 }
350 else{ // Reference PMs
351 index = (detector-1)/3+22;
352 }
9d0b658a 353 //
354 if(index==-1) printf("ERROR in ZDCLASERda.cxx -> det %d quad %d res %d index %d ADC %d\n",
355 detector, sector, rawStreamZDC->GetADCGain(), index, rawStreamZDC->GetADCValue());
bd94dbdd 356
9d0b658a 357 Float_t Pedestal=0.;
19419b67 358 if(rawStreamZDC->GetADCGain()==0) Pedestal = MeanPedhg[index];
9d0b658a 359 else if(rawStreamZDC->GetADCGain()==1) Pedestal = MeanPedlg[index];
360 //
ead118d8 361 Float_t CorrADC = rawStreamZDC->GetADCValue() - Pedestal;
9d0b658a 362 //
363 //printf("\tdet %d sec %d res %d index %d ped %1.0f ADCcorr %1.0f\n",
364 // detector, sector, rawStreamZDC->GetADCGain(), index, Pedestal,CorrADC);
78beff0d 365
bd94dbdd 366 // **** Detector PMs
9d0b658a 367 if(sector!=5){
368 if(rawStreamZDC->GetADCGain()==0){ // --- High gain chain ---
369 // ---- side C
19419b67 370 if(detector==1) hZNChg[sector]->Fill(CorrADC);
9d0b658a 371 else if(detector==2) hZPChg[sector]->Fill(CorrADC);
372 // ---- side A
373 else if(detector==4) hZNAhg[sector]->Fill(CorrADC);
374 else if(detector==5) hZPAhg[sector]->Fill(CorrADC);
375 // ---- ZEM
376 else if(detector==3) hZEMhg[sector-1]->Fill(CorrADC);
377 }
378 else if(rawStreamZDC->GetADCGain()==1){ // --- Low gain chain ---
379 // ---- side C
19419b67 380 if(detector==1) hZNClg[sector]->Fill(CorrADC);
9d0b658a 381 else if(detector==2) hZPClg[sector]->Fill(CorrADC);
382 // ---- side A
383 else if(detector==4) hZNAlg[sector]->Fill(CorrADC);
384 else if(detector==5) hZPAlg[sector]->Fill(CorrADC);
385 // ---- ZEM
386 else if(detector==3) hZEMlg[sector-1]->Fill(CorrADC);
387 }
ead118d8 388 }
bd94dbdd 389 // **** Reference PMs
9d0b658a 390 else if(sector==5){
391 if(rawStreamZDC->GetADCGain()==0){ // --- High gain chain ---
392 // ---- PMRef chain side C
393 if(detector==1) hPMRefChg->Fill(CorrADC);
394 // ---- PMRef side A
395 else if(detector==4) hPMRefAhg->Fill(CorrADC);
396 }
397 else if(rawStreamZDC->GetADCGain()==1){ // --- Low gain chain ---
398 // ---- PMRef chain side C
399 if(detector==1) hPMRefClg->Fill(CorrADC);
400 // ---- PMRef side A
401 else if(detector==4) hPMRefAlg->Fill(CorrADC);
402 }
403 }
bd94dbdd 404 }//IsADCDataWord()+NOunderflow+NOoverflow
ead118d8 405 //
406 }
407 //
408 nevents_physics++;
409 //
410 delete reader;
411 delete rawStreamZDC;
412
413 }//(if PHYSICS_EVENT)
ead118d8 414
d3f256ca 415 /* exit when last event received, no need to wait for TERM signal */
416 else if(eventT==END_OF_RUN) {
417 printf(" -> EOR event detected\n");
418 break;
419 }
420
421
422 nevents_total++;
ead118d8 423
424 }
d3f256ca 425
426 /* free resources */
427 free(event);
ead118d8 428 }
429
430 /* Analysis of the histograms */
431 //
19419b67 432 Int_t detector[2*kNChannels], quad[2*kNChannels];
9d0b658a 433 Int_t maxBin[2*kNChannels], nBin[2*kNChannels];
434 Float_t xMax[2*kNChannels], maxXval[2*kNChannels], xlow[2*kNChannels];
435 Float_t mean[2*kNChannels], sigma[2*kNChannels];
19419b67 436 for(Int_t t=0; t<2*kNChannels; t++){
437 detector[t] = quad[t] = 0;
438 maxBin[t] = nBin[t] = 0;
439 xMax[t] = maxXval[t] = xlow[t] = 0.;
440 mean[t] = sigma[t] = 0.;
441 }
9d0b658a 442 TF1 *fun[2*kNChannels];
19419b67 443 Int_t atLeastOneHisto=0;
6f427255 444
9d0b658a 445 // ******** High gain chain ********
bd94dbdd 446 for(Int_t k=0; k<5; k++){
447 // --- ZNC
19419b67 448 detector[k] = 1;
9d0b658a 449 quad[k] = k;
450 maxBin[k] = hZNChg[k]->GetMaximumBin();
451 nBin[k] = (hZNChg[k]->GetXaxis())->GetNbins();
452 xMax[k] = (hZNChg[k]->GetXaxis())->GetXmax();
bd94dbdd 453 if(nBin[k]!=0) maxXval[k] = maxBin[k]*xMax[k]/nBin[k];
bd94dbdd 454 if(maxXval[k]-150.<0.) xlow[k]=0.;
455 else xlow[k] = maxXval[k]-150.;
19419b67 456 // checking if at least one histo is fitted
457 if(hZNChg[k]->GetEntries() != 0){
458 atLeastOneHisto=1;
459 //
460 hZNChg[k]->Fit("gaus","Q","",xlow[k],maxXval[k]+150.);
461 fun[k] = hZNChg[k]->GetFunction("gaus");
462 mean[k] = (Float_t) (fun[k]->GetParameter(1));
463 sigma[k] = (Float_t) (fun[k]->GetParameter(2));
464 }
bd94dbdd 465 // --- ZPC
19419b67 466 detector[k+5] = 2;
9d0b658a 467 quad[k+5] = k;
468 maxBin[k+5] = hZPChg[k]->GetMaximumBin();
469 nBin[k+5] = (hZPChg[k]->GetXaxis())->GetNbins();
470 xMax[k+5] = (hZPChg[k]->GetXaxis())->GetXmax();
bd94dbdd 471 if(nBin[k+5]!=0) maxXval[k+5] = maxBin[k+5]*xMax[k+5]/nBin[k+5];
bd94dbdd 472 if(maxXval[k+5]-150.<0.) xlow[k+5]=0.;
473 else xlow[k+5] = maxXval[k+5]-150.;
19419b67 474 if(hZPChg[k]->GetEntries() != 0){
475 atLeastOneHisto=1;
476 //
477 hZPChg[k]->Fit("gaus","Q","",xlow[k+5],maxXval[k+5]+150.);
478 fun[k+5] = hZPChg[k]->GetFunction("gaus");
479 mean[k+5] = (Float_t) (fun[k+5]->GetParameter(1));
480 sigma[k+5] = (Float_t) (fun[k+5]->GetParameter(2));
481 }
9d0b658a 482 // --- ZEM1
483 if(k<2){
19419b67 484 detector[k+10] = 3;
9d0b658a 485 quad[k+10] = k+1;
486 maxBin[k+10] = hZEMhg[k]->GetMaximumBin();
487 nBin[k+10] = (hZEMhg[k]->GetXaxis())->GetNbins();
488 xMax[k+10] = (hZEMhg[k]->GetXaxis())->GetXmax();
489 if(nBin[k+10]!=0) maxXval[k+10] = maxBin[k+10]*xMax[k+10]/nBin[k+10];
490 if(maxXval[k+10]-150.<0.) xlow[k+10]=0.;
491 else xlow[k+10] = maxXval[k+10]-150.;
19419b67 492 if(hZEMhg[k]->GetEntries() != 0){
493 atLeastOneHisto=1;
494 //
495 hZEMhg[k]->Fit("gaus","Q","",xlow[k+10],maxXval[k+10]+150.);
496 fun[k+10] = hZEMhg[k]->GetFunction("gaus");
497 mean[k+10] = (Float_t) (fun[k+10]->GetParameter(1));
498 sigma[k+10] = (Float_t) (fun[k+10]->GetParameter(2));
499 }
9d0b658a 500 }
bd94dbdd 501 // --- ZNA
19419b67 502 detector[k+12] = 4;
9d0b658a 503 quad[k+12] = k;
504 maxBin[k+12] = hZNAhg[k]->GetMaximumBin();
505 nBin[k+12] = (hZNAhg[k]->GetXaxis())->GetNbins();
506 xMax[k+12] = (hZNAhg[k]->GetXaxis())->GetXmax();
507 if(nBin[k+12]!=0) maxXval[k+12] = maxBin[k+12]*xMax[k+12]/nBin[k+12];
508 if(maxXval[k+12]-150.<0.) xlow[k+12]=0.;
509 else xlow[k+12] = maxXval[k+12]-150.;
19419b67 510 if(hZNAhg[k]->GetEntries() != 0){
511 atLeastOneHisto=1;
512 //
513 hZNAhg[k]->Fit("gaus","Q","",xlow[k+12],maxXval[k+12]+150.);
514 fun[k+12] = hZNAhg[k]->GetFunction("gaus");
515 mean[k+12] = (Float_t) (fun[k+12]->GetParameter(1));
516 sigma[k+12] = (Float_t) (fun[k+12]->GetParameter(2));
517 }
bd94dbdd 518 // --- ZPA
19419b67 519 detector[k+17] = 4;
9d0b658a 520 quad[k+17] = 5;
521 maxBin[k+17] = hZPAhg[k]->GetMaximumBin();
522 nBin[k+17] = (hZPAhg[k]->GetXaxis())->GetNbins();
523 xMax[k+17] = (hZPAhg[k]->GetXaxis())->GetXmax();
524 if(nBin[k+17]!=0) maxXval[k+17] = maxBin[k+17]*xMax[k+17]/nBin[k+17];
525 if(maxXval[k+17]-150.<0.) xlow[k+17]=0.;
526 else xlow[k+17] = maxXval[k+17]-150.;
19419b67 527 if(hZPAhg[k]->GetEntries() != 0){
528 atLeastOneHisto=1;
529 //
530 hZPAhg[k]->Fit("gaus","Q","",xlow[k+17],maxXval[k+17]+150.);
531 fun[k+17] = hZPAhg[k]->GetFunction("gaus");
532 mean[k+17] = (Float_t) (fun[k+17]->GetParameter(1));
533 sigma[k+17] = (Float_t) (fun[k+17]->GetParameter(2));
534 }
bd94dbdd 535 }
bd94dbdd 536 // ~~~~~~~~ PM Ref side C ~~~~~~~~
19419b67 537 detector[22] = 1;
9d0b658a 538 quad[22] = 5;
539 maxBin[22] = hPMRefChg->GetMaximumBin();
540 nBin[22] = (hPMRefChg->GetXaxis())->GetNbins();
541 xMax[22] = (hPMRefChg->GetXaxis())->GetXmax();
542 if(nBin[22]!=0) maxXval[22] = maxBin[22]*xMax[22]/nBin[22];
543 if(maxXval[22]-150.<0.) xlow[22]=0.;
19419b67 544 else xlow[22] = maxXval[22]-150.;
545 if(hPMRefChg->GetEntries() != 0){
546 atLeastOneHisto=1;
547 //
548 hPMRefChg->Fit("gaus","Q","",xlow[22],maxXval[22]+150.);
549 fun[22] = hPMRefChg->GetFunction("gaus");
550 mean[22] = (Float_t) (fun[22]->GetParameter(1));
551 sigma[22] = (Float_t) (fun[22]->GetParameter(2));
552 }
9d0b658a 553 // ~~~~~~~~ PM Ref side A ~~~~~~~~
19419b67 554 detector[23] = 4;
9d0b658a 555 quad[23] = 5;
556 maxBin[23] = hPMRefAhg->GetMaximumBin();
557 nBin[23] = (hPMRefAhg->GetXaxis())->GetNbins();
558 xMax[23] = (hPMRefAhg->GetXaxis())->GetXmax();
559 if(nBin[23]!=0) maxXval[23] = maxBin[23]*xMax[23]/nBin[23];
560 if(maxXval[23]-100.<0.) xlow[23]=0.;
19419b67 561 else xlow[23] = maxXval[23]-150.;
562 if(hPMRefAhg->GetEntries() != 0){
563 atLeastOneHisto=1;
564 //
565 hPMRefAhg->Fit("gaus","Q","",xlow[23],maxXval[23]+100.);
566 fun[23] = hPMRefAhg->GetFunction("gaus");
567 mean[23] = (Float_t) (fun[23]->GetParameter(1));
568 sigma[23] = (Float_t) (fun[23]->GetParameter(2));
569 }
6f427255 570
9d0b658a 571 // ******** Low gain chain ********
572 Int_t kOffset = 24;
573 for(Int_t k=0; k<5; k++){
574 // --- ZNC
19419b67 575 detector[k+kOffset] = 1;
9d0b658a 576 quad[k+kOffset] = k;
577 maxBin[k+kOffset] = hZNClg[k]->GetMaximumBin();
578 nBin[k+kOffset] = (hZNClg[k]->GetXaxis())->GetNbins();
579 xMax[k+kOffset] = (hZNClg[k]->GetXaxis())->GetXmax();
580 if(nBin[k+kOffset]!=0) maxXval[k+kOffset] = maxBin[k+kOffset]*xMax[k+kOffset]/nBin[k+kOffset];
581 if(maxXval[k+kOffset]-150.<0.) xlow[k+kOffset]=0.;
582 else xlow[k+kOffset] = maxXval[k+kOffset]-150.;
19419b67 583 if(hZNClg[k]->GetEntries() != 0){
584 atLeastOneHisto=1;
585 //
586 hZNClg[k]->Fit("gaus","Q","",xlow[k+kOffset],maxXval[k+kOffset]+150.);
587 fun[k+kOffset] = hZNClg[k]->GetFunction("gaus");
588 mean[k+kOffset] = (Float_t) (fun[k+kOffset]->GetParameter(1));
589 sigma[k+kOffset] = (Float_t) (fun[k+kOffset]->GetParameter(2));
590 }
9d0b658a 591 // --- ZPC
19419b67 592 detector[k+kOffset+5] = 2;
9d0b658a 593 quad[k+kOffset+5] = k;
594 maxBin[k+kOffset+5] = hZPClg[k]->GetMaximumBin();
595 nBin[k+kOffset+5] = (hZPClg[k]->GetXaxis())->GetNbins();
596 xMax[k+kOffset+5] = (hZPClg[k]->GetXaxis())->GetXmax();
597 if(nBin[k+kOffset+5]!=0) maxXval[k+kOffset+5] = maxBin[k+kOffset+5]*xMax[k+kOffset+5]/nBin[k+kOffset+5];
598 if(maxXval[k+kOffset+5]-150.<0.) xlow[k+kOffset+5]=0.;
599 else xlow[k+kOffset+5] = maxXval[k+kOffset+5]-150.;
19419b67 600 if(hZPClg[k]->GetEntries() != 0){
601 atLeastOneHisto=1;
602 //
603 hZPClg[k]->Fit("gaus","Q","",xlow[k+kOffset+5],maxXval[k+kOffset+5]+150.);
604 fun[k+kOffset+5] = hZPClg[k]->GetFunction("gaus");
605 mean[k+kOffset+5] = (Float_t) (fun[k+kOffset+5]->GetParameter(1));
606 sigma[k+kOffset+5] = (Float_t) (fun[k+kOffset+5]->GetParameter(2));
607 }
9d0b658a 608 // --- ZEM1
609 if(k+kOffset<2){
19419b67 610 detector[k+kOffset+10] = 3;
9d0b658a 611 quad[k+kOffset+10] = k+1;
612 maxBin[k+kOffset+10] = hZEMlg[k]->GetMaximumBin();
613 nBin[k+kOffset+10] = (hZEMlg[k]->GetXaxis())->GetNbins();
614 xMax[k+kOffset+10] = (hZEMlg[k]->GetXaxis())->GetXmax();
615 if(nBin[k+kOffset+10]!=0) maxXval[k+kOffset+10] = maxBin[k+kOffset+10]*xMax[k+kOffset+10]/nBin[k+kOffset+10];
616 if(maxXval[k+kOffset+10]-150.<0.) xlow[k+kOffset+10]=0.;
617 else xlow[k+kOffset+10] = maxXval[k+kOffset+10]-150.;
19419b67 618 if(hZEMlg[k]->GetEntries() != 0){
619 atLeastOneHisto=1;
620 //
621 hZEMlg[k]->Fit("gaus","Q","",xlow[k+kOffset+10],maxXval[k+kOffset+10]+150.);
622 fun[k+kOffset+10] = hZEMlg[k]->GetFunction("gaus");
623 mean[k+kOffset+10] = (Float_t) (fun[k+kOffset+10]->GetParameter(1));
624 sigma[k+kOffset+10] = (Float_t) (fun[k+kOffset+10]->GetParameter(2));
625 }
9d0b658a 626 }
627 // --- ZNA
19419b67 628 detector[k+kOffset+12] = 4;
9d0b658a 629 quad[k+kOffset+12] = k;
630 maxBin[k+kOffset+12] = hZNAlg[k]->GetMaximumBin();
631 nBin[k+kOffset+12] = (hZNAlg[k]->GetXaxis())->GetNbins();
632 xMax[k+kOffset+12] = (hZNAlg[k]->GetXaxis())->GetXmax();
633 if(nBin[k+kOffset+12]!=0) maxXval[k+kOffset+12] = maxBin[k+kOffset+12]*xMax[k+kOffset+12]/nBin[k+kOffset+12];
634 if(maxXval[k+kOffset+12]-150.<0.) xlow[k+kOffset+12]=0.;
635 else xlow[k+kOffset+12] = maxXval[k+kOffset+12]-150.;
19419b67 636 if(hZNAlg[k]->GetEntries() != 0){
637 atLeastOneHisto=1;
638 //
639 hZNAlg[k]->Fit("gaus","Q","",xlow[k+kOffset+12],maxXval[k+kOffset+12]+150.);
640 fun[k+kOffset+12] = hZNAlg[k]->GetFunction("gaus");
641 mean[k+kOffset+12] = (Float_t) (fun[k+kOffset+12]->GetParameter(1));
642 sigma[k+kOffset+12] = (Float_t) (fun[k+kOffset+12]->GetParameter(2));
643 }
9d0b658a 644 // --- ZPA
19419b67 645 detector[k+kOffset+17] = 5;
9d0b658a 646 quad[k+kOffset+17] = k;
647 maxBin[k+kOffset+17] = hZPAlg[k]->GetMaximumBin();
648 nBin[k+kOffset+17] = (hZPAlg[k]->GetXaxis())->GetNbins();
649 xMax[k+kOffset+17] = (hZPAlg[k]->GetXaxis())->GetXmax();
650 if(nBin[k+kOffset+17]!=0) maxXval[k+kOffset+17] = maxBin[k+kOffset+17]*xMax[k+kOffset+17]/nBin[k+kOffset+17];
651 if(maxXval[k+kOffset+17]-150.<0.) xlow[k+kOffset+17]=0.;
652 else xlow[k+kOffset+17] = maxXval[k+kOffset+17]-150.;
19419b67 653 if(hZPAlg[k]->GetEntries() != 0){
654 atLeastOneHisto=1;
655 //
656 hZPAlg[k]->Fit("gaus","Q","",xlow[k+kOffset+17],maxXval[k+kOffset+17]+150.);
657 fun[k+kOffset+17] = hZPAlg[k]->GetFunction("gaus");
658 mean[k+kOffset+17] = (Float_t) (fun[k+kOffset+17]->GetParameter(1));
659 sigma[k+kOffset+17] = (Float_t) (fun[k+kOffset+17]->GetParameter(2));
660 }
9d0b658a 661 }
662 // ~~~~~~~~ PM Ref side C ~~~~~~~~
19419b67 663 detector[46] = 1;
9d0b658a 664 quad[46] = 5;
665 maxBin[46] = hPMRefClg->GetMaximumBin();
666 nBin[46] = (hPMRefClg->GetXaxis())->GetNbins();
667 xMax[46] = (hPMRefClg->GetXaxis())->GetXmax();
668 if(nBin[46]!=0) maxXval[46] = maxBin[46]*xMax[46]/nBin[46];
669 if(maxXval[46]-150.<0.) xlow[46]=0.;
19419b67 670 else xlow[46] = maxXval[46]-150.;
671 if(hPMRefClg->GetEntries() != 0){
672 atLeastOneHisto=1;
673 //
674 hPMRefClg->Fit("gaus","Q","",xlow[46],maxXval[46]+150.);
675 fun[46] = hPMRefClg->GetFunction("gaus");
676 mean[46] = (Float_t) (fun[46]->GetParameter(1));
677 sigma[46] = (Float_t) (fun[46]->GetParameter(2));
678 }
bd94dbdd 679 // ~~~~~~~~ PM Ref side A ~~~~~~~~
19419b67 680 detector[47] = 4;
9d0b658a 681 quad[47] = 5;
682 maxBin[47] = hPMRefAlg->GetMaximumBin();
683 nBin[47] = (hPMRefAlg->GetXaxis())->GetNbins();
684 xMax[47] = (hPMRefAlg->GetXaxis())->GetXmax();
685 if(nBin[47]!=0) maxXval[47] = maxBin[47]*xMax[47]/nBin[47];
686 if(maxXval[47]-100.<0.) xlow[47]=0.;
19419b67 687 else xlow[47] = maxXval[47]-150.;
688 if(hPMRefAlg->GetEntries() != 0){
689 atLeastOneHisto=1;
690 //
691 hPMRefAlg->Fit("gaus","Q","",xlow[47],maxXval[47]+100.);
692 fun[47] = hPMRefAlg->GetFunction("gaus");
693 mean[47] = (Float_t) (fun[47]->GetParameter(1));
694 sigma[47] = (Float_t) (fun[47]->GetParameter(2));
695 }
696
697 if(atLeastOneHisto==0){
698 printf("\n WARNING! Empty LASER histos -> ending DA WITHOUT writing output\n\n");
699 return -1;
700 }
bd94dbdd 701
ead118d8 702 FILE *fileShuttle;
218f916a 703 fileShuttle = fopen(LASDATA_FILE,"w");
9d0b658a 704 for(Int_t i=0; i<2*kNChannels; i++){
19419b67 705 fprintf(fileShuttle,"\t%d\t%d\t%f\t%f\n",detector[i],quad[i],mean[i], sigma[i]);
bd94dbdd 706 }
ead118d8 707 //
708 fclose(fileShuttle);
65448375 709 /* report progress */
710 daqDA_progressReport(80);
711 //
712 TFile *histofile = new TFile(LASHISTO_FILE,"RECREATE");
713 histofile->cd();
714 for(int j=0; j<5; j++){
715 hZNChg[j]->Write();
716 hZPChg[j]->Write();
717 hZNAhg[j]->Write();
718 hZPAhg[j]->Write();
719 hZNClg[j]->Write();
720 hZPClg[j]->Write();
721 hZNAlg[j]->Write();
722 hZPAlg[j]->Write();
723 if(j<2){
724 hZEMhg[j]->Write();
725 hZEMlg[j]->Write();
726 }
727 }
728 hPMRefChg->Write();
729 hPMRefAhg->Write();
730 hPMRefClg->Write();
731 hPMRefAlg->Write();
732 //
733 histofile->Close();
ead118d8 734 //
bd94dbdd 735 for(Int_t j=0; j<5; j++){
9d0b658a 736 delete hZNChg[j];
737 delete hZPChg[j];
738 delete hZNAhg[j];
739 delete hZPAhg[j];
740 delete hZNClg[j];
741 delete hZPClg[j];
742 delete hZNAlg[j];
743 delete hZPAlg[j];
744 if(j<2){
745 delete hZEMhg[j];
746 delete hZEMlg[j];
747 }
bd94dbdd 748 }
9d0b658a 749 delete hPMRefChg;
750 delete hPMRefAhg;
751 delete hPMRefClg;
752 delete hPMRefAlg;
442e1b18 753
ead118d8 754 /* write report */
755 fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
756
757 /* close result file */
758 fclose(fp);
759
760 /* report progress */
761 daqDA_progressReport(90);
9d0b658a 762
ead118d8 763 /* store the result file on FES */
65448375 764 // [1] File with mapping
a8eaff64 765 status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
442e1b18 766 if(status){
767 printf("Failed to export file : %d\n",status);
768 return -1;
769 }
770 //
65448375 771 // [2] File with laser data
a8eaff64 772 status = daqDA_FES_storeFile(LASDATA_FILE, "LASERDATA");
ead118d8 773 if(status){
774 printf("Failed to export file : %d\n",status);
775 return -1;
776 }
65448375 777 // [3] File with laser histos
a8eaff64 778 status = daqDA_FES_storeFile(LASHISTO_FILE, "LASERHISTOS");
65448375 779 if(status){
780 printf("Failed to export pedestal histos file to DAQ FES\n");
781 return -1;
782 }
ead118d8 783
784 /* report progress */
785 daqDA_progressReport(100);
786
ead118d8 787 return status;
788}