/* This program reads the DAQ data files passed as argument using the monitoring library. It computes the average event size and populates local "./result.txt" file with the result. The program reports about its processing progress. Messages on stdout are exported to DAQ log system. DA for ZDC standalone pedestal runs Contact: Chiara.Oppedisano@to.infn.it Link: Run Type: STANDALONE_PEDESTAL_RUN DA Type: LDC Number of events needed: no constraint (tipically ~10^3) Input Files: Output Files: ZDCPedestal.dat, ZDCChMapping.dat Trigger Types Used: Standalone Trigger */ #include #include #include // DATE #include #include #include //ROOT #include #include #include #include #include #include #include //AliRoot #include #include #include /* Main routine Arguments: list of DATE raw data files */ int main(int argc, char **argv) { TFitter *minuitFit = new TFitter(4); TVirtualFitter::SetFitter(minuitFit); int status = 0; /* log start of process */ printf("\n ZDC PEDESTAL program started\n"); /* check that we got some arguments = list of files */ if (argc<2) { printf("Wrong number of arguments\n"); return -1; } // --- Histograms for ADC pedestals // [22 signal channels + 2 reference PTMs] x 2 gain chains // TH1F::AddDirectory(0); int const kNChannels = 24; TH1F *hPedhg[kNChannels], *hPedOutOfTimehg[kNChannels]; TH2F *hPedCorrhg[kNChannels]; TH1F *hPedlg[kNChannels], *hPedOutOfTimelg[kNChannels]; TH2F *hPedCorrlg[kNChannels]; // char namhist1hg[50], namhist2hg[50], namhist3hg[50]; char namhist1lg[50], namhist2lg[50], namhist3lg[50]; for(Int_t j=0; j=5 && j<10){ // ZP1 sprintf(namhist1hg,"PedZP1hg_%d",j-5); sprintf(namhist2hg,"PedZP1hgOutOfTime_%d",j-5); sprintf(namhist3hg,"PedCorrZP1hg_%d",j-5); // sprintf(namhist1lg,"PedZP1lg_%d",j-5); sprintf(namhist2lg,"PedZP1lgOutOfTime_%d",j-5); sprintf(namhist3lg,"PedCorrZP1lg_%d",j-5); } else if(j>=10 && j<12){ // ZEM sprintf(namhist1hg,"PedZEMhg_%d",j-10); sprintf(namhist2hg,"PedZEMhgOutOfTime_%d",j-10); sprintf(namhist3hg,"PedCorrZEMhg_%d",j-10); // sprintf(namhist1lg,"PedZEMlg_%d",j-10); sprintf(namhist2lg,"PedZEMlgOutOfTime_%d",j-10); sprintf(namhist3lg,"PedCorrZEMlg_%d",j-10); } else if(j>=12 && j<17){ // ZN2 sprintf(namhist1hg,"PedZN2hg_%d",j-12); sprintf(namhist2hg,"PedZN2hgOutOfTime_%d",j-12); sprintf(namhist3hg,"PedCorrZN2hg_%d",j-12); // sprintf(namhist1lg,"PedZN2lg_%d",j-12); sprintf(namhist2lg,"PedZN2lgOutOfTime_%d",j-12); sprintf(namhist3lg,"PedCorrZN2lg_%d",j-12); } else if(j>=17 && j<22){ // ZP2 sprintf(namhist1hg,"PedZP2hg_%d",j-17); sprintf(namhist2hg,"PedZP2hgOutOfTime_%d",j-17); sprintf(namhist3hg,"PedCorrZP2hg_%d",j-17); // sprintf(namhist1lg,"PedZP2lg_%d",j-17); sprintf(namhist2lg,"PedZP2lgOutOfTime_%d",j-17); sprintf(namhist3lg,"PedCorrZP2lg_%d",j-17); } else if(j>=22 && j<24){ //Reference PMs sprintf(namhist1hg,"PedRefhg_%d",j-22); sprintf(namhist2hg,"PedRefhgOutOfTime_%d",j-22); sprintf(namhist3hg,"PedCorrRefhg_%d",j-22); // sprintf(namhist1lg,"PedReflg_%d",j-22); sprintf(namhist2lg,"PedReflgOutOfTime_%d",j-22); sprintf(namhist3lg,"PedCorrReflg_%d",j-22); } // --- High gain chain histos hPedhg[j] = new TH1F(namhist1hg, namhist1hg, 200,0., 200.); hPedOutOfTimehg[j] = new TH1F(namhist2hg, namhist2hg, 200,0., 200.); hPedCorrhg[j] = new TH2F(namhist3hg,namhist3hg,100,0.,200.,100,0.,200.); // --- Low gain chain histos hPedlg[j] = new TH1F(namhist1lg, namhist1lg, 100,0., 600.); hPedOutOfTimelg[j] = new TH1F(namhist2lg, namhist2lg, 100,0., 600.); hPedCorrlg[j] = new TH2F(namhist3lg,namhist3lg,100,0.,600.,100,0.,600.); } /* open result file */ FILE *fp=NULL; fp=fopen("./result.txt","a"); if (fp==NULL) { printf("Failed to open file\n"); return -1; } FILE *mapFile4Shuttle; const char *mapfName = "ZDCChMapping.dat"; /* report progress */ daqDA_progressReport(10); /* init some counters */ int nevents_physics=0; int nevents_total=0; /* read the data files */ int n; for(n=1;nSelect("ZDC"); // --- Reading event header //UInt_t evtype = reader->GetType(); //printf("\n\t ZDCPEDESTALda -> ev. type %d\n",evtype); //printf("\t ZDCPEDESTALda -> run # %d\n",reader->GetRunNumber()); // AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader); /* use event - here, just write event id to result file */ eventT=event->eventType; Int_t ich=0, adcMod[48], adcCh[48], sigCode[48], det[48], sec[48]; if(eventT==START_OF_DATA){ if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n"); else{ while(rawStreamZDC->Next()){ if(rawStreamZDC->IsChMapping()){ adcMod[ich] = rawStreamZDC->GetADCModFromMap(ich); adcCh[ich] = rawStreamZDC->GetADCChFromMap(ich); sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich); det[ich] = rawStreamZDC->GetDetectorFromMap(ich); sec[ich] = rawStreamZDC->GetTowerFromMap(ich); ich++; } } } // -------------------------------------------------------- // --- Writing ascii data file for the Shuttle preprocessor mapFile4Shuttle = fopen(mapfName,"w"); for(Int_t i=0; i ch.%d mod %d, ch %d, code %d det %d, sec %d\n", // i,adcMod[i],adcCh[i],sigCode[i],det[i],sec[i]); } fclose(mapFile4Shuttle); } if(eventT==PHYSICS_EVENT){ fprintf(fp,"Run #%lu, event size: %lu, BC:%u, Orbit:%u, Period:%u\n", (unsigned long)event->eventRunNb, (unsigned long)event->eventSize, EVENT_ID_GET_BUNCH_CROSSING(event->eventId), EVENT_ID_GET_ORBIT(event->eventId), EVENT_ID_GET_PERIOD(event->eventId)); // --- Reading data header reader->ReadHeader(); const AliRawDataHeader* header = reader->GetDataHeader(); if(header){ UChar_t message = header->GetAttributes(); if(message & 0x20){ // PEDESTAL RUN //printf("\t STANDALONE_PEDESTAL RUN raw data found\n"); } else{ printf("\t NO STANDALONE_PEDESTAL RUN raw data found\n"); return -1; } } else{ printf("\t ATTENTION! No Raw Data Header found!!!\n"); return -1; } if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n"); // // ----- Setting ch. mapping ----- for(Int_t jk=0; jk<48; jk++){ rawStreamZDC->SetMapADCMod(jk, adcMod[jk]); rawStreamZDC->SetMapADCCh(jk, adcCh[jk]); rawStreamZDC->SetMapADCSig(jk, sigCode[jk]); rawStreamZDC->SetMapDet(jk, det[jk]); rawStreamZDC->SetMapTow(jk, sec[jk]); } // Int_t iraw=0; Int_t RawADChg[kNChannels], RawADCoothg[kNChannels]; Int_t RawADClg[kNChannels], RawADCootlg[kNChannels]; for(Int_t j=0; jNext()){ Int_t index=-1; if(rawStreamZDC->IsADCDataWord()){ if(rawStreamZDC->GetSector(1)!=5){ // Physics signals if(rawStreamZDC->GetSector(0)==1) index = rawStreamZDC->GetSector(1); // *** ZNC else if(rawStreamZDC->GetSector(0)==2) index = rawStreamZDC->GetSector(1)+5; // *** ZPC else if(rawStreamZDC->GetSector(0)==3) index = rawStreamZDC->GetSector(1)+9; // *** ZEM else if(rawStreamZDC->GetSector(0)==4) index = rawStreamZDC->GetSector(1)+12; // *** ZNA else if(rawStreamZDC->GetSector(0)==5) index = rawStreamZDC->GetSector(1)+17; // *** ZPA } else{ // Reference PMs index = (rawStreamZDC->GetSector(0)-1)/3+22; } // /*printf("\t iraw %d index %d det %d quad %d res %d ADC %d\n", iraw, index, rawStreamZDC->GetSector(0), rawStreamZDC->GetSector(1), rawStreamZDC->GetADCGain(), rawStreamZDC->GetADCValue()); */ // if(iraw<2*kNChannels){ // --- In-time pedestals (1st 48 raw data) if(rawStreamZDC->GetADCGain()==0){ hPedhg[index]->Fill(rawStreamZDC->GetADCValue()); RawADChg[index] = rawStreamZDC->GetADCValue(); // //printf("\t filling histo hPedhg[%d]\n",index); } else{ hPedlg[index]->Fill(rawStreamZDC->GetADCValue()); RawADClg[index] = rawStreamZDC->GetADCValue(); // //printf("\t filling histo hPedlg[%d]\n",index); } } else{ // --- Out-of-time pedestals if(rawStreamZDC->GetADCGain()==0){ hPedOutOfTimehg[index]->Fill(rawStreamZDC->GetADCValue()); RawADCoothg[index] = rawStreamZDC->GetADCValue(); // //printf("\t filling histo hPedOutOfTimehg[%d]\n",index); } else{ hPedOutOfTimelg[index]->Fill(rawStreamZDC->GetADCValue()); RawADCootlg[index] = rawStreamZDC->GetADCValue(); // //printf("\t filling histo hPedOutOfTimelg[%d]\n",index); } } iraw++; }//IsADCDataWord() // if(iraw == 4*kNChannels){ // Last ADC channel -> Filling correlation histos for(Int_t k=0; kFill(RawADCoothg[k], RawADChg[k]); hPedCorrlg[k]->Fill(RawADCootlg[k], RawADClg[k]); } } } // nevents_physics++; // delete reader; delete rawStreamZDC; }//(if PHYSICS_EVENT) nevents_total++; /* free resources */ free(event); } } /* Analysis of the histograms */ // FILE *fileShuttle; const char *fName = "ZDCPedestal.dat"; fileShuttle = fopen(fName,"w"); // Float_t MeanPed[2*kNChannels], MeanPedWidth[2*kNChannels], MeanPedOOT[2*kNChannels], MeanPedWidthOOT[2*kNChannels]; // --- In-time pedestals TF1 *ADCfunchg[kNChannels]; for(Int_t i=0; iFit("gaus","Q"); ADCfunchg[i] = hPedhg[i]->GetFunction("gaus"); MeanPed[i] = (Double_t) ADCfunchg[i]->GetParameter(1); MeanPedWidth[i] = (Double_t) ADCfunchg[i]->GetParameter(2); fprintf(fileShuttle,"\t%d\t%f\t%f\n",i,MeanPed[i],MeanPedWidth[i]); //printf("\t MeanPed[%d] = %f\n",i, MeanPed[i]); } TF1 *ADCfunclg[kNChannels]; for(Int_t i=0; iFit("gaus","Q"); ADCfunclg[i] = hPedlg[i]->GetFunction("gaus"); MeanPed[i+kNChannels] = (Double_t) ADCfunclg[i]->GetParameter(1); MeanPedWidth[i+kNChannels] = (Double_t) ADCfunclg[i]->GetParameter(2); fprintf(fileShuttle,"\t%d\t%f\t%f\n",i+kNChannels,MeanPed[i+kNChannels],MeanPedWidth[i+kNChannels]); //printf("\t MeanPed[%d] = %f\n",i+kNChannels, MeanPed[i+kNChannels]); } // --- Out-of-time pedestals TF1 *ADCootfunchg[kNChannels]; for(Int_t i=0; iFit("gaus","Q"); ADCootfunchg[i] = hPedOutOfTimehg[i]->GetFunction("gaus"); MeanPedOOT[i] = (Double_t) ADCootfunchg[i]->GetParameter(1); MeanPedWidthOOT[i] = (Double_t) ADCootfunchg[i]->GetParameter(2); fprintf(fileShuttle,"\t%d\t%f\t%f\n",i,MeanPedOOT[i],MeanPedWidthOOT[i]); //printf("\t MeanPedOOT[%d] = %f\n",i, MeanPedOOT[i]); } TF1 *ADCootfunclg[kNChannels]; for(Int_t i=0; iFit("gaus","Q"); ADCootfunclg[i] = hPedOutOfTimelg[i]->GetFunction("gaus"); MeanPedOOT[i+kNChannels] = (Double_t) ADCootfunclg[i]->GetParameter(1); MeanPedWidthOOT[i+kNChannels] = (Double_t) ADCootfunclg[i]->GetParameter(2); fprintf(fileShuttle,"\t%d\t%f\t%f\n",i+kNChannels,MeanPedOOT[i+kNChannels],MeanPedWidthOOT[i+kNChannels]); //printf("\t MeanPedOOT[%d] = %f\n",i+kNChannels, MeanPedOOT[i+kNChannels]); } // *************************************************** // Unless we have a narrow correlation to fit we // don't fit and store in-time vs. out-of-time // histograms -> mean pedstal subtracted!!!!!! // *************************************************** // --- Correlations /* Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels]; TProfile *hPedCorrProfhg[kNChannels], *hPedCorrProflg[kNChannels]; TF1 *ffunchg[kNChannels], *ffunclg[kNChannels]; char namhist4[50]; for(int i=0;iProfileX(namhist4,-1,-1,"S"); hPedCorrProfhg[i]->SetName(namhist4); hPedCorrProfhg[i]->Fit("pol1","Q"); ffunchg[i] = hPedCorrProfhg[i]->GetFunction("pol1"); CorrCoeff0[i] = (Double_t) ffunchg[i]->GetParameter(0); CorrCoeff1[i] = (Double_t) ffunchg[i]->GetParameter(1); fprintf(fileShuttle,"\t%d\t%f\t%f\n",i,CorrCoeff0[i],CorrCoeff1[i]); //printf("\t CorrCoeff0[%d] = %f, CorrCoeff1[%d] = %f\n",i, CorrCoeff0[i], i, CorrCoeff1[i]); } for(int i=0;iProfileX(namhist4,-1,-1,"S"); hPedCorrProflg[i]->SetName(namhist4); hPedCorrProflg[i]->Fit("pol1","Q"); ffunclg[i] = hPedCorrProflg[i]->GetFunction("pol1"); CorrCoeff0[i+kNChannels] = (Double_t) ffunclg[i]->GetParameter(0); CorrCoeff1[i+kNChannels] = (Double_t) ffunclg[i]->GetParameter(1); fprintf(fileShuttle,"\t%d\t%f\t%f\n",i+kNChannels,CorrCoeff0[i+kNChannels],CorrCoeff1[i+kNChannels]); //printf("\t CorrCoeff0[%d] = %f, CorrCoeff1[%d] = %f\n", // i+kNChannels, CorrCoeff0[i+kNChannels], i+kNChannels, CorrCoeff1[i+kNChannels]); } */ // fclose(fileShuttle); // for(Int_t j=0; j