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