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