]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/ITSSDDBASda.cxx
Fix problem in simulated raw data + code cleanup (F.Prino)
[u/mrichter/AliRoot.git] / ITS / ITSSDDBASda.cxx
1 //////////////////////////////////////////////////////////////////////////////
2 // Detector Algorithm for analysis of SDD baseline runs.                    //
3 //                                                                          //
4 // Produces ASCII and ROOT output files with:                               //
5 // 1. anode quality bit                                                     //
6 // 1. Baseline values                                                       //
7 // 2. Raw noise                                                             //
8 // 3. Common mode coefficients                                              //
9 // 4. Noise corrected for common mode                                       //
10 // Files are stored locally.                                                //
11 // The next DAQ-DA step on Test Pulse run writes files to FXS               //
12 //                                                                          //
13 // Author: F. Prino (prino@to.infn.it)                                      //
14 //                                                                          //
15 //////////////////////////////////////////////////////////////////////////////
16
17
18
19 extern "C" {
20 #include "daqDA.h"
21 }
22
23 #include "event.h"
24 #include "monitor.h"
25
26 #include <stdio.h>
27 #include <stdlib.h>
28
29 // ROOT includes
30 #include <TFile.h>
31 #include <TH1F.h>
32 #include <TH2F.h>
33
34 // AliRoot includes
35 #include "AliRawReaderDate.h"
36 #include "AliITSOnlineSDDBase.h"
37 #include "AliITSOnlineSDDCMN.h"
38 #include "AliITSRawStreamSDD.h"
39 /* Main routine
40       Arguments: list of DATE raw data files
41 */
42 int main(int argc, char **argv) {
43
44   int status = 0;
45
46
47   /* log start of process */
48   printf("ITS SDD BASE algorithm program started\n");  
49
50
51   /* check that we got some arguments = list of files */
52   if (argc<2) {
53     printf("Wrong number of arguments\n");
54     return -1;
55   }
56
57
58   Int_t eqOffset = 256;
59   Int_t DDLrange = 24;
60   Int_t maxNEvents=10; // maximum number of events to be analyzed
61   const Int_t nSDDmodules=12;  // temp for test raw data
62   AliITSOnlineSDDBase **base=new AliITSOnlineSDDBase*[2*nSDDmodules];
63   AliITSOnlineSDDCMN **corr=new AliITSOnlineSDDCMN*[2*nSDDmodules];
64   TH2F **histo=new TH2F*[2*nSDDmodules];
65
66   Char_t hisnam[20];
67   for(Int_t imod=0; imod<nSDDmodules;imod++){
68     for(Int_t isid=0;isid<2;isid++){
69       Int_t index=2*imod+isid;
70       base[index]=new AliITSOnlineSDDBase(imod,isid);
71       corr[index]=new AliITSOnlineSDDCMN(imod,isid);
72       sprintf(hisnam,"his%03ds%d",imod,isid);
73       histo[index]=new TH2F(hisnam,"",256,-0.5,255.5,256,-0.5,255.5);
74     }
75   }
76   
77   /* report progress */
78   daqDA_progressReport(8);
79   Int_t iev;
80   for(Int_t iStep=0;iStep<2;iStep++){
81     printf("Start Analysis Step %d\n",iStep);
82     iev=0;
83     if(iStep==1){
84       for(Int_t imod=0; imod<nSDDmodules;imod++){
85         for(Int_t isid=0;isid<2;isid++){
86           Int_t index=2*imod+isid;
87           corr[index]->Reset();
88         }
89       }
90     }
91
92     /* read the data files */
93     int n;
94     for (n=1;n<argc;n++) {
95    
96       status=monitorSetDataSource( argv[n] );
97       if (status!=0) {
98         printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
99         return -1;
100       }
101
102       /* report progress */
103       /* in this example, indexed on the number of files */
104       // Progress report inside the event loop as well?
105       daqDA_progressReport(10+40*iStep*n/argc);
106
107       /* read the file */
108       for(;;) {
109         struct eventHeaderStruct *event;
110         eventTypeType eventT;
111
112         /* get next event */
113         status=monitorGetEventDynamic((void **)&event);
114         if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
115         if (status!=0) {
116           printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
117           return -1;
118         }
119
120         /* retry if got no event */
121         if (event==NULL) {
122           break;
123         }
124
125         iev++; 
126         if(iev>maxNEvents) break;
127
128         /* use event - here, just write event id to result file */
129         eventT=event->eventType;
130         switch (event->eventType){
131       
132           /* START OF RUN */
133         case START_OF_RUN:
134           break;
135       
136           /* END OF RUN */
137         case END_OF_RUN:
138           break;
139
140           //      case PHYSICS_EVENT:  // comment this line for test raw data
141           //    break;               // comment this line for test raw data
142
143
144         case CALIBRATION_EVENT:
145           break;  // uncomment this line for test raw data
146         case PHYSICS_EVENT: // uncomment this line for test raw data
147           printf(" event number = %i \n",iev);
148           AliRawReader *rawReader = new AliRawReaderDate((void*)event);
149           rawReader->RequireHeader(kFALSE);
150           rawReader->SelectEquipment(17,eqOffset+1,eqOffset+DDLrange);
151
152
153           Int_t evtyp=0;
154           while(rawReader->ReadHeader()){
155             const UInt_t *subev = rawReader->GetSubEventAttributes();
156             if(subev[0]==0 && subev[1]==0 && subev[2]==0) evtyp=1; 
157           }
158
159           rawReader->Reset();
160           for(Int_t imod=0; imod<nSDDmodules;imod++){
161             for(Int_t isid=0; isid<2;isid++){
162               Int_t index=2*imod+isid;
163               histo[index]->Reset();
164             }
165           }
166           AliITSRawStreamSDD s(rawReader);
167
168           while(s.Next()){
169             Int_t iddl=rawReader->GetDDLID();
170             iddl=0; // temporary for test raw data
171             Int_t isddmod=s.GetModuleNumber(iddl,s.GetCarlosId()); 
172             isddmod-=240;  // to have SDD modules from 0 to 259
173             isddmod=s.GetCarlosId(); // temporary for test raw data
174             if(isddmod<nSDDmodules&&s.IsCompletedModule()==kFALSE){ 
175               Int_t index=2*isddmod+s.GetChannel(); 
176               histo[index]->Fill(s.GetCoord2(),s.GetCoord1(),s.GetSignal());
177             }
178           }
179           delete rawReader;
180           for(Int_t imod=0; imod<nSDDmodules;imod++){
181             for(Int_t isid=0; isid<2;isid++){
182               Int_t index=2*imod+isid;
183               if(iStep==0) base[index]->AddEvent(histo[index]);    
184               if(iStep==1) corr[index]->AddEvent(histo[index]);    
185             }
186           }
187
188           /* free resources */
189           free(event);
190         }
191       }
192     
193     }
194     
195     for(Int_t imod=0; imod<nSDDmodules;imod++){
196       for(Int_t isid=0; isid<2;isid++){
197         Int_t index=2*imod+isid;
198         if(iStep==0){
199           base[index]->ValidateAnodes();
200           base[index]->WriteToASCII();
201         }
202         if(iStep==1){
203           corr[index]->ValidateAnodes();
204           corr[index]->WriteToASCII();
205         }
206       }
207     }  
208   }
209
210   /* write report */
211   printf("Run #%s, received %d calibration events\n",getenv("DATE_RUN_NUMBER"),iev);
212
213   /* report progress */
214   daqDA_progressReport(90);
215
216
217
218   TFile *fh=new TFile("SDDbaseHistos.root","RECREATE");
219   for(Int_t imod=0; imod<nSDDmodules;imod++){
220     for(Int_t isid=0; isid<2;isid++){
221       Int_t index=2*imod+isid;
222       corr[index]->WriteToROOT(fh);
223     }
224   }
225   fh->Close();
226
227
228   /* report progress */
229   daqDA_progressReport(100);
230
231
232
233   return status;
234 }