]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/ITSSDDINJda.cxx
correct treatment of DDL header; modification needed to solve I/O problems; added...
[u/mrichter/AliRoot.git] / ITS / ITSSDDINJda.cxx
CommitLineData
00252d8f 1/*
2- Contact: - prino@to.infn.it
3- Link: -
4- Run Type: - PHYSICS
5- DA Type: - LDC
6- Number of events needed:
7- Input Files: -
8- Output Files: - SDDinj_mod*_sid*.data
9- Trigger types used:
10*/
11
12
4c82df4c 13//////////////////////////////////////////////////////////////////////////////
14// Detector Algorithm for analysis of SDD injector events. //
15// //
16// Produces ASCII output files with: //
17// 1. event number //
18// 2. event timestamp //
19// 3. Fit parameters of drift vel. vs. anode //
20// Tar Files are written to FXS //
21// //
22// Author: F. Prino (prino@to.infn.it) //
23// //
24//////////////////////////////////////////////////////////////////////////////
25
26
27
28extern "C" {
29#include "daqDA.h"
30}
31
32#include "event.h"
33#include "monitor.h"
34
35#include <stdio.h>
36#include <stdlib.h>
37
38// ROOT includes
39#include <TFile.h>
40#include <TH1F.h>
41#include <TH2F.h>
42#include <TSystem.h>
00252d8f 43#include <TROOT.h>
44#include <TPluginManager.h>
45
4c82df4c 46
47// AliRoot includes
48#include "AliRawReaderDate.h"
49#include "AliITSOnlineSDDInjectors.h"
50#include "AliITSRawStreamSDD.h"
51/* Main routine
52 Arguments: list of DATE raw data files
53*/
54int main(int argc, char **argv) {
55
56 int status = 0;
57
58
00252d8f 59 // line added to solve IO problems
60 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
61 "*",
62 "TStreamerInfo",
63 "RIO",
64 "TStreamerInfo()");
65
4c82df4c 66 /* log start of process */
67 printf("ITS SDD INJ algorithm program started\n");
68
69
70 /* check that we got some arguments = list of files */
71 if (argc<2) {
72 printf("Wrong number of arguments\n");
73 return -1;
74 }
75
76
767934e8 77 Int_t eqOffset = 256;
78 Int_t DDLrange = 24;
4c82df4c 79 Int_t maxNEvents=10; // maximum number of events to be analyzed
80 const Int_t nSDDmodules=12; // temp for test raw data
81 AliITSOnlineSDDInjectors **injan=new AliITSOnlineSDDInjectors*[2*nSDDmodules];
82 TH2F **histo=new TH2F*[2*nSDDmodules];
83 Int_t nWrittenEv[2*nSDDmodules];
84
85 Char_t hisnam[20];
86 for(Int_t imod=0; imod<nSDDmodules;imod++){
87 for(Int_t isid=0;isid<2;isid++){
88 Int_t index=2*imod+isid;
89 injan[index]=new AliITSOnlineSDDInjectors(imod,isid);
90 sprintf(hisnam,"his%03ds%d",imod,isid);
91 histo[index]=new TH2F(hisnam,"",256,-0.5,255.5,256,-0.5,255.5);
92 nWrittenEv[index]=0;
93 }
94 }
95
96 /* report progress */
97 daqDA_progressReport(10);
98 Int_t iev=0;
99 /* read the data files */
100 int n;
101 for (n=1;n<argc;n++) {
102
103 status=monitorSetDataSource( argv[n] );
104 if (status!=0) {
105 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
106 return -1;
107 }
108
109 /* report progress */
110 /* in this example, indexed on the number of files */
111 // Progress report inside the event loop as well?
112 daqDA_progressReport(10+80*n/argc);
113
114 /* read the file */
115 for(;;) {
116 struct eventHeaderStruct *event;
117 eventTypeType eventT;
118
119 /* get next event */
120 status=monitorGetEventDynamic((void **)&event);
121 if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
122 if (status!=0) {
123 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
124 return -1;
125 }
126
127 /* retry if got no event */
128 if (event==NULL) {
129 break;
130 }
131
132 iev++;
133 if(iev>maxNEvents) break;
134
135 /* use event - here, just write event id to result file */
136 eventT=event->eventType;
137 switch (event->eventType){
138
139 /* START OF RUN */
140 case START_OF_RUN:
141 break;
142
143 /* END OF RUN */
144 case END_OF_RUN:
145 break;
146
147 // case PHYSICS_EVENT: // comment this line for test raw data
148 // break; // comment this line for test raw data
149
150 case CALIBRATION_EVENT:
151 break; // uncomment this line for test raw data
152 case PHYSICS_EVENT: // uncomment this line for test raw data
153 printf(" event number = %i \n",iev);
154 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
767934e8 155 rawReader->SelectEquipment(17,eqOffset+1,eqOffset+DDLrange);
4c82df4c 156
cd452b3b 157 UInt_t timeSt=rawReader->GetTimestamp();
158 //UInt_t timeSt=iev*5000+12; // fake timestamp for test
4c82df4c 159 Int_t evtyp=0;
160 while(rawReader->ReadHeader()){
161 const UInt_t *subev = rawReader->GetSubEventAttributes();
162 if(subev[0]==0 && subev[1]==0 && subev[2]==0) evtyp=1;
163 }
164
165 rawReader->Reset();
166 for(Int_t imod=0; imod<nSDDmodules;imod++){
167 for(Int_t isid=0; isid<2;isid++){
168 Int_t index=2*imod+isid;
169 histo[index]->Reset();
170 }
171 }
172 AliITSRawStreamSDD s(rawReader);
4c82df4c 173
174 while(s.Next()){
175 Int_t iddl=rawReader->GetDDLID();
176 iddl=0; // temporary for test raw data
177 Int_t isddmod=s.GetModuleNumber(iddl,s.GetCarlosId());
178 isddmod-=240; // to have SDD modules from 0 to 259
179 isddmod=s.GetCarlosId(); // temporary for test raw data
5dfa68c5 180 if(isddmod<nSDDmodules&&s.IsCompletedModule()==kFALSE){
4c82df4c 181 Int_t index=2*isddmod+s.GetChannel();
182 histo[index]->Fill(s.GetCoord2(),s.GetCoord1(),s.GetSignal());
183 }
184 }
185 delete rawReader;
4c82df4c 186 for(Int_t imod=0; imod<nSDDmodules;imod++){
187 for(Int_t isid=0; isid<2;isid++){
188 Int_t index=2*imod+isid;
189 injan[index]->Reset();
190 injan[index]->AnalyzeEvent(histo[index]);
191 injan[index]->WriteToASCII(iev,timeSt,nWrittenEv[index]);
192 nWrittenEv[index]++;
193 }
194 }
195
196 /* free resources */
197 free(event);
198 }
199 }
200
201 }
202
203 Char_t filnam[100],command[120];
204 for(Int_t imod=0; imod<nSDDmodules;imod++){
205 for(Int_t isid=0; isid<2;isid++){
206 sprintf(filnam,"SDDinj_mod%03d_sid%d.data",imod,isid);
207 sprintf(command,"tar -rf SDDinj_LDC1.tar %s",filnam);
208 gSystem->Exec(command);
209 }
210 }
211
212 /* write report */
213 printf("Run #%s, received %d injector events\n",getenv("DATE_RUN_NUMBER"),iev);
214
215 /* report progress */
216 daqDA_progressReport(90);
217
218
219 status=daqDA_FES_storeFile("./SDDinj_LDC1.tar","SDD_Calib");
220
221 /* report progress */
222 daqDA_progressReport(100);
223
224
225
226 return status;
227}