]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/ITSSDDBASda.cxx
6993725f28b8613f0ff77d6a5f322b7a03873784
[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 maxNEvents=10; // maximum number of events to be analyzed
59   const Int_t nSDDmodules=12;  // temp for test raw data
60   AliITSOnlineSDDBase **base=new AliITSOnlineSDDBase*[2*nSDDmodules];
61   AliITSOnlineSDDCMN **corr=new AliITSOnlineSDDCMN*[2*nSDDmodules];
62   TH2F **histo=new TH2F*[2*nSDDmodules];
63
64   Char_t hisnam[20];
65   for(Int_t imod=0; imod<nSDDmodules;imod++){
66     for(Int_t isid=0;isid<2;isid++){
67       Int_t index=2*imod+isid;
68       base[index]=new AliITSOnlineSDDBase(imod,isid);
69       corr[index]=new AliITSOnlineSDDCMN(imod,isid);
70       sprintf(hisnam,"his%03ds%d",imod,isid);
71       histo[index]=new TH2F(hisnam,"",256,-0.5,255.5,256,-0.5,255.5);
72     }
73   }
74   
75   /* report progress */
76   daqDA_progressReport(8);
77   Int_t iev;
78   for(Int_t iStep=0;iStep<2;iStep++){
79     printf("Start Analysis Step %d\n",iStep);
80     iev=0;
81     if(iStep==1){
82       for(Int_t imod=0; imod<nSDDmodules;imod++){
83         for(Int_t isid=0;isid<2;isid++){
84           Int_t index=2*imod+isid;
85           corr[index]->Reset();
86         }
87       }
88     }
89
90     /* read the data files */
91     int n;
92     for (n=1;n<argc;n++) {
93    
94       status=monitorSetDataSource( argv[n] );
95       if (status!=0) {
96         printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
97         return -1;
98       }
99
100       /* report progress */
101       /* in this example, indexed on the number of files */
102       // Progress report inside the event loop as well?
103       daqDA_progressReport(10+40*iStep*n/argc);
104
105       /* read the file */
106       for(;;) {
107         struct eventHeaderStruct *event;
108         eventTypeType eventT;
109
110         /* get next event */
111         status=monitorGetEventDynamic((void **)&event);
112         if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
113         if (status!=0) {
114           printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
115           return -1;
116         }
117
118         /* retry if got no event */
119         if (event==NULL) {
120           break;
121         }
122
123         iev++; 
124         if(iev>maxNEvents) break;
125
126         /* use event - here, just write event id to result file */
127         eventT=event->eventType;
128         switch (event->eventType){
129       
130           /* START OF RUN */
131         case START_OF_RUN:
132           break;
133       
134           /* END OF RUN */
135         case END_OF_RUN:
136           break;
137
138           //      case PHYSICS_EVENT:  // comment this line for test raw data
139           //    break;               // comment this line for test raw data
140
141
142         case CALIBRATION_EVENT:
143           break;  // uncomment this line for test raw data
144         case PHYSICS_EVENT: // uncomment this line for test raw data
145           printf(" event number = %i \n",iev);
146           AliRawReader *rawReader = new AliRawReaderDate((void*)event);
147           rawReader->RequireHeader(kFALSE);
148           rawReader->SelectEquipment(17,101,101);  // temp for test raw data
149
150
151           Int_t evtyp=0;
152           while(rawReader->ReadHeader()){
153             const UInt_t *subev = rawReader->GetSubEventAttributes();
154             if(subev[0]==0 && subev[1]==0 && subev[2]==0) evtyp=1; 
155           }
156
157           rawReader->Reset();
158           for(Int_t imod=0; imod<nSDDmodules;imod++){
159             for(Int_t isid=0; isid<2;isid++){
160               Int_t index=2*imod+isid;
161               histo[index]->Reset();
162             }
163           }
164           AliITSRawStreamSDD s(rawReader);
165           rawReader->SelectEquipment(17,101,101);  // temp for test raw data
166
167           while(s.Next()){
168             Int_t iddl=rawReader->GetDDLID();
169             iddl=0; // temporary for test raw data
170             Int_t isddmod=s.GetModuleNumber(iddl,s.GetCarlosId()); 
171             isddmod-=240;  // to have SDD modules from 0 to 259
172             isddmod=s.GetCarlosId(); // temporary for test raw data
173             if(isddmod<nSDDmodules){ 
174               Int_t index=2*isddmod+s.GetChannel(); 
175               histo[index]->Fill(s.GetCoord2(),s.GetCoord1(),s.GetSignal());
176             }
177           }
178           delete rawReader;
179           for(Int_t imod=0; imod<nSDDmodules;imod++){
180             for(Int_t isid=0; isid<2;isid++){
181               Int_t index=2*imod+isid;
182               if(iStep==0) base[index]->AddEvent(histo[index]);    
183               if(iStep==1) corr[index]->AddEvent(histo[index]);    
184             }
185           }
186
187           /* free resources */
188           free(event);
189         }
190       }
191     
192     }
193     
194     for(Int_t imod=0; imod<nSDDmodules;imod++){
195       for(Int_t isid=0; isid<2;isid++){
196         Int_t index=2*imod+isid;
197         if(iStep==0){
198           base[index]->ValidateAnodes();
199           base[index]->WriteToASCII();
200         }
201         if(iStep==1){
202           corr[index]->ValidateAnodes();
203           corr[index]->WriteToASCII();
204         }
205       }
206     }  
207   }
208
209   /* write report */
210   printf("Run #%s, received %d calibration events\n",getenv("DATE_RUN_NUMBER"),iev);
211
212   /* report progress */
213   daqDA_progressReport(90);
214
215
216
217   TFile *fh=new TFile("SDDbaseHistos.root","RECREATE");
218   for(Int_t imod=0; imod<nSDDmodules;imod++){
219     for(Int_t isid=0; isid<2;isid++){
220       Int_t index=2*imod+isid;
221       corr[index]->WriteToROOT(fh);
222     }
223   }
224   fh->Close();
225
226
227   /* report progress */
228   daqDA_progressReport(100);
229
230
231
232   return status;
233 }