]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/DA/PHOSCPVPEDda.cxx
PHOS: CPV PED DA added
[u/mrichter/AliRoot.git] / PHOS / DA / PHOSCPVPEDda.cxx
1 /*
2 CPV PED DA for processing pedestal runs and producing pedestal tables for loading to CPV FEE.
3
4 Contact: Sergey Evdokimov <sevdokim@cern.ch>
5 Link: https://twiki.cern.ch/twiki/bin/view/ALICE/CPVda
6 Reference run: 211758 (/afs/cern.ch/user/s/sevdokim/public/CPV_run211758_pedestal.raw)
7 Run Type:  PHYSICS
8 DA Type: PED
9 Number of events needed: 2000 events
10 Input files: raw data file
11 Output files: thr?_??.dat CpvPeds.root
12 Trigger types used: PHYSICS_EVENT
13 */
14
15 #include "event.h"
16 #include "monitor.h"
17 #include "daqDA.h"
18
19 #include <Riostream.h>
20 #include <stdlib.h>
21 #include <fstream>
22 #include <string>
23
24 //AliRoot
25 #include "AliPHOSCpvRawStream.h"
26 #include "AliPHOSCpvPedProducer.h"
27 #include "AliPHOSCpvParam.h"
28 #include "AliRawReaderDate.h"
29 #include "AliBitPacking.h"
30 #include "TMath.h"
31
32 //ROOT
33 #include "TROOT.h"
34 #include "TPluginManager.h"
35 #include "TSAXParser.h"
36 #include "TTree.h"
37
38 #include "TString.h"
39 #include "TFile.h"
40 #include "TSystem.h"
41 #include "TKey.h"
42 #include "TH2S.h"
43 #include "TH2F.h"
44 #include "TObject.h"
45 #include "TBenchmark.h"
46 #include "TMath.h"
47 #include "TRandom.h"
48
49
50 int main( int argc, char **argv )
51 {
52   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
53                                         "*",
54                                         "TStreamerInfo",
55                                         "RIO",
56                                         "TStreamerInfo()");
57   Int_t status,print;
58   Int_t sigcut=3;
59   Bool_t turbo = kTRUE;
60
61   if (argc!=2) {
62     printf("Wrong number of arguments\n");
63     return -1;
64   }
65
66   // log start of process
67   printf("Cpv DA program started\n");
68
69   /* report progress */
70   daqDA_progressReport(0);
71
72   status=monitorSetDataSource( argv[1] );
73   if (status!=0) {
74     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
75     return -1;
76   }
77
78   /* declare monitoring program */
79   status=monitorDeclareMp( __FILE__ );
80   if (status!=0) {
81     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
82     return -1;
83   }
84
85   /* define wait event timeout - 1s max */
86   monitorSetNowait();
87   monitorSetNoWaitNetworkTimeout(1000);
88
89     /* report progress */
90   daqDA_progressReport(5);
91
92
93   // init event counter
94   Int_t iPhysEvnt=0;
95   Int_t iTotEvnt =0;
96
97   AliPHOSCpvPedProducer * pedProducer = new AliPHOSCpvPedProducer(sigcut);
98
99   // Reader
100   AliRawReader * reader;
101
102   /* report progress */
103   daqDA_progressReport(10);
104
105   /* main loop (infinite) */
106   for(;;) { // infinite loop
107     struct eventHeaderStruct *event;
108     eventTypeType eventT;
109
110     /* check shutdown condition */
111     if (daqDA_checkShutdown()) {break;}
112
113     // get next event
114     status=monitorGetEventDynamic((void **)&event);
115     if (status==MON_ERR_EOF) { // end of monitoring file has been reached
116       printf("End of monitoring file has been reached! \n");
117       break;
118     }
119
120     if (status!=0) {
121       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
122       break;
123     }
124
125     // retry if got no event
126     if (event==NULL) continue;
127
128     // use event - here, just write event id to result file
129     eventT=event->eventType;
130     if (eventT==PHYSICS_EVENT) { //we use PHYSICS_EVENT for pedestal not CALIBRATION_EVENT???
131       iTotEvnt++;
132       reader = new AliRawReaderDate((void*)event);
133       pedProducer->LoadNewEvent(reader);
134       if(pedProducer->FillPedestal()) iPhysEvnt++;
135       delete reader;
136     } // if PHYSICS_EVENT
137
138     free(event);
139
140     /* exit when last event received, no need to wait for TERM signal */
141     if (eventT==END_OF_RUN) {
142       printf("EOR event detected\n");
143       break;
144     }
145   }
146
147   Printf(" Received %d events, %d good events",iTotEvnt,iPhysEvnt);
148   /* report progress */
149   daqDA_progressReport(80);
150
151   for(int iDDL = 0; iDDL<AliPHOSCpvParam::kNDDL; iDDL++){
152     if(pedProducer -> CalcPedestal(iDDL)){
153       pedProducer -> WritePedFiles(iDDL);
154       //for (int iCC = 0; iCC<AliPHOSCpvParam::kNRows){
155       //        status=daqDA_DB_storeFile(Form("thr%d_%02d.dat", iDDL, iCC));
156       //        if(status) printf("Failed to store thr%d_%02d.dat in DAQ DB!\n",iDDL, iCC);
157       //        status=daqDA_FES_storeFile(Form("thr%d_%02d.dat", iDDL, iCC));
158       //        if(status) printf("Failed to export thr%d_%02d.dat to DAQ FES!\n",iDDL, iCC);
159       //}
160     }
161   }
162
163   pedProducer->WriteAllHistsToFile("CpvPeds.root");
164   status = daqDA_DB_storeFile("CpvPeds.root","CpvPeds.root");
165   if(status) printf("Failed to store CpvPeds.root in DAQ DB!\n");
166
167   /* report progress */
168   daqDA_progressReport(95);
169
170
171   return status;
172 }