]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TOF/TOFpulserda.cxx
Do not get main runloader from gAlice
[u/mrichter/AliRoot.git] / TOF / TOFpulserda.cxx
CommitLineData
6b034754 1/*
2
3TOF DA for online calibration from pulser data
4
5Contact: Chiara.Zampolli@bo.infn.it
6Link: www.bo.infn.it/~zampolli
7fffa85b 7Run Type: PULSER
6b034754 8DA Type: LDC
9Number of events needed: 10000
10Input Files: TOF<nrun>.raw, where <nrun> is the run number
11Output 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
12Trigger types used: PHYSICS_EVENT, for the time being
13
14*/
15
16// DATE
17#include "event.h"
18#include "monitor.h"
89b3ea8c 19#include "daqDA.h"
6b034754 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>
3bc66139 41#include "TROOT.h"
42#include "TPluginManager.h"
6b034754 43
44/* Main routine
45 Arguments: list of DATE raw data files
46*/
47int main(int argc, char **argv) {
48
3bc66139 49/* magic line from Rene */
50 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
51 "*",
52 "TStreamerInfo",
53 "RIO",
54 "TStreamerInfo()");
55
6b034754 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);
9512a964 130 const AliRawDataHeader *currentCDH;
131 AliTOFDecoder * decoderTOF = new AliTOFDecoder();
6b034754 132 AliTOFHitDataBuffer *DataBuffer[72];
133 AliTOFHitDataBuffer *PackedDataBuffer[72];
134 for (Int_t i=0;i<AliDAQ::NumberOfDdls("TOF");i++){
135 DataBuffer[i]=new AliTOFHitDataBuffer();
136 PackedDataBuffer[i]=new AliTOFHitDataBuffer();
137 }
138 Int_t currentEquipment;
139 Int_t currentDDL;
140 UChar_t *data = 0x0;
141 while(rawReader->ReadHeader()){
142 ldcId = rawReader->GetLDCId();
143 //printf ("ldcId = %i \n",ldcId);
144
145 //memory leak prevention (actually data should be always 0x0 here)
146 if (data != 0x0)delete [] data;
147
148 //get equipment infos
149 currentEquipment = rawReader->GetEquipmentId();
150 currentDDL = rawReader->GetDDLID();
9512a964 151 currentCDH = rawReader->GetDataHeader();
6b034754 152 Int_t nchDDL = 0;
153 if (currentDDL%2==0) {
154 nchDDL = 2160;
155 }
156 else {
157 nchDDL = 2208;
158 }
159 Int_t * array = new Int_t[nchDDL];
160 decoderTOF->GetArrayDDL(array, currentDDL);
161
162 for (Int_t i=0;i<nchDDL;i++){
163 if (htofPulser->GetBinContent(array[i]+1)<0) htofPulser->SetBinContent(array[i]+1,0);
164 }
165
166 //printf(" Equipment = %i, and DDL = %i \n", currentEquipment,currentDDL);
167 const Int_t kDataSize = rawReader->GetDataSize();
168 const Int_t kDataWords = kDataSize / 4;
169 data = new UChar_t[kDataSize];
170 decoderTOF->SetDataBuffer(DataBuffer[currentDDL]);
171 decoderTOF->SetPackedDataBuffer(PackedDataBuffer[currentDDL]);
172 //start decoding
173 if (!rawReader->ReadNext(data, kDataSize)){
174 rawReader->AddMajorErrorLog(AliTOFRawStream::kDDLdataReading);
175 printf("Error while reading DDL data. Go to next equipment \n");
176 delete [] data;
177 data = 0x0;
178 continue;
179 }
9512a964 180 if (decoderTOF->Decode((UInt_t *)data, kDataWords, currentCDH) == kTRUE) {
6b034754 181 rawReader->AddMajorErrorLog(AliTOFRawStream::kDDLDecoder,Form("DDL # = %d",currentDDL));
182 printf("Error while decoding DDL # %d: decoder returned with errors \n", currentDDL);
183 }
184
185 Int_t nDBEntries = DataBuffer[currentDDL]->GetEntries();
186 Int_t nPDBEntries = PackedDataBuffer[currentDDL]->GetEntries();
187 nPDBEntriesToT+=nPDBEntries;
188 nDBEntriesToT+=nDBEntries;
189 /* reset buffer */
190 DataBuffer[currentDDL]->Reset();
191
192 /* read data buffer hits */
193 for (Int_t iHit = 0; iHit < nPDBEntries; iHit++){
194 HitData = PackedDataBuffer[currentDDL]->GetHit(iHit);
195 /* add volume information */
196 HitData->SetDDLID(currentDDL);
197 rawStreamTOF->EquipmentId2VolumeId(HitData, Volume);
198 if (Volume[0]==-1 ||
199 Volume[1]==-1 ||
200 Volume[2]==-1 ||
201 Volume[3]==-1 ||
202 Volume[4]==-1) continue;
203 else {
204 dummy = Volume[3];
205 Volume[3] = Volume[4];
206 Volume[4] = dummy;
207 Int_t index = geom->GetIndex(Volume);
208 //printf ("index = %i \n",index);
209 htofPulser->Fill(index); //channel index start from 0, bin index from 1
210 //debugging printings
211 //printf("sector %i, plate %i, strip %i, padz %i, padx %i \n",Volume[0],Volume[1],Volume[2],Volume[3],Volume[4]);
212 }
213 }
214 /* reset buffer */
215 PackedDataBuffer[currentDDL]->Reset();
216 delete [] data;
217 data = 0x0;
218 delete [] array;
219 }
220
221 //printf(" Packed Hit Buffer Entries = %i \n",nPDBEntriesToT);
222 //printf(" Hit Buffer Entries = %i \n",nDBEntriesToT);
223
224 delete decoderTOF;
225 decoderTOF=0x0;
226
227 for (Int_t i=0;i<72;i++){
228 delete DataBuffer[i];
229 delete PackedDataBuffer[i];
230 }
231
232 delete rawStreamTOF;
233 rawStreamTOF = 0x0;
234
235 delete rawReader;
236 rawReader = 0x0;
237 nevents_physics++;
238 }
239 nevents_total++;
240
241 /* free resources */
242 free(event);
243 }
9512a964 244
6b034754 245 }
246
247 delete geom;
248 geom = 0x0;
249
250 //write the Run level file
251 char filename[100];
252 sprintf(filename,"TOFoutPulserLDC_%02i.root",ldcId);
253 TFile * fileRun = new TFile (filename,"RECREATE");
254 htofPulser->Write();
255 fileRun->Close();
256
257 /* write report */
258 fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
259
260 /* close result file */
261 fclose(fp);
262
263 /* store the result file on FES */
98a837bf 264 status=daqDA_FES_storeFile(filename,"PULSER");
6b034754 265 if (status) {
266 status = -2;
267 }
268
269 return status;
270}