]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/ITSSDDINJda.cxx
increased number of channels used in the TOFFEElight data structure, from 157248...
[u/mrichter/AliRoot.git] / ITS / ITSSDDINJda.cxx
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_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
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 */
54 int main(int argc, char **argv) {
55
56   int status = 0;
57
58
59   // line added to solve IO problems
60   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
61                                         "*",
62                                         "TStreamerInfo",
63                                         "RIO",
64                                         "TStreamerInfo()");
65
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
78   const Int_t kTotDDL=24;
79   const Int_t kModPerDDL=12;
80   const Int_t kSides=2;
81   Int_t adcSamplFreq=40;
82   Bool_t readfeeconf=kFALSE;
83   gSystem->Exec("rm -f SDDinj_ddl*.data");
84   if(gSystem->Getenv("DAQ_DETDB_LOCAL")!=NULL){
85     const char* dir=gSystem->Getenv("DAQ_DETDB_LOCAL");    
86     TString filnam=Form("%s/fee.conf",dir);    
87     FILE* feefil=fopen(filnam.Data(),"r"); 
88     if(feefil){
89       fscanf(feefil,"%d \n",&adcSamplFreq);
90       fclose(feefil);
91       readfeeconf=kTRUE;
92       printf("ADC sampling frequency = %d MHz\n",adcSamplFreq);
93     }
94   }
95   if(!readfeeconf) printf("File fee.conf not found, sampling frequency set to 40 MHz\n");
96
97
98
99   AliITSOnlineSDDInjectors **injan=new AliITSOnlineSDDInjectors*[kTotDDL*kModPerDDL*kSides];
100   TH2F **histo=new TH2F*[kTotDDL*kModPerDDL*kSides];
101   Int_t nWrittenEv[kTotDDL*kModPerDDL*kSides];
102   Bool_t isFilled[kTotDDL*kModPerDDL*kSides];
103
104   Char_t hisnam[20];
105   for(Int_t iddl=0; iddl<kTotDDL;iddl++){
106     for(Int_t imod=0; imod<kModPerDDL;imod++){
107       for(Int_t isid=0;isid<kSides;isid++){
108         Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
109         injan[index]=new AliITSOnlineSDDInjectors(iddl,imod,isid);
110         if(adcSamplFreq==20) injan[index]->Set20MHzConfig();
111         else injan[index]->Set40MHzConfig();
112         sprintf(hisnam,"h%02dc%02ds%d",iddl,imod,isid);
113         histo[index]=new TH2F(hisnam,"",256,-0.5,255.5,256,-0.5,255.5);
114         nWrittenEv[index]=0;
115         isFilled[index]=0;
116       }
117     }
118   }
119   
120   /* report progress */
121   daqDA_progressReport(10);
122   Int_t iev=0;
123   Int_t ievInj=0;
124   /* read the data files */
125   int n;
126   for (n=1;n<argc;n++) {
127    
128     status=monitorSetDataSource( argv[n] );
129     if (status!=0) {
130       printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
131       return -1;
132     }
133
134     /* report progress */
135     /* in this example, indexed on the number of files */
136     // Progress report inside the event loop as well?
137     daqDA_progressReport(10+80*n/argc);
138
139     /* read the file */
140     for(;;) {
141       struct eventHeaderStruct *event;
142       eventTypeType eventT;
143
144       /* get next event */
145       status=monitorGetEventDynamic((void **)&event);
146       if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
147       if (status!=0) {
148         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
149         return -1;
150       }
151       /* retry if got no event */
152       if (event==NULL) {
153         break;
154       }
155       iev++;
156
157       
158       /* use event - here, just write event id to result file */
159       eventT=event->eventType;
160       switch (event->eventType){
161         
162         /* START OF RUN */
163       case START_OF_RUN:
164         break;
165         
166         /* END OF RUN */
167       case END_OF_RUN:
168         break;
169         
170         //      case PHYSICS_EVENT:  // comment this line for test raw data
171         //      break;               // comment this line for test raw data
172
173       case CALIBRATION_EVENT:
174         break;  // uncomment this line for test raw data
175       case PHYSICS_EVENT: // uncomment this line for test raw data
176         printf(" event number = %i \n",iev);
177         ievInj++; 
178         AliRawReader *rawReader = new AliRawReaderDate((void*)event);
179
180         UInt_t timeSt=rawReader->GetTimestamp();
181
182         Int_t evtyp=0;
183         while(rawReader->ReadHeader()){
184           const UInt_t *subev = rawReader->GetSubEventAttributes();
185           if(subev[0]==0 && subev[1]==0 && subev[2]==0) evtyp=1; 
186         }
187
188         rawReader->Reset();
189         for(Int_t iddl=0; iddl<kTotDDL;iddl++){
190           for(Int_t imod=0; imod<kModPerDDL;imod++){
191             for(Int_t isid=0;isid<kSides;isid++){
192               Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
193               histo[index]->Reset();
194             }
195           }
196         }
197         AliITSRawStreamSDD s(rawReader);
198         
199         while(s.Next()){
200           Int_t iDDL=rawReader->GetDDLID();
201           Int_t iCarlos=s.GetCarlosId();
202           if(s.IsCompletedModule()) continue;
203           if(s.IsCompletedDDL()) continue;
204           if(iDDL>=0 && iDDL<kTotDDL){ 
205             Int_t index=kSides*(kModPerDDL*iDDL+iCarlos)+s.GetChannel(); 
206             histo[index]->Fill(s.GetCoord2(),s.GetCoord1(),s.GetSignal());
207             isFilled[index]=1;
208           }
209         }
210         delete rawReader;
211         for(Int_t iddl=0; iddl<kTotDDL;iddl++){
212           for(Int_t imod=0; imod<kModPerDDL;imod++){
213             for(Int_t isid=0;isid<kSides;isid++){
214               Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
215               if(isFilled[index]){
216                 injan[index]->Reset();
217                 injan[index]->AnalyzeEvent(histo[index]);    
218                 injan[index]->WriteToASCII(iev,timeSt,nWrittenEv[index]);
219                 nWrittenEv[index]++;
220               }
221             }
222           }
223         }
224         /* free resources */
225         free(event);
226       }
227       if(ievInj>=maxNEvents) break;
228     }
229   }
230   /* write report */
231   printf("Run #%s, received %d injector events\n",getenv("DATE_RUN_NUMBER"),ievInj);
232
233   gSystem->Exec("rm -f  SDDinj_LDC.tar");
234   Char_t filnam[100],command[120];
235   for(Int_t iddl=0; iddl<kTotDDL;iddl++){
236     for(Int_t imod=0; imod<kModPerDDL;imod++){
237       for(Int_t isid=0;isid<kSides;isid++){
238         Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
239         if(nWrittenEv[index]>0){
240           sprintf(filnam,"SDDinj_ddl%02dc%02d_sid%d.data",iddl,imod,isid);
241           sprintf(command,"tar -rf SDDinj_LDC.tar %s",filnam);
242           gSystem->Exec(command);
243         }
244       }  
245     }
246   }
247
248
249   /* report progress */
250   daqDA_progressReport(90);
251
252
253   status=daqDA_FES_storeFile("./SDDinj_LDC.tar","SDD_Injec");
254
255   /* report progress */
256   daqDA_progressReport(100);
257
258
259
260   return status;
261 }