]>
Commit | Line | Data |
---|---|---|
1 | /* | |
2 | ||
3 | TOF DA for pulser data | |
4 | ||
5 | Contact: Chiara.Zampolli@bo.infn.it | |
6 | Link: www.bo.infn.it/~zampolli | |
7 | Run Type: PULSER | |
8 | DA Type: LDC | |
9 | Number of events needed: 10000 | |
10 | Input Files: TOF<nrun>.raw, where <nrun> is the run number | |
11 | Output Files: TOFoutPulserLDC_<LDCid>.root, where <LDCid> is the id of the LDC which is read (2 characters field, e.g. TOFoutPulserLDC_03.root),to be exported to the DAQ FXS | |
12 | Trigger types used: PHYSICS_EVENT (for the time being) | |
13 | ||
14 | */ | |
15 | ||
16 | // DATE | |
17 | #include "event.h" | |
18 | #include "monitor.h" | |
19 | #include "daqDA.h" | |
20 | ||
21 | #include <stdio.h> | |
22 | #include <stdlib.h> | |
23 | ||
24 | //AliRoot | |
25 | #include <AliTOFRawStream.h> | |
26 | #include <AliRawReaderDate.h> | |
27 | #include <AliRawReader.h> | |
28 | #include <AliTOFGeometry.h> | |
29 | #include <AliDAQ.h> | |
30 | #include <AliTOFHitData.h> | |
31 | #include <AliTOFHitDataBuffer.h> | |
32 | #include <AliTOFDecoder.h> | |
33 | #include <AliTOFNoiseConfigHandler.h> | |
34 | ||
35 | //ROOT | |
36 | #include <TFile.h> | |
37 | #include <TKey.h> | |
38 | #include <TH1I.h> | |
39 | #include <TObject.h> | |
40 | #include <TMath.h> | |
41 | #include <TSystem.h> | |
42 | #include "TROOT.h" | |
43 | #include "TPluginManager.h" | |
44 | #include "TSAXParser.h" | |
45 | ||
46 | /* Main routine | |
47 | Arguments: list of DATE raw data files | |
48 | */ | |
49 | int main(int argc, char **argv) { | |
50 | ||
51 | /* magic line from Rene */ | |
52 | gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo", | |
53 | "*", | |
54 | "TStreamerInfo", | |
55 | "RIO", | |
56 | "TStreamerInfo()"); | |
57 | ||
58 | AliTOFGeometry * geom = new AliTOFGeometry(); | |
59 | ||
60 | static const Int_t size = AliTOFGeometry::NPadXSector()*AliTOFGeometry::NSectors(); | |
61 | TH1F::AddDirectory(0); | |
62 | TH1I * htofPulser = new TH1I("hTOFpulser","histo with signals on TOF during pulser", size,-0.5,size-0.5); | |
63 | for (Int_t ibin =1;ibin<=size;ibin++) | |
64 | htofPulser->SetBinContent(ibin,-1); | |
65 | ||
66 | UInt_t ldcId=99; | |
67 | UInt_t ldcIdOLD=99; | |
68 | ||
69 | int status; | |
70 | ||
71 | /* log start of process */ | |
72 | printf("TOF DA started\n"); | |
73 | ||
74 | /* check that we got some arguments = list of files */ | |
75 | if (argc<2) { | |
76 | printf("Wrong number of arguments\n"); | |
77 | return -1; | |
78 | } | |
79 | ||
80 | /* retrieve config file */ | |
81 | int getConfigFile = daqDA_DB_getFile("TOFNoiseConfig.xml","TOFNoiseConfig.xml"); | |
82 | if (getConfigFile != 0){ | |
83 | printf("Failed to retrieve config file from DB! returning...\n"); | |
84 | return -1; | |
85 | } | |
86 | ||
87 | AliTOFNoiseConfigHandler* tofHandler = new AliTOFNoiseConfigHandler(); | |
88 | TSAXParser *parser = new TSAXParser(); | |
89 | parser->ConnectToHandler("AliTOFNoiseConfigHandler", tofHandler); | |
90 | if (parser->ParseFile("./TOFNoiseConfig.xml") != 0) { | |
91 | printf("Failed parsing config file! retunring... \n"); | |
92 | return -1; | |
93 | } | |
94 | ||
95 | Int_t debugFlag = tofHandler->GetDebugFlag(); | |
96 | printf("the debug flag is %i\n",debugFlag); | |
97 | ||
98 | /* init some counters */ | |
99 | int nevents_physics=0; | |
100 | int nevents_total=0; | |
101 | ||
102 | Int_t nPDBEntriesToT = 0; | |
103 | Int_t nDBEntriesToT = 0; | |
104 | AliTOFHitData *HitData = 0; | |
105 | Int_t dummy = -1; | |
106 | Int_t Volume[5]; | |
107 | for (Int_t i=0;i<5;i++) Volume[i]=-1; | |
108 | AliTOFRawStream *rawStreamTOF = new AliTOFRawStream(); | |
109 | AliTOFDecoder * decoderTOF = new AliTOFDecoder(); | |
110 | AliTOFHitDataBuffer *DataBuffer[AliDAQ::NumberOfDdls("TOF")]; | |
111 | AliTOFHitDataBuffer *PackedDataBuffer[AliDAQ::NumberOfDdls("TOF")]; | |
112 | for (Int_t i=0;i<AliDAQ::NumberOfDdls("TOF");i++) { | |
113 | DataBuffer[i]=new AliTOFHitDataBuffer(); | |
114 | PackedDataBuffer[i]=new AliTOFHitDataBuffer(); | |
115 | } | |
116 | Int_t currentEquipment; | |
117 | Int_t currentDDL; | |
118 | UChar_t *data = 0x0; | |
119 | Int_t nchDDL = 0; | |
120 | Int_t nDBEntries = 0; | |
121 | Int_t nPDBEntries = 0; | |
122 | ||
123 | struct eventHeaderStruct *event; | |
124 | eventTypeType eventT; | |
125 | ||
126 | /* read the data files */ | |
127 | int n; | |
128 | for (n=1;n<argc;n++) { | |
129 | ||
130 | status=monitorSetDataSource( argv[n] ); | |
131 | if (status!=0) { | |
132 | printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status)); | |
133 | return -1; | |
134 | } | |
135 | ||
136 | /* read the file */ | |
137 | for(;;) { | |
138 | ||
139 | rawStreamTOF->Clear(); | |
140 | currentEquipment = 0; | |
141 | currentDDL = 0; | |
142 | data = 0x0; | |
143 | nPDBEntriesToT = 0; | |
144 | nDBEntriesToT = 0; | |
145 | dummy = -1; | |
146 | nchDDL = 0; | |
147 | nDBEntries = 0; | |
148 | nPDBEntries = 0; | |
149 | HitData = 0; | |
150 | ||
151 | /* get next event */ | |
152 | status=monitorGetEventDynamic((void **)&event); | |
153 | if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */ | |
154 | if (status!=0) { | |
155 | printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status)); | |
156 | return -1; | |
157 | } | |
158 | ||
159 | /* retry if got no event */ | |
160 | if (event==NULL) | |
161 | break; | |
162 | ||
163 | /* use event - here, just write event id to result file */ | |
164 | eventT=event->eventType; | |
165 | ||
166 | if (eventT==PHYSICS_EVENT) { | |
167 | //printf ("event %i \n", nevents_physics); | |
168 | ||
169 | AliRawReader *rawReader = new AliRawReaderDate((void*)event); | |
170 | rawStreamTOF->SetRawReader(rawReader); | |
171 | ||
172 | while (rawReader->ReadHeader()) { | |
173 | ldcId = rawReader->GetLDCId(); | |
174 | //debugging printings | |
175 | if (debugFlag && ldcId!=ldcIdOLD) { | |
176 | ldcIdOLD = ldcId; | |
177 | printf ("ldcId = %i \n",ldcId); | |
178 | } | |
179 | ||
180 | //memory leak prevention (actually data should be always 0x0 here) | |
181 | if (data != 0x0) delete [] data; | |
182 | ||
183 | //get equipment infos | |
184 | currentEquipment = rawReader->GetEquipmentId(); | |
185 | currentDDL = rawReader->GetDDLID(); | |
186 | const AliRawDataHeader *currentCDH = (AliRawDataHeader*)rawReader->GetDataHeader(); | |
187 | if (currentDDL%2==0) | |
188 | nchDDL = 2160; | |
189 | else | |
190 | nchDDL = 2208; | |
191 | ||
192 | //Int_t * array = new Int_t[nchDDL]; | |
193 | Int_t array[nchDDL]; | |
194 | for (Int_t ii=0; ii<nchDDL; ii++) array[ii] = 0; | |
195 | decoderTOF->GetArrayDDL(array, currentDDL); | |
196 | ||
197 | for (Int_t i=0;i<nchDDL;i++) | |
198 | if (htofPulser->GetBinContent(array[i]+1)<0) htofPulser->SetBinContent(array[i]+1,0); | |
199 | ||
200 | //debugging printings | |
201 | //if (debugFlag) printf(" Equipment = %i, and DDL = %i \n", currentEquipment,currentDDL); | |
202 | const Int_t kDataSize = rawReader->GetDataSize(); | |
203 | const Int_t kDataWords = kDataSize / 4; | |
204 | data = new UChar_t[kDataSize]; | |
205 | decoderTOF->SetDataBuffer(DataBuffer[currentDDL]); | |
206 | decoderTOF->SetPackedDataBuffer(PackedDataBuffer[currentDDL]); | |
207 | //start decoding | |
208 | if (!rawReader->ReadNext(data, kDataSize)) { | |
209 | rawReader->AddMajorErrorLog(AliTOFRawStream::kDDLdataReading); | |
210 | printf("Error while reading DDL data. Go to next equipment \n"); | |
211 | delete [] data; | |
212 | data = 0x0; | |
213 | continue; | |
214 | } | |
215 | if (decoderTOF->Decode((UInt_t *)data, kDataWords, currentCDH) == kTRUE) { | |
216 | rawReader->AddMajorErrorLog(AliTOFRawStream::kDDLDecoder,Form("DDL # = %d",currentDDL)); | |
217 | printf("Error while decoding DDL # %d: decoder returned with errors \n", currentDDL); | |
218 | } | |
219 | ||
220 | nDBEntries = DataBuffer[currentDDL]->GetEntries(); | |
221 | nPDBEntries = PackedDataBuffer[currentDDL]->GetEntries(); | |
222 | nPDBEntriesToT+=nPDBEntries; | |
223 | nDBEntriesToT+=nDBEntries; | |
224 | /* reset buffer */ | |
225 | DataBuffer[currentDDL]->Reset(); | |
226 | ||
227 | /* read data buffer hits */ | |
228 | for (Int_t iHit = 0; iHit < nPDBEntries; iHit++){ | |
229 | HitData = PackedDataBuffer[currentDDL]->GetHit(iHit); | |
230 | /* add volume information */ | |
231 | HitData->SetDDLID(currentDDL); | |
232 | for (Int_t i=0;i<5;i++) Volume[i]=-1; | |
233 | rawStreamTOF->EquipmentId2VolumeId(HitData, Volume); | |
234 | if (Volume[0]==-1 || | |
235 | Volume[1]==-1 || | |
236 | Volume[2]==-1 || | |
237 | Volume[3]==-1 || | |
238 | Volume[4]==-1) continue; | |
239 | else { | |
240 | dummy = Volume[3]; | |
241 | Volume[3] = Volume[4]; | |
242 | Volume[4] = dummy; | |
243 | Int_t index = geom->GetIndex(Volume); | |
244 | // to check array indexes | |
245 | /* | |
246 | Bool_t found =kFALSE; | |
247 | for (Int_t j=0;j<nchDDL;j++){ | |
248 | if (index==array[j]) { | |
249 | found = kTRUE; | |
250 | break; | |
251 | } | |
252 | } | |
253 | printf ("index = %6d, found = %6d\n",index, (Int_t)found); | |
254 | */ | |
255 | //printf ("index = %6d \n",index); | |
256 | htofPulser->Fill(index); //channel index start from 0, bin index from 1 | |
257 | //debugging printings | |
258 | //if (debugFlag) printf("sector %2d, plate %1d, strip %2d, padz %1d, padx %2d \n",Volume[0],Volume[1],Volume[2],Volume[3],Volume[4]); // too verbose | |
259 | } | |
260 | } | |
261 | /* reset buffer */ | |
262 | PackedDataBuffer[currentDDL]->Reset(); | |
263 | delete [] data; | |
264 | data = 0x0; | |
265 | //delete [] array; | |
266 | } | |
267 | ||
268 | //debugging printings | |
269 | //if (debugFlag) { | |
270 | // printf(" Packed Hit Buffer Entries = %i \n",nPDBEntriesToT); // too verbose | |
271 | // printf(" Hit Buffer Entries = %i \n",nDBEntriesToT); // too verbose | |
272 | //} | |
273 | ||
274 | delete rawReader; | |
275 | rawReader = 0x0; | |
276 | nevents_physics++; | |
277 | } | |
278 | nevents_total++; | |
279 | ||
280 | /* free resources */ | |
281 | free(event); | |
282 | } | |
283 | ||
284 | } | |
285 | ||
286 | delete decoderTOF; | |
287 | decoderTOF=0x0; | |
288 | ||
289 | for (Int_t i=0;i<AliDAQ::NumberOfDdls("TOF");i++){ | |
290 | delete DataBuffer[i]; | |
291 | delete PackedDataBuffer[i]; | |
292 | } | |
293 | ||
294 | delete rawStreamTOF; | |
295 | rawStreamTOF = 0x0; | |
296 | ||
297 | delete geom; | |
298 | geom = 0x0; | |
299 | ||
300 | //write the Run level file | |
301 | char filename[100]; | |
302 | sprintf(filename,"TOFoutPulserLDC_%02i.root",ldcId); | |
303 | TFile * fileRun = new TFile (filename,"RECREATE"); | |
304 | htofPulser->Write(); | |
305 | fileRun->Close(); | |
306 | ||
307 | /* write report */ | |
308 | printf("Run #%s, received %d physics events out of %d\n", | |
309 | getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total); | |
310 | ||
311 | status = 0; | |
312 | ||
313 | /* store the result file on FES */ | |
314 | status=daqDA_FES_storeFile(filename,"PULSER"); | |
315 | if (status) | |
316 | status = -2; | |
317 | ||
318 | return status; | |
319 | } |