]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/DA/CPVPEDda.cxx
CMake: DA rpm creation, enable more detectors from the same folder
[u/mrichter/AliRoot.git] / PHOS / DA / CPVPEDda.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 //daqDA
15 #include "event.h"
16 #include "monitor.h"
17 #include "daqDA.h"
18 //AMORE monitoring framework
19 #include <AmoreDA.h>
20
21 //system
22 #include <Riostream.h>
23 #include <stdlib.h>
24 #include <fstream>
25 #include <string>
26
27 //AliRoot
28 #include "AliPHOSCpvRawStream.h"
29 #include "AliPHOSCpvPedProducer.h"
30 #include "AliPHOSCpvParam.h"
31 #include "AliRawReaderDate.h"
32 #include "AliBitPacking.h"
33
34 //ROOT
35 #include "TROOT.h"
36 #include "TPluginManager.h"
37 #include "TSAXParser.h"
38 #include "TTree.h"
39 #include "TMath.h"
40 #include "TString.h"
41 #include "TFile.h"
42 #include "TSystem.h"
43 #include "TKey.h"
44 #include "TH2S.h"
45 #include "TH2F.h"
46 #include "TObject.h"
47 #include "TBenchmark.h"
48 #include "TMath.h"
49 #include "TRandom.h"
50
51
52 int main( int argc, char **argv )
53 {
54   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
55                                         "*",
56                                         "TStreamerInfo",
57                                         "RIO",
58                                         "TStreamerInfo()");
59   Int_t status,print;
60   Int_t sigcut=3;
61   Bool_t turbo = kTRUE;
62
63   if (argc!=2) {
64     printf("Wrong number of arguments\n");
65     return -1;
66   }
67
68   // log start of process
69   printf("Cpv pedestal DA program started\n");
70
71   /* report progress */
72   daqDA_progressReport(0);
73
74   status=monitorSetDataSource( argv[1] );
75   if (status!=0) {
76     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
77     return -1;
78   }
79
80   /* declare monitoring program */
81   status=monitorDeclareMp( __FILE__ );
82   if (status!=0) {
83     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
84     return -1;
85   }
86
87   /* define wait event timeout - 1s max */
88   monitorSetNowait();
89   monitorSetNoWaitNetworkTimeout(1000);
90
91     /* report progress */
92   daqDA_progressReport(5);
93
94
95   // init event counter
96   Int_t iPhysEvnt=0;
97   Int_t iTotEvnt =0;
98
99   AliPHOSCpvPedProducer * pedProducer = new AliPHOSCpvPedProducer(sigcut);
100
101   // Reader
102   AliRawReader * reader;
103
104   /* report progress */
105   daqDA_progressReport(10);
106
107   /* main loop (infinite) */
108   for(;;) { // infinite loop
109     struct eventHeaderStruct *event;
110     eventTypeType eventT;
111
112     /* check shutdown condition */
113     if (daqDA_checkShutdown()) {break;}
114
115     // get next event
116     status=monitorGetEventDynamic((void **)&event);
117     if (status==MON_ERR_EOF) { // end of monitoring file has been reached
118       printf("End of monitoring file has been reached! \n");
119       break;
120     }
121
122     if (status!=0) {
123       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
124       break;
125     }
126
127     // retry if got no event
128     if (event==NULL) continue;
129
130     // use event - here, just write event id to result file
131     eventT=event->eventType;
132     if (eventT==PHYSICS_EVENT) { //we use PHYSICS_EVENT for pedestal not CALIBRATION_EVENT???
133       iTotEvnt++;
134       reader = new AliRawReaderDate((void*)event);
135       pedProducer->LoadNewEvent(reader);
136       if(pedProducer->FillPedestal()) iPhysEvnt++;
137       delete reader;
138     } // if PHYSICS_EVENT
139
140     free(event);
141
142     /* exit when last event received, no need to wait for TERM signal */
143     if (eventT==END_OF_RUN) {
144       printf("EOR event detected\n");
145       break;
146     }
147   }
148
149   Printf(" Received %d events, %d good events",iTotEvnt,iPhysEvnt);
150   /* report progress */
151   daqDA_progressReport(95);
152
153   pedProducer->WriteAllHistsToFile("CpvPeds.root");
154
155   //send file with histos to amore
156   amore::da::AmoreDA* myAmore = new amore::da::AmoreDA(amore::da::AmoreDA::kSender);
157   TFile* fSend = TFile::Open("CpvPeds.root");
158   TH2F *pedsMean = (TH2F*)fSend->Get("fPedMeanMap4");
159   pedsMean->SetDrawOption("colz");
160   myAmore->Send("fPedMeanMap4",pedsMean);
161   TH2F *pedsSig = (TH2F*)fSend->Get("fPedSigMap4");
162   pedsSig->SetDrawOption("colz");
163   myAmore->Send("fPedSigMap4",pedsSig);
164   TH1F *h1DPedMean=(TH1F*)fSend->Get("f1DPedMean4");
165   h1DPedMean->GetXaxis()->SetRangeUser(0.,750.);
166   myAmore->Send("f1DPedMean4",h1DPedMean);
167   TH1F *h1DPedSig=(TH1F*)fSend->Get("f1DPedSig4");
168   h1DPedSig->GetXaxis()->SetRangeUser(0.,10.);
169   myAmore->Send("f1DPedSig4",h1DPedSig);
170   fSend->Close();
171
172   if(iPhysEvnt>=1000){//we have enough events to publish data
173     status = daqDA_DB_storeFile("CpvPeds.root","CpvPeds.root");
174     if(status) printf("Failed to store CpvPeds.root in DAQ DB!\n");
175   
176     status = daqDA_FES_storeFile("CpvPeds.root","CpvPeds.root");
177     if(status) printf("Failed to store CpvPeds.root in DAQ FXS!\n");
178
179     for(int iDDL = 0; iDDL<2*AliPHOSCpvParam::kNDDL; iDDL++){
180       if(pedProducer -> CalcPedestal(iDDL)){
181         pedProducer -> WritePedFiles(iDDL);
182         for (int iCC = 0; iCC<AliPHOSCpvParam::kNRows; iCC++){
183           status=daqDA_DB_storeFile(Form("thr%d_%02d.dat", iDDL, iCC),Form("thr%d_%02d.dat", iDDL, iCC));
184           if(status) printf("Failed to store thr%d_%02d.dat in DAQ DB!\n",iDDL, iCC);
185         }
186       }
187     }
188   }else return 10;//error code 10 (not enough events!)
189
190
191   /* report progress */
192   daqDA_progressReport(100);
193
194
195   return 0;
196 }