]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/ITSSDDINJda.cxx
modifications to satisfy the coding conventions
[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: -
979b5a5f 8- Output Files: - SDDinj_ddl*c*_sid*.data
00252d8f 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
77 Int_t maxNEvents=10; // maximum number of events to be analyzed
979b5a5f 78 const Int_t kTotDDL=24;
79 const Int_t kModPerDDL=12;
80 const Int_t kSides=2;
81 AliITSOnlineSDDInjectors **injan=new AliITSOnlineSDDInjectors*[kTotDDL*kModPerDDL*kSides];
82 TH2F **histo=new TH2F*[kTotDDL*kModPerDDL*kSides];
83 Int_t nWrittenEv[kTotDDL*kModPerDDL*kSides];
84 Bool_t isFilled[kTotDDL*kModPerDDL*kSides];
4c82df4c 85
86 Char_t hisnam[20];
979b5a5f 87 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
88 for(Int_t imod=0; imod<kModPerDDL;imod++){
89 for(Int_t isid=0;isid<kSides;isid++){
90 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
91 injan[index]=new AliITSOnlineSDDInjectors(iddl,imod,isid);
92 sprintf(hisnam,"h%02dc%02ds%d",iddl,imod,isid);
93 histo[index]=new TH2F(hisnam,"",256,-0.5,255.5,256,-0.5,255.5);
94 nWrittenEv[index]=0;
95 isFilled[index]=0;
96 }
4c82df4c 97 }
98 }
99
100 /* report progress */
101 daqDA_progressReport(10);
102 Int_t iev=0;
103 /* read the data files */
104 int n;
105 for (n=1;n<argc;n++) {
106
107 status=monitorSetDataSource( argv[n] );
108 if (status!=0) {
109 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
110 return -1;
111 }
112
113 /* report progress */
114 /* in this example, indexed on the number of files */
115 // Progress report inside the event loop as well?
116 daqDA_progressReport(10+80*n/argc);
117
118 /* read the file */
119 for(;;) {
120 struct eventHeaderStruct *event;
121 eventTypeType eventT;
122
123 /* get next event */
124 status=monitorGetEventDynamic((void **)&event);
125 if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
126 if (status!=0) {
127 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
128 return -1;
129 }
130
131 /* retry if got no event */
132 if (event==NULL) {
133 break;
134 }
135
136 iev++;
137 if(iev>maxNEvents) break;
138
139 /* use event - here, just write event id to result file */
140 eventT=event->eventType;
141 switch (event->eventType){
142
143 /* START OF RUN */
144 case START_OF_RUN:
145 break;
146
147 /* END OF RUN */
148 case END_OF_RUN:
149 break;
150
151 // case PHYSICS_EVENT: // comment this line for test raw data
152 // break; // comment this line for test raw data
153
154 case CALIBRATION_EVENT:
155 break; // uncomment this line for test raw data
156 case PHYSICS_EVENT: // uncomment this line for test raw data
157 printf(" event number = %i \n",iev);
158 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
4c82df4c 159
cd452b3b 160 UInt_t timeSt=rawReader->GetTimestamp();
161 //UInt_t timeSt=iev*5000+12; // fake timestamp for test
4c82df4c 162 Int_t evtyp=0;
163 while(rawReader->ReadHeader()){
164 const UInt_t *subev = rawReader->GetSubEventAttributes();
165 if(subev[0]==0 && subev[1]==0 && subev[2]==0) evtyp=1;
166 }
167
168 rawReader->Reset();
979b5a5f 169 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
170 for(Int_t imod=0; imod<kModPerDDL;imod++){
171 for(Int_t isid=0;isid<kSides;isid++){
172 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
173 histo[index]->Reset();
174 }
4c82df4c 175 }
176 }
177 AliITSRawStreamSDD s(rawReader);
4c82df4c 178
179 while(s.Next()){
979b5a5f 180 Int_t iDDL=rawReader->GetDDLID();
181 Int_t iCarlos=s.GetCarlosId();
de075dae 182 if(s.IsCompletedModule()) continue;
183 if(s.IsCompletedDDL()) continue;
184 if(iDDL>=0 && iDDL<kTotDDL){
979b5a5f 185 Int_t index=kSides*(kModPerDDL*iDDL+iCarlos)+s.GetChannel();
4c82df4c 186 histo[index]->Fill(s.GetCoord2(),s.GetCoord1(),s.GetSignal());
979b5a5f 187 isFilled[index]=1;
4c82df4c 188 }
189 }
190 delete rawReader;
979b5a5f 191 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
192 for(Int_t imod=0; imod<kModPerDDL;imod++){
193 for(Int_t isid=0;isid<kSides;isid++){
194 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
195 if(isFilled[index]){
196 injan[index]->Reset();
197 injan[index]->AnalyzeEvent(histo[index]);
198 injan[index]->WriteToASCII(iev,timeSt,nWrittenEv[index]);
199 nWrittenEv[index]++;
200 }
201 }
4c82df4c 202 }
203 }
4c82df4c 204 /* free resources */
205 free(event);
206 }
207 }
208
209 }
210
211 Char_t filnam[100],command[120];
979b5a5f 212 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
213 for(Int_t imod=0; imod<kModPerDDL;imod++){
214 for(Int_t isid=0;isid<kSides;isid++){
215 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
216 if(nWrittenEv[index]>0){
217 sprintf(filnam,"SDDinj_ddl%02dc%02d_sid%d.data",iddl,imod,isid);
218 sprintf(command,"tar -rf SDDinj_LDC.tar %s",filnam);
219 gSystem->Exec(command);
220 }
221 }
222 }
4c82df4c 223 }
224
225 /* write report */
226 printf("Run #%s, received %d injector events\n",getenv("DATE_RUN_NUMBER"),iev);
227
228 /* report progress */
229 daqDA_progressReport(90);
230
231
2b066e71 232 status=daqDA_FES_storeFile("./SDDinj_LDC.tar","SDD_Injec");
4c82df4c 233
234 /* report progress */
235 daqDA_progressReport(100);
236
237
238
239 return status;
240}