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