]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/ITSSDDBASda.cxx
### files: AliTPCTempMap.h (.cxx)
[u/mrichter/AliRoot.git] / ITS / ITSSDDBASda.cxx
1 /*
2 - Contact: - prino@to.infn.it
3 - Link: - http://www.to.infn.it/~prino/alice/RawData/run11161.date
4 - Run Type: - PEDESTAL_RUN
5 - DA Type: - LDC
6 - Number of events needed: 100
7 - Input Files: - 
8 - Output Files: - SDDbase_step1_ddl*c*_sid*.data SDDbase_step2_ddl*c*_sid*.data
9 - Trigger types used: 
10 */
11
12
13 //////////////////////////////////////////////////////////////////////////////
14 // Detector Algorithm for analysis of SDD baseline runs.                    //
15 //                                                                          //
16 // Produces ASCII and ROOT output files with:                               //
17 // 1. anode quality bit                                                     //
18 // 1. Baseline values                                                       //
19 // 2. Raw noise                                                             //
20 // 3. Common mode coefficients                                              //
21 // 4. Noise corrected for common mode                                       //
22 // Files are stored locally.                                                //
23 // The next DAQ-DA step on Test Pulse run writes files to FXS               //
24 //                                                                          //
25 // Author: F. Prino (prino@to.infn.it)                                      //
26 //                                                                          //
27 //////////////////////////////////////////////////////////////////////////////
28
29
30
31 extern "C" {
32 #include "daqDA.h"
33 }
34
35 #include "event.h"
36 #include "monitor.h"
37
38 #include <stdio.h>
39 #include <stdlib.h>
40
41 // ROOT includes
42 #include <TFile.h>
43 #include <TH1F.h>
44 #include <TH2F.h>
45 #include <TROOT.h>
46 #include <TPluginManager.h>
47
48 // AliRoot includes
49 #include "AliRawReaderDate.h"
50 #include "AliITSOnlineSDDBase.h"
51 #include "AliITSOnlineSDDCMN.h"
52 #include "AliITSRawStreamSDD.h"
53 /* Main routine
54       Arguments: list of DATE raw data files
55 */
56 int main(int argc, char **argv) {
57
58   int status = 0;
59
60   // line added to solve IO problems
61   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
62                                         "*",
63                                         "TStreamerInfo",
64                                         "RIO",
65                                         "TStreamerInfo()");
66
67   /* log start of process */
68   printf("ITS SDD BASELINE+NOISE algorithm program started\n");  
69
70
71   /* check that we got some arguments = list of files */
72   if (argc<2) {
73     printf("Wrong number of arguments\n");
74     return -1;
75   }
76
77
78   Int_t maxNEvents=18; // maximum number of events to be analyzed
79   const Int_t kTotDDL=24;
80   const Int_t kModPerDDL=12;
81   const Int_t kSides=2;
82   AliITSOnlineSDDBase **base=new AliITSOnlineSDDBase*[kTotDDL*kModPerDDL*kSides];
83   AliITSOnlineSDDCMN **corr=new AliITSOnlineSDDCMN*[kTotDDL*kModPerDDL*kSides];
84   TH2F **histo=new TH2F*[kTotDDL*kModPerDDL*kSides];
85
86   Char_t hisnam[20];
87   for(Int_t iddl=0; iddl<kTotDDL;iddl++){
88     for(Int_t imod=0; imod<kModPerDDL;imod++){
89       for(Int_t isid=0;isid<kSides;isid++){
90         Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
91         base[index]=new AliITSOnlineSDDBase(iddl,imod,isid);
92         sprintf(hisnam,"h%02dc%02ds%d",iddl,imod,isid);
93         histo[index]=new TH2F(hisnam,"",256,-0.5,255.5,256,-0.5,255.5);
94       }
95     }
96   }
97   
98   /* report progress */
99   daqDA_progressReport(8);
100   Int_t iev;
101   for(Int_t iStep=0;iStep<2;iStep++){
102     printf("Start Analysis Step %d\n",iStep);
103     iev=0;
104     if(iStep==1){
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             corr[index]=new AliITSOnlineSDDCMN(iddl,imod,isid);
110           }
111         }
112       }
113     }
114
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+40*iStep*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
167         case CALIBRATION_EVENT:
168           break;  // uncomment this line for test raw data
169         case PHYSICS_EVENT: // uncomment this line for test raw data
170           printf(" event number = %i \n",iev);
171           AliRawReader *rawReader = new AliRawReaderDate((void*)event);
172
173
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(iDDL>=0 && iDDL<kTotDDL && s.IsCompletedModule()==kFALSE){ 
195               Int_t index=kSides*(kModPerDDL*iDDL+iCarlos)+s.GetChannel(); 
196               histo[index]->Fill(s.GetCoord2(),s.GetCoord1(),s.GetSignal());
197             }
198           }
199           delete rawReader;
200           for(Int_t iddl=0; iddl<kTotDDL;iddl++){
201             for(Int_t imod=0; imod<kModPerDDL;imod++){
202               for(Int_t isid=0;isid<kSides;isid++){
203                 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
204                 if(iStep==0) base[index]->AddEvent(histo[index]);
205                 if(iStep==1) corr[index]->AddEvent(histo[index]);
206               }
207             }
208           }
209
210           /* free resources */
211           free(event);
212         }
213       }
214     }
215
216     if(iStep==0){
217       for(Int_t iddl=0; iddl<kTotDDL;iddl++){
218         for(Int_t imod=0; imod<kModPerDDL;imod++){
219           for(Int_t isid=0;isid<kSides;isid++){
220             Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
221             base[index]->ValidateAnodes();
222             base[index]->WriteToASCII();
223           }
224         }
225       }  
226     }
227   }
228
229   /* write report */
230   printf("Run #%s, received %d calibration events\n",getenv("DATE_RUN_NUMBER"),iev);
231
232   /* report progress */
233   daqDA_progressReport(90);
234
235
236
237   TFile *fh=new TFile("SDDbaseHistos.root","RECREATE");
238   for(Int_t iddl=0; iddl<kTotDDL;iddl++){
239     for(Int_t imod=0; imod<kModPerDDL;imod++){
240       for(Int_t isid=0;isid<kSides;isid++){
241         Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
242         corr[index]->ValidateAnodes();
243         corr[index]->WriteToASCII();
244         corr[index]->WriteToROOT(fh);
245       }
246     }
247   }
248   fh->Close();
249
250
251   /* report progress */
252   daqDA_progressReport(100);
253
254
255
256   return status;
257 }