]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/TOFpulserda.cxx
Adding the track fit residuals as a consequence of the ExB distortions (Marian)
[u/mrichter/AliRoot.git] / TOF / TOFpulserda.cxx
1 /*
2
3 TOF DA for online calibration from 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
34 //ROOT
35 #include <TFile.h>
36 #include <TKey.h>
37 #include <TH1I.h>
38 #include <TObject.h>
39 #include <TMath.h>
40 #include <TSystem.h>
41 #include "TROOT.h"
42 #include "TPluginManager.h"
43
44 /* Main routine
45       Arguments: list of DATE raw data files
46 */
47 int main(int argc, char **argv) {
48
49 /* magic line from Rene */
50   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
51                     "*",
52                     "TStreamerInfo",
53                     "RIO",
54                     "TStreamerInfo()");
55
56   AliTOFGeometry * geom = new AliTOFGeometry();
57
58   static const Int_t size = AliTOFGeometry::NPadXSector()*AliTOFGeometry::NSectors();
59   TH1F::AddDirectory(0);
60   TH1I * htofPulser = new TH1I("hTOFpulser","histo with signals on TOF during pulser", size,-0.5,size-0.5);
61   for (Int_t ibin =1;ibin<=size;ibin++){
62     htofPulser->SetBinContent(ibin,-1);
63   }
64
65   UInt_t ldcId=0;
66
67   int status;
68
69   /* log start of process */
70   printf("TOF DA started\n");  
71
72   /* check that we got some arguments = list of files */
73   if (argc<2) {
74     printf("Wrong number of arguments\n");
75     return -1;
76   }
77
78   /* open result file */
79   FILE *fp=NULL;
80   fp=fopen("./result.txt","a");
81   if (fp==NULL) {
82     printf("Failed to open file\n");
83     return -1;
84   }
85   
86   /* init some counters */
87   int nevents_physics=0;
88   int nevents_total=0;
89
90   /* read the data files */
91   int n;
92   for (n=1;n<argc;n++) {
93    
94     status=monitorSetDataSource( argv[n] );
95     if (status!=0) {
96       printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
97       return -1;
98     }
99
100     /* read the file */
101     for(;;) {
102       struct eventHeaderStruct *event;
103       eventTypeType eventT;
104
105       /* get next event */
106       status=monitorGetEventDynamic((void **)&event);
107       if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
108       if (status!=0) {
109         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
110         return -1;
111       }
112
113       /* retry if got no event */
114       if (event==NULL) {
115         break;
116       }
117
118       /* use event - here, just write event id to result file */
119       eventT=event->eventType;
120
121       if (eventT==PHYSICS_EVENT) {
122         //printf ("event %i \n", nevents_physics);
123         Int_t nPDBEntriesToT = 0;
124         Int_t nDBEntriesToT = 0;
125         AliTOFHitData *HitData;
126         Int_t dummy = -1;
127         Int_t Volume[5];
128         AliRawReader *rawReader = new AliRawReaderDate((void*)event);
129         AliTOFRawStream *rawStreamTOF = new AliTOFRawStream(rawReader);
130         AliTOFDecoder * decoderTOF = new AliTOFDecoder();
131         AliTOFHitDataBuffer *DataBuffer[72]; 
132         AliTOFHitDataBuffer *PackedDataBuffer[72];      
133         for (Int_t i=0;i<AliDAQ::NumberOfDdls("TOF");i++){
134           DataBuffer[i]=new AliTOFHitDataBuffer();
135           PackedDataBuffer[i]=new AliTOFHitDataBuffer();
136         }
137         Int_t currentEquipment;
138         Int_t currentDDL;
139         UChar_t *data = 0x0;
140         while(rawReader->ReadHeader()){
141           ldcId = rawReader->GetLDCId(); 
142           //printf ("ldcId = %i \n",ldcId);
143
144           //memory leak prevention (actually data should be always 0x0 here)
145           if (data != 0x0)delete [] data;
146           
147           //get equipment infos
148           currentEquipment = rawReader->GetEquipmentId();
149           currentDDL = rawReader->GetDDLID();
150           Int_t nchDDL = 0;
151           if (currentDDL%2==0) {
152             nchDDL = 2160;
153           }
154           else {
155             nchDDL = 2208;
156           }
157           Int_t * array = new Int_t[nchDDL];
158           decoderTOF->GetArrayDDL(array, currentDDL);
159           
160           for (Int_t i=0;i<nchDDL;i++){
161             if (htofPulser->GetBinContent(array[i]+1)<0) htofPulser->SetBinContent(array[i]+1,0);
162           }
163           
164           //printf(" Equipment = %i, and DDL = %i \n", currentEquipment,currentDDL); 
165           const Int_t kDataSize = rawReader->GetDataSize();
166           const Int_t kDataWords = kDataSize / 4;
167           data = new UChar_t[kDataSize];
168           decoderTOF->SetDataBuffer(DataBuffer[currentDDL]);
169           decoderTOF->SetPackedDataBuffer(PackedDataBuffer[currentDDL]);
170           //start decoding
171           if (!rawReader->ReadNext(data, kDataSize)){
172             rawReader->AddMajorErrorLog(AliTOFRawStream::kDDLdataReading);
173             printf("Error while reading DDL data. Go to next equipment \n");
174             delete [] data;
175             data = 0x0;
176             continue;
177           }
178           if (decoderTOF->Decode((UInt_t *)data, kDataWords) == kTRUE) {
179             rawReader->AddMajorErrorLog(AliTOFRawStream::kDDLDecoder,Form("DDL # = %d",currentDDL));
180             printf("Error while decoding DDL # %d: decoder returned with errors \n", currentDDL);
181           }
182           
183           Int_t nDBEntries = DataBuffer[currentDDL]->GetEntries();
184           Int_t nPDBEntries = PackedDataBuffer[currentDDL]->GetEntries();
185           nPDBEntriesToT+=nPDBEntries;
186           nDBEntriesToT+=nDBEntries;
187           /* reset buffer */
188           DataBuffer[currentDDL]->Reset();
189           
190           /* read data buffer hits */
191           for (Int_t iHit = 0; iHit < nPDBEntries; iHit++){
192             HitData = PackedDataBuffer[currentDDL]->GetHit(iHit);
193             /* add volume information */
194             HitData->SetDDLID(currentDDL);
195             rawStreamTOF->EquipmentId2VolumeId(HitData, Volume);
196             if (Volume[0]==-1 ||
197                 Volume[1]==-1 ||
198                 Volume[2]==-1 ||
199                 Volume[3]==-1 ||
200                 Volume[4]==-1) continue;
201             else {
202               dummy = Volume[3];
203               Volume[3] = Volume[4];
204               Volume[4] = dummy;
205               Int_t index = geom->GetIndex(Volume);
206               //printf ("index = %i \n",index);
207               htofPulser->Fill(index); //channel index start from 0, bin index from 1
208               //debugging printings
209               //printf("sector %i, plate %i, strip %i, padz %i, padx %i \n",Volume[0],Volume[1],Volume[2],Volume[3],Volume[4]);
210             }
211           }
212           /* reset buffer */
213           PackedDataBuffer[currentDDL]->Reset();
214           delete [] data;
215           data = 0x0;
216           delete [] array;
217         }
218         
219         //printf(" Packed Hit Buffer Entries = %i \n",nPDBEntriesToT);
220         //printf(" Hit Buffer Entries = %i \n",nDBEntriesToT);
221         
222         delete decoderTOF;
223         decoderTOF=0x0;
224
225         for (Int_t i=0;i<72;i++){ 
226           delete DataBuffer[i];
227           delete PackedDataBuffer[i];
228         }
229         
230         delete rawStreamTOF;
231         rawStreamTOF = 0x0;
232
233         delete rawReader;
234         rawReader = 0x0;
235         nevents_physics++;
236       }
237       nevents_total++;
238
239       /* free resources */
240       free(event);
241     }
242   }
243
244   delete geom;
245   geom = 0x0;
246
247   //write the Run level file   
248   char filename[100];
249   sprintf(filename,"TOFoutPulserLDC_%02i.root",ldcId);
250   TFile * fileRun = new TFile (filename,"RECREATE"); 
251   htofPulser->Write();
252   fileRun->Close();
253   
254   /* write report */
255   fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
256
257   /* close result file */
258   fclose(fp);
259   
260   /* store the result file on FES */
261   status=daqDA_FES_storeFile(filename,"PULSER");
262   if (status) {
263     status = -2;
264   }
265
266   return status;
267 }