]>
Commit | Line | Data |
---|---|---|
1 | /* | |
2 | ||
3 | HMPID DA for online calibration | |
4 | ||
5 | Contact: Levente.Molnar@ba.infn.it, Giacomo.Volpe@ba.infn.it | |
6 | Link: http://richpc1.ba.infn.it/~levente/Files4Public/ValidateHmpidDA/ | |
7 | Run Type: PEDESTAL -- but we select on the PHYSICS_EVENTS in th HMPIDda.cxx | |
8 | DA Type: LDC | |
9 | Number of events needed: 1000 events | |
10 | Input Files: Raw pedestal file, EXTERNAL config file: HmpidSigmaCut.txt on both HMPID LDCs | |
11 | Output Files: 14 txt files including pedestal values | |
12 | Trigger types used: PEDESTAL RUN (selecting on PHYSICS_EVENT) | |
13 | ||
14 | */ | |
15 | ||
16 | extern "C" { | |
17 | #include <daqDA.h> | |
18 | } | |
19 | ||
20 | #include "event.h" | |
21 | #include "monitor.h" | |
22 | ||
23 | #include <Riostream.h> | |
24 | #include "fstream.h" | |
25 | #include <stdio.h> | |
26 | #include <stdlib.h> | |
27 | ||
28 | //AliRoot | |
29 | #include "AliHMPIDRawStream.h" | |
30 | #include "AliHMPIDCalib.h" | |
31 | #include "AliRawReaderDate.h" | |
32 | #include "AliBitPacking.h" | |
33 | ||
34 | //ROOT | |
35 | #include "TROOT.h" | |
36 | #include "TFile.h" | |
37 | #include "TSystem.h" | |
38 | #include "TString.h" | |
39 | #include "TObject.h" | |
40 | #include "TPluginManager.h" | |
41 | ||
42 | ||
43 | int main(int argc, char **argv){ | |
44 | ||
45 | int status; | |
46 | ||
47 | /* log start of process */ | |
48 | printf("HMPID DA program started\n"); | |
49 | ||
50 | gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo", | |
51 | "*", | |
52 | "TStreamerInfo", | |
53 | "RIO", | |
54 | "TStreamerInfo()"); | |
55 | ||
56 | /* check that we got some arguments = list of files */ | |
57 | if (argc<2) { | |
58 | printf("Wrong number of arguments\n"); | |
59 | return -1; | |
60 | } | |
61 | ||
62 | /* copy locally a file from daq detector config db | |
63 | status=daqDA_DB_getFile("myconfig","./myconfig.txt"); | |
64 | if (status) { | |
65 | printf("Failed to get config file : %d\n",status); | |
66 | return -1; | |
67 | } | |
68 | and possibly use it */ | |
69 | ||
70 | /* report progress */ | |
71 | daqDA_progressReport(10); | |
72 | ||
73 | ||
74 | ||
75 | /* define wait event timeout - 1s max */ | |
76 | monitorSetNowait(); | |
77 | monitorSetNoWaitNetworkTimeout(1000); | |
78 | ||
79 | /* set local storage of ped files for Fe2C */ | |
80 | ||
81 | ||
82 | /* init the pedestal calculation */ | |
83 | AliHMPIDCalib *pCal=new AliHMPIDCalib(); | |
84 | /* Set the number of sigma cuts inside the file HmpidSigmaCut.txt on BOTH LDCs! */ | |
85 | /* If the file is NOT present then the default cut 3 will be used!*/ | |
86 | ||
87 | //pCal->SetSigCutFromFile("HmpidSigmaCut.txt"); | |
88 | pCal->SetSigCutFromShell("HMPID_SIGMA_CUT"); //decision to make later: wether to use file or env variable... | |
89 | pCal->SetDaOutFromShell("HMPID_DA_OUT"); //decision to make later: wether to use file or env variable... | |
90 | pCal->SetFeeInFromShell("HMPID_FEE_IN"); //decision to make later: wether to use file or env variable... | |
91 | ||
92 | ||
93 | /* ONLY set this option to kTRUE if you want to create the ADC dsitributions for all 161280 pads!!!!*/ | |
94 | /* kTRUE is not suggested for production mode b/c of the memory consumption! */ | |
95 | pCal->SetWriteHistoPads(kFALSE); //use this option for default production useage!!! | |
96 | //pCal->SetWriteHistoPads(kTRUE); //only for expert debug | |
97 | //pCal->SetWriteHistoPads(kTRUE,kTRUE,13); //DO NOT USE THIS OPTION! | |
98 | ||
99 | ||
100 | /* init event counter */ | |
101 | Int_t firstEvt=0; | |
102 | Int_t iEvtNcal=firstEvt; //Start from 1 not 0! | |
103 | ULong_t runNum=0; | |
104 | ULong_t ldcId=0; | |
105 | ||
106 | int n; | |
107 | for (n=1;n<argc;n++) { | |
108 | ||
109 | status=monitorSetDataSource( argv[n] ); | |
110 | if (status!=0) { | |
111 | printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status)); | |
112 | return -1; | |
113 | } | |
114 | ||
115 | /* report progress */ | |
116 | /* in this example, indexed on the number of files */ | |
117 | daqDA_progressReport(10+80*n/argc); | |
118 | ||
119 | for(;;) { // infinite loop | |
120 | ||
121 | /* check shutdown condition */ | |
122 | if (daqDA_checkShutdown()) {break;} | |
123 | ||
124 | struct eventHeaderStruct *event; | |
125 | eventTypeType eventT; | |
126 | ||
127 | /* get next event */ | |
128 | status=monitorGetEventDynamic((void **)&event); | |
129 | if (status==MON_ERR_EOF) /* end of monitoring file has been reached */ | |
130 | { | |
131 | printf("End of monitoring file has been reached! \n"); | |
132 | break; | |
133 | } | |
134 | ||
135 | ||
136 | if (status!=0) { | |
137 | printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status)); | |
138 | return -1; | |
139 | } | |
140 | ||
141 | /* retry if got no event */ | |
142 | if (event==NULL) { | |
143 | //break; | |
144 | continue; | |
145 | } | |
146 | ||
147 | /* use event - here, just write event id to result file */ | |
148 | eventT=event->eventType; | |
149 | ||
150 | if (eventT==PHYSICS_EVENT || eventT==CALIBRATION_EVENT ) { //updated: 18/02/2008 based on http://alice-ecs.web.cern.ch/alice-ecs/runtypes_3.16.html | |
151 | runNum=(unsigned long)event->eventRunNb; //assuming that only one run is processed at a time | |
152 | ldcId=(unsigned long)event->eventLdcId; | |
153 | if(iEvtNcal==firstEvt && pCal->GetWritePads()==kTRUE) | |
154 | { | |
155 | if(pCal->GetLargePads()==kFALSE) pCal->InitFile((Int_t)ldcId); //The number for iEvtNcal should be the same as for the first value | |
156 | else pCal->InitFile((Int_t)runNum); //The number for iEvtNcal should be the same as for the first value | |
157 | } | |
158 | ||
159 | iEvtNcal++; | |
160 | AliRawReader *reader = new AliRawReaderDate((void*)event); | |
161 | AliHMPIDRawStream stream(reader); | |
162 | while(stream.Next()) | |
163 | { | |
164 | for(Int_t iPad=0;iPad<stream.GetNPads();iPad++) { | |
165 | pCal->FillPedestal(stream.GetPadArray()[iPad],stream.GetChargeArray()[iPad]); | |
166 | } //pads | |
167 | }//while -- loop on det load in one event | |
168 | ||
169 | for(Int_t iddl=0;iddl<stream.GetNDDL();iddl++){ | |
170 | pCal->FillDDLCnt(iddl,stream.GetnDDLInStream()[iddl],stream.GetnDDLOutStream()[iddl]); | |
171 | for(Int_t ierr=0; ierr < stream.GetNErrors(); ierr++) { | |
172 | pCal->FillErrors(iddl,ierr,stream.GetErrors(iddl,ierr)); | |
173 | } | |
174 | }//err | |
175 | ||
176 | pCal->SetRunParams(runNum,stream.GetTimeStamp(),stream.GetLDCNumber()); //Get the last TimeStam read and the LDC ID | |
177 | stream.Delete(); | |
178 | }// if CALIBRATION_EVENT | |
179 | ||
180 | /* exit when last event received, no need to wait for TERM signal */ | |
181 | if (eventT==END_OF_RUN) { | |
182 | printf("EOR event detected\n"); | |
183 | break; | |
184 | ||
185 | } // events loop | |
186 | ||
187 | free(event); | |
188 | } | |
189 | ||
190 | }//arg | |
191 | ||
192 | /* write report */ | |
193 | printf("HMPID DA processed RUN #%s on LDC#%d, with %d calibration events\n",getenv("DATE_RUN_NUMBER"),ldcId,iEvtNcal); | |
194 | ||
195 | if (!iEvtNcal) { | |
196 | printf("No calibration events have been read. Exiting\n"); | |
197 | return -1; | |
198 | } | |
199 | ||
200 | /* report progress */ | |
201 | daqDA_progressReport(90); | |
202 | TString sDaOut=pCal->GetDaOutFromShell(); | |
203 | ||
204 | for(Int_t nDDL=0; nDDL < AliHMPIDRawStream::kNDDL; nDDL++) { | |
205 | ||
206 | /* Calculate pedestal for the given ddl, if there is no ddl go t next */ | |
207 | if(!pCal->CalcPedestal(nDDL,Form("%sHmpidPedDdl%02i.txt",sDaOut.Data(),nDDL),iEvtNcal)) continue; | |
208 | /* to create pedestal file as Paolo uncomment the line */ | |
209 | // if(!pCal->CalcPedestalPaolo(nDDL,Form("%sHmpidPedDdl%02i.txt",sDaOut.Data(),nDDL),iEvtNcal)) continue; | |
210 | if(!pCal->WriteErrors(nDDL,Form("%sHmpidErrorsDdl%02i.txt",sDaOut.Data(),nDDL),iEvtNcal)) continue; | |
211 | ||
212 | /* store the result file on FES */ | |
213 | ||
214 | status=daqDA_FES_storeFile(Form("%sHmpidPedDdl%02i.txt",sDaOut.Data(),nDDL),Form("HMPID_DA_Pedestals_ddl=%02i",nDDL)); | |
215 | if (status) { printf("Failed to export file : %d\n",status); } | |
216 | status=daqDA_FES_storeFile(Form("%sHmpidErrorsDdl%02i.txt",sDaOut.Data(),nDDL),Form("HMPID_DA_Errors_ddl=%02i",nDDL)); | |
217 | if (status) { printf("Failed to export file : %d\n",status); } | |
218 | ||
219 | }//nDDL | |
220 | ||
221 | if(pCal->GetWritePads()==kTRUE) pCal->CloseFile(); | |
222 | delete pCal; | |
223 | if (status) return -1; | |
224 | ||
225 | /* report progress */ | |
226 | daqDA_progressReport(100); | |
227 | ||
228 | ||
229 | return status; | |
230 | } |