]>
Commit | Line | Data |
---|---|---|
1 | /* | |
2 | - Contact: - prino@to.infn.it | |
3 | - Link: - alien:///alice/data/2009/LHC09c_SDD/000079098/raw/09000079098024.10.root | |
4 | - Run Type: - INJECTOR | |
5 | - DA Type: - LDC | |
6 | - Number of events needed: >20 | |
7 | - Input Files: - | |
8 | - Output Files: - SDDinj_ddl*c*_sid*.data | |
9 | - Trigger types used: | |
10 | */ | |
11 | ||
12 | ||
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 | ||
28 | extern "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> | |
43 | #include <TROOT.h> | |
44 | #include <TPluginManager.h> | |
45 | #include <TObjArray.h> | |
46 | #include <TObjString.h> | |
47 | #include <TDatime.h> | |
48 | ||
49 | ||
50 | // AliRoot includes | |
51 | #include "AliRawReaderDate.h" | |
52 | #include "AliITSOnlineSDDInjectors.h" | |
53 | #include "AliITSRawStreamSDD.h" | |
54 | #include "AliITSRawStreamSDDCompressed.h" | |
55 | ||
56 | #ifdef ALI_AMORE | |
57 | #include <AmoreDA.h> | |
58 | #endif | |
59 | ||
60 | /* Main routine | |
61 | Arguments: list of DATE raw data files | |
62 | */ | |
63 | int main(int argc, char **argv) { | |
64 | // main - Arguments: list of DATE raw data files | |
65 | int status = 0; | |
66 | ||
67 | ||
68 | // line added to solve IO problems | |
69 | gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo", | |
70 | "*", | |
71 | "TStreamerInfo", | |
72 | "RIO", | |
73 | "TStreamerInfo()"); | |
74 | ||
75 | /* log start of process */ | |
76 | printf("ITS SDD INJ algorithm program started\n"); | |
77 | ||
78 | ||
79 | /* check that we got some arguments = list of files */ | |
80 | if (argc<2) { | |
81 | printf("Wrong number of arguments\n"); | |
82 | return -1; | |
83 | } | |
84 | ||
85 | ||
86 | Int_t maxNEvents=20; // maximum number of events to be analyzed | |
87 | const Int_t kTotDDL=24; | |
88 | const Int_t kModPerDDL=12; | |
89 | const Int_t kSides=2; | |
90 | Int_t amSamplFreq=40; | |
91 | ||
92 | gSystem->Exec("rm -f SDDinj_ddl*.data"); | |
93 | gSystem->Exec("rm -f SDDinj_LDC.tar"); | |
94 | ||
95 | ||
96 | ||
97 | AliITSOnlineSDDInjectors **injan=new AliITSOnlineSDDInjectors*[kTotDDL*kModPerDDL*kSides]; | |
98 | TH2F **histo=new TH2F*[kTotDDL*kModPerDDL*kSides]; | |
99 | Int_t nWrittenEv[kTotDDL*kModPerDDL*kSides]; | |
100 | Bool_t isFilled[kTotDDL*kModPerDDL*kSides]; | |
101 | Bool_t writtenoutput=kFALSE; | |
102 | UInt_t timeSt=0; | |
103 | Char_t hisnam[20]; | |
104 | for(Int_t iddl=0; iddl<kTotDDL;iddl++){ | |
105 | for(Int_t imod=0; imod<kModPerDDL;imod++){ | |
106 | for(Int_t isid=0;isid<kSides;isid++){ | |
107 | Int_t index=kSides*(kModPerDDL*iddl+imod)+isid; | |
108 | injan[index]=new AliITSOnlineSDDInjectors(iddl,imod,isid); | |
109 | sprintf(hisnam,"h%02dc%02ds%d",iddl,imod,isid); | |
110 | histo[index]=new TH2F(hisnam,"",256,-0.5,255.5,256,-0.5,255.5); | |
111 | nWrittenEv[index]=0; | |
112 | isFilled[index]=0; | |
113 | } | |
114 | } | |
115 | } | |
116 | ||
117 | /* report progress */ | |
118 | daqDA_progressReport(10); | |
119 | ||
120 | Int_t iev=0; | |
121 | Int_t ievInj=0; | |
122 | Int_t ievUsed=0; | |
123 | Int_t nEvToBeSkipped=5; | |
124 | ||
125 | /* read the data files */ | |
126 | int n; | |
127 | for (n=1;n<argc;n++) { | |
128 | ||
129 | status=monitorSetDataSource( argv[n] ); | |
130 | if (status!=0) { | |
131 | printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status)); | |
132 | return -1; | |
133 | } | |
134 | ||
135 | /* report progress */ | |
136 | /* in this example, indexed on the number of files */ | |
137 | // Progress report inside the event loop as well? | |
138 | daqDA_progressReport(10+80*n/argc); | |
139 | ||
140 | /* read the file */ | |
141 | for(;;) { | |
142 | struct eventHeaderStruct *event; | |
143 | eventTypeType eventT; | |
144 | ||
145 | /* get next event */ | |
146 | status=monitorGetEventDynamic((void **)&event); | |
147 | if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */ | |
148 | if (status!=0) { | |
149 | printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status)); | |
150 | return -1; | |
151 | } | |
152 | /* retry if got no event */ | |
153 | if (event==NULL) { | |
154 | break; | |
155 | } | |
156 | iev++; | |
157 | ||
158 | ||
159 | /* use event - here, just write event id to result file */ | |
160 | eventT=event->eventType; | |
161 | switch (event->eventType){ | |
162 | ||
163 | /* START OF RUN */ | |
164 | case START_OF_RUN: | |
165 | break; | |
166 | ||
167 | /* END OF RUN */ | |
168 | case END_OF_RUN: | |
169 | break; | |
170 | ||
171 | // case PHYSICS_EVENT: // comment this line for test raw data | |
172 | // break; // comment this line for test raw data | |
173 | ||
174 | case CALIBRATION_EVENT: | |
175 | break; // uncomment this line for test raw data | |
176 | case PHYSICS_EVENT: // uncomment this line for test raw data | |
177 | printf(" Event number = %i ",iev); | |
178 | ievInj++; | |
179 | if(ievInj<=nEvToBeSkipped){ | |
180 | printf(" -> SKIP\n"); | |
181 | break; | |
182 | } | |
183 | printf(" -> Analyze\n"); | |
184 | ievUsed++; | |
185 | ||
186 | AliRawReader *rawReader = new AliRawReaderDate((void*)event); | |
187 | timeSt=rawReader->GetTimestamp(); | |
188 | rawReader->Reset(); | |
189 | UChar_t cdhAttr=AliITSRawStreamSDD::ReadBlockAttributes(rawReader); | |
190 | amSamplFreq=AliITSRawStreamSDD::ReadAMSamplFreqFromCDH(cdhAttr); | |
191 | AliITSRawStream* s=AliITSRawStreamSDD::CreateRawStreamSDD(rawReader,cdhAttr); | |
192 | if(!writtenoutput){ | |
193 | printf("Use %s raw stream, sampling frequency %d MHz\n",s->ClassName(),amSamplFreq); | |
194 | writtenoutput=kTRUE; | |
195 | } | |
196 | for(Int_t iddl=0; iddl<kTotDDL;iddl++){ | |
197 | for(Int_t imod=0; imod<kModPerDDL;imod++){ | |
198 | for(Int_t isid=0;isid<kSides;isid++){ | |
199 | Int_t index=kSides*(kModPerDDL*iddl+imod)+isid; | |
200 | histo[index]->Reset(); | |
201 | } | |
202 | } | |
203 | } | |
204 | ||
205 | while(s->Next()){ | |
206 | Int_t iDDL=rawReader->GetDDLID(); | |
207 | Int_t iCarlos=s->GetCarlosId(); | |
208 | if(s->IsCompletedModule()) continue; | |
209 | if(s->IsCompletedDDL()) continue; | |
210 | if(iDDL>=0 && iDDL<kTotDDL){ | |
211 | Int_t index=kSides*(kModPerDDL*iDDL+iCarlos)+s->GetChannel(); | |
212 | histo[index]->Fill(s->GetCoord2(),s->GetCoord1(),s->GetSignal()); | |
213 | isFilled[index]=1; | |
214 | } | |
215 | } | |
216 | delete s; | |
217 | delete rawReader; | |
218 | ||
219 | for(Int_t iddl=0; iddl<kTotDDL;iddl++){ | |
220 | for(Int_t imod=0; imod<kModPerDDL;imod++){ | |
221 | for(Int_t isid=0;isid<kSides;isid++){ | |
222 | Int_t index=kSides*(kModPerDDL*iddl+imod)+isid; | |
223 | if(isFilled[index]){ | |
224 | if(amSamplFreq==20) injan[index]->Set20MHzConfig(); | |
225 | else injan[index]->Set40MHzConfig(); | |
226 | injan[index]->AddEvent(histo[index]); | |
227 | // injan[index]->WriteToASCII(iev,timeSt,nWrittenEv[index]); | |
228 | nWrittenEv[index]++; | |
229 | } | |
230 | } | |
231 | } | |
232 | } | |
233 | /* free resources */ | |
234 | free(event); | |
235 | } | |
236 | if(ievUsed>=maxNEvents) break; | |
237 | } | |
238 | } | |
239 | ||
240 | /* write report */ | |
241 | TDatime time; | |
242 | TObjString timeinfo(Form("%02d%02d%02d%02d%02d%02d",time.GetYear()-2000,time.GetMonth(),time.GetDay(),time.GetHour(),time.GetMinute(),time.GetSecond())); | |
243 | printf("Run #%s, received %d calibration events, time %s\n",getenv("DATE_RUN_NUMBER"),ievUsed,timeinfo.GetString().Data()); | |
244 | ||
245 | /* report progress */ | |
246 | daqDA_progressReport(90); | |
247 | TObjArray* dspHistos=new TObjArray(); | |
248 | ||
249 | ||
250 | Char_t filnam[100],command[120]; | |
251 | for(Int_t iddl=0; iddl<kTotDDL;iddl++){ | |
252 | for(Int_t imod=0; imod<kModPerDDL;imod++){ | |
253 | for(Int_t isid=0;isid<kSides;isid++){ | |
254 | Int_t index=kSides*(kModPerDDL*iddl+imod)+isid; | |
255 | if(nWrittenEv[index]>0){ | |
256 | injan[index]->FitMeanDriftSpeedVsAnode(); | |
257 | injan[index]->WriteToASCII(0,timeSt,0); | |
258 | injan[index]->WriteInjectorStatusToASCII(); | |
259 | dspHistos->AddLast(injan[index]->GetMeanDriftSpeedVsPadHisto()); | |
260 | sprintf(filnam,"SDDinj_ddl%02dc%02d_sid%d.data",iddl,imod,isid); | |
261 | sprintf(command,"tar -rf SDDinj_LDC.tar %s",filnam); | |
262 | gSystem->Exec(command); | |
263 | } | |
264 | } | |
265 | } | |
266 | } | |
267 | status=daqDA_FES_storeFile("./SDDinj_LDC.tar","SDD_Injec"); | |
268 | ||
269 | ||
270 | #ifdef ALI_AMORE | |
271 | amore::da::AmoreDA amoreDA(amore::da::AmoreDA::kSender); | |
272 | Int_t statusamore =0; | |
273 | statusamore += amoreDA.Send("TimeInfoInject",&timeinfo); | |
274 | statusamore += amoreDA.Send("DriftSpeed",dspHistos); | |
275 | if ( statusamore ) | |
276 | printf("Warning: Failed to write Arrays in the AMORE database\n"); | |
277 | else | |
278 | printf("amoreDA.Send() OK\n"); | |
279 | #else | |
280 | printf("Warning: SDDINJ DA not compiled with AMORE support\n"); | |
281 | #endif | |
282 | ||
283 | TFile *fh=new TFile("SDDinjectHistos.root","RECREATE"); | |
284 | dspHistos->Write(); | |
285 | fh->Close(); | |
286 | ||
287 | ||
288 | /* report progress */ | |
289 | daqDA_progressReport(100); | |
290 | ||
291 | ||
292 | ||
293 | return status; | |
294 | } |