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