]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/TOFcalibda.cxx
Coverity 10113
[u/mrichter/AliRoot.git] / TOF / TOFcalibda.cxx
1 /*
2
3 TOF DA for online calibration
4
5 Version: "1.0"
6 Contact: Roberto.Preghenella@bo.infn.it
7 Reference Run: 115401
8 Run Type: PHYSICS
9 DA Type: MON
10 Number of events needed:
11 Input Files: no input
12 Output Files: TOFdaCalib.root
13 Event types used: CALIBRATION_EVENT
14
15 */
16
17 #define FILE_CALIB "TOFdaCalib.root"
18 #define EVENT_DETECTOR_TOF 5
19
20 // DATE
21 #include "event.h"
22 #include "monitor.h"
23 #include "daqDA.h"
24
25 #include <stdio.h>
26 #include <stdlib.h>
27 #include <errno.h>
28
29 //ROOT
30 #include "TROOT.h"
31 #include "TFile.h"
32 #include "TH1F.h"
33 #include "TPluginManager.h"
34
35 //AliRoot
36 #include "AliLog.h"
37 #include "AliTOFRawStream.h"
38 #include "AliRawReaderDate.h"
39 #include "AliRawReader.h"
40 #include "AliDAQ.h"
41 #include "AliTOFGeometry.h"
42 #include "AliTOFDecoderV2.h"
43 #include "AliTOFDecoderSummaryData.h"
44 #include "AliTOFDRMSummaryData.h"
45 #include "AliTOFTRMSummaryData.h"
46 #include "AliTOFChainSummaryData.h"
47 #include "AliTOFTDCHitBuffer.h"
48 #include "AliTOFTDCHit.h"
49
50 /* Main routine
51       Arguments: 
52       1- monitoring data source
53 */
54 int 
55 main(int argc, char **argv) 
56 {
57   
58   /* magic line from Rene */
59   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
60                                         "*",
61                                         "TStreamerInfo",
62                                         "RIO",
63                                         "TStreamerInfo()");
64   
65
66   /* log start of process */
67   printf("TOF DA started\n");  
68   
69   /* check that we got some arguments = list of files */
70   if (argc!=2) {
71     printf("Wrong number of arguments\n");
72     return -1;
73   }
74
75
76   /*
77    * INIT 
78    */
79
80   /* constants */
81   const Int_t nChannels = 157248;
82   /* counters and flags */
83   Int_t nCalibEvents;
84   /* variables */
85   Int_t ddl, slot, trm, chain, tdc, channel, index, det[5], dummy;
86   /* TOF raw data handling */
87   AliTOFRawStream *rawStream = new AliTOFRawStream();
88   AliTOFDecoderV2 *decoder = rawStream->GetDecoderV2();
89   AliTOFDecoderSummaryData *decodersd;
90   AliTOFDRMSummaryData *drmsd;
91   AliTOFTRMSummaryData *trmsd;  
92   AliTOFChainSummaryData *chainsd;
93   AliTOFTDCHitBuffer *hitBuffer;
94   AliTOFTDCHit *hit;
95   UChar_t *data = 0x0;
96   Int_t dataSize;
97   Int_t dataWords;
98   Int_t currentDDL;
99
100   /* init counters and flags */
101   nCalibEvents = 0;
102   
103   /* open CALIB output file */
104   TFile *fileOutCalib = new TFile(FILE_CALIB, "RECREATE"); 
105   /* create calib hit histo */
106   TH1F *hCalibHit = new TH1F("hCalibHit", "Calibration events;index;N_{hits}/N_{events}", nChannels, 0., nChannels);
107
108   /*
109    * ONLINE MONITOR
110    */
111
112   AliLog::SetGlobalLogLevel(AliLog::kFatal);
113   struct eventHeaderStruct *event;
114   int ret;
115   /* define monitoring table */
116   char *monTable[5] = {
117     "ALL", "no",
118     "CAL", "yes",
119     NULL
120   };
121   ret = monitorDeclareTable(monTable);
122   if (ret != 0) {
123     printf("monitorDeclareTable() failed: %s\n", monitorDecodeError(ret));
124     return -1;
125   }
126   /* define data source : this is argument 1 */  
127   ret = monitorSetDataSource(argv[1]);
128   if (ret != 0) {
129     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(ret));
130     return -1;
131   }
132   /* declare monitoring program */
133   ret = monitorDeclareMp("TOFdaCalib");
134   if (ret != 0) {
135     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(ret));
136     return -1;
137   }
138   /* define wait event timeout - 1s max */
139   monitorSetNowait();
140   monitorSetNoWaitNetworkTimeout(1000);
141   
142   /* loop over events */
143   while (1) {
144     
145     /* check shutdown condition */
146     if (daqDA_checkShutdown()) break;
147     
148     /*
149      * GET EVENT
150      */
151
152     /* get next event (blocking call until timeout) */
153     ret = monitorGetEventDynamic((void **)&event);
154     if (ret == MON_ERR_EOF) {
155       printf ("End of File detected\n");
156       break; /* end of monitoring file has been reached */
157     }
158     if (ret != 0) {
159       printf("monitorGetEventDynamic() failed (ret=%d errno=%d): %s\n", ret, errno, monitorDecodeError(ret));
160       break;
161     }
162     /* retry if got no event */
163     if (event==NULL) continue;
164     /* check TOF in partecipating detectors */
165     if (!TEST_DETECTOR_IN_PATTERN(event->eventDetectorPattern, EVENT_DETECTOR_TOF)) {
166       free(event);
167       continue;
168     }
169     /* check event type */
170     if (event->eventType != CALIBRATION_EVENT) {
171       printf("not a calibration event: %d\n", event->eventType);
172       free(event);
173       continue;
174     }
175     /* increment number of calib events */
176     nCalibEvents++;
177
178     /*
179      * DECODE EVENT
180      */
181
182     /* create and setup raw reader */
183     AliRawReader *rawReader = new AliRawReaderDate((void *)event);
184     rawReader->Reset();
185     rawReader->Select("TOF", 0, AliDAQ::NumberOfDdls("TOF") - 1);
186     /* setup raw stream */
187     rawStream->SetRawReader(rawReader);
188
189     /* loop over DDLs - rawReader->ReadHeader() */
190     while (rawReader->ReadHeader()) {
191
192       /* read equipment data */
193       dataSize = rawReader->GetDataSize();
194       data = new UChar_t[dataSize];
195       if (!rawReader->ReadNext(data, dataSize)){
196         delete [] data;
197         continue;
198       }
199
200       /* decode data */
201       dataWords = dataSize / 4;
202       decoder->Decode((UInt_t *)data, dataWords);
203       delete [] data;
204
205       /* read equipment info */
206       currentDDL = rawReader->GetDDLID();
207       /* read decoder summary data */
208       decodersd = decoder->GetDecoderSummaryData();
209       /* check DRM header/trailer */
210       drmsd = decodersd->GetDRMSummaryData();
211       if (!drmsd->GetHeader() || !drmsd->GetTrailer()) continue;
212       /* loop over TRM to get hits */
213       for (Int_t itrm = 0; itrm < 10; itrm++) {
214         trmsd = drmsd->GetTRMSummaryData(itrm);
215         /* check header/trailer */
216         if (!trmsd->GetHeader() || !trmsd->GetTrailer()) continue;
217         /* loop over chains */
218         for (Int_t ichain = 0; ichain < 2; ichain++) {
219           chainsd = trmsd->GetChainSummaryData(ichain);
220           /* check header/trailer */
221           if (!chainsd->GetHeader() || !chainsd->GetTrailer()) continue;
222           hitBuffer = chainsd->GetTDCPackedHitBuffer();
223
224           /*
225            * HIT MANIPULATION
226            */
227           
228           /* loop over hits in buffer */
229           for (Int_t ihit = 0; ihit < hitBuffer->GetEntries(); ihit++) {
230             
231             /* get hit */
232             hit = hitBuffer->GetHit(ihit);
233             
234             /* get channel info */
235             ddl = currentDDL;
236             slot = trmsd->GetSlotID();
237             trm = slot - 3;
238             chain = chainsd->GetChain();
239             tdc = hit->GetTDCID();
240             channel = hit->GetChan();
241             /* get index */
242             rawStream->EquipmentId2VolumeId(ddl, slot, chain, tdc, channel, det);
243             dummy = det[4];
244             det[4] = det[3];
245             det[3] = dummy;
246             /* check valid index */
247             if (det[0] < 0 || det[0] > 17 ||
248                 det[1] < 0 || det[1] > 5 ||
249                 det[2] < 0 || det[2] > 18 ||
250                 det[3] < 0 || det[3] > 1 ||
251                 det[4] < 0 || det[4] > 47) continue;
252             index = AliTOFGeometry::GetIndex(det);
253             
254             /* fill calib hit histo */
255             hCalibHit->Fill(index);
256             
257           } /* end of loop over hits in buffer */
258         } /* end of loop over chains */
259       } /* end of loop over TRMs */
260     } /* end of loop over DDLs - rawReader->ReadHeader() */
261     
262     /* delete raw reader */
263     delete rawReader;
264     /* free event */
265     free(event);
266     
267   } /* end of loop over events */
268
269   /* scale calib hit histo by number of calib events */
270   printf("found %d calibration events\n", nCalibEvents);
271   hCalibHit->Sumw2();
272   if (nCalibEvents > 0)
273     hCalibHit->Scale(1. / nCalibEvents);
274   
275   /* write calib hit histo on CALIB file */
276   fileOutCalib->cd();
277   hCalibHit->Write();
278   fileOutCalib->Close();
279   /* export file to FXS */
280   if (daqDA_FES_storeFile(FILE_CALIB, "CALIB"))
281     return -2;
282
283   return 0;
284 }