]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/ITSSDDBASda.cxx
Added check for wrong records in the config file (R. Shaoyan)
[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 <TSystem.h>
44 #include <TH1F.h>
45 #include <TH2F.h>
46 #include <TROOT.h>
47 #include <TPluginManager.h>
48
49 // AliRoot includes
50 #include "AliRawReaderDate.h"
51 #include "AliITSOnlineSDDBase.h"
52 #include "AliITSOnlineSDDCMN.h"
53 #include "AliITSRawStreamSDD.h"
54 /* Main routine
55       Arguments: list of DATE raw data files
56 */
57 int main(int argc, char **argv) {
58
59   int status = 0;
60
61   // line added to solve IO problems
62   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
63                                         "*",
64                                         "TStreamerInfo",
65                                         "RIO",
66                                         "TStreamerInfo()");
67
68   /* log start of process */
69   printf("ITS SDD BASELINE+NOISE algorithm program started\n");  
70
71
72   /* check that we got some arguments = list of files */
73   if (argc<2) {
74     printf("Wrong number of arguments\n");
75     return -1;
76   }
77
78
79   Int_t maxNEvents=18; // maximum number of events to be analyzed
80   const Int_t kTotDDL=24;
81   const Int_t kModPerDDL=12;
82   const Int_t kSides=2;
83   Int_t adcSamplFreq=40;
84   if(gSystem->Getenv("DAQ_DETDB_LOCAL")!=NULL){
85     const char* dir=gSystem->Getenv("DAQ_DETDB_LOCAL");    
86     TString filnam=Form("%s/fee.conf",dir); 
87     FILE* feefil=fopen(filnam.Data(),"r"); 
88     fscanf(feefil,"%d \n",&adcSamplFreq);
89     fclose(feefil);
90   }
91   
92   AliITSOnlineSDDBase **base=new AliITSOnlineSDDBase*[kTotDDL*kModPerDDL*kSides];
93   AliITSOnlineSDDCMN **corr=new AliITSOnlineSDDCMN*[kTotDDL*kModPerDDL*kSides];
94   TH2F **histo=new TH2F*[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         base[index]=new AliITSOnlineSDDBase(iddl,imod,isid);
102         if(adcSamplFreq==20) base[index]->SetLastGoodTB(126);
103         else base[index]->SetLastGoodTB(254);
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       }
107     }
108   }
109   
110   /* report progress */
111   daqDA_progressReport(8);
112   Int_t iev;
113   for(Int_t iStep=0;iStep<2;iStep++){
114     printf("Start Analysis Step %d\n",iStep);
115     iev=0;
116     if(iStep==1){
117       for(Int_t iddl=0; iddl<kTotDDL;iddl++){
118         for(Int_t imod=0; imod<kModPerDDL;imod++){
119           for(Int_t isid=0;isid<kSides;isid++){
120             Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
121             corr[index]=new AliITSOnlineSDDCMN(iddl,imod,isid);
122             if(adcSamplFreq==20) corr[index]->SetLastGoodTB(126);
123             else corr[index]->SetLastGoodTB(254);
124           }
125         }
126       }
127     }
128
129     /* read the data files */
130     int n;
131     for (n=1;n<argc;n++) {
132    
133       status=monitorSetDataSource( argv[n] );
134       if (status!=0) {
135         printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
136         return -1;
137       }
138
139       /* report progress */
140       /* in this example, indexed on the number of files */
141       // Progress report inside the event loop as well?
142       daqDA_progressReport(10+40*iStep*n/argc);
143
144       /* read the file */
145       for(;;) {
146         struct eventHeaderStruct *event;
147         eventTypeType eventT;
148
149         /* get next event */
150         status=monitorGetEventDynamic((void **)&event);
151         if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
152         if (status!=0) {
153           printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
154           return -1;
155         }
156
157         /* retry if got no event */
158         if (event==NULL) {
159           break;
160         }
161
162         iev++; 
163         if(iev>maxNEvents) break;
164
165         /* use event - here, just write event id to result file */
166         eventT=event->eventType;
167         switch (event->eventType){
168       
169           /* START OF RUN */
170         case START_OF_RUN:
171           break;
172       
173           /* END OF RUN */
174         case END_OF_RUN:
175           break;
176
177           //      case PHYSICS_EVENT:  // comment this line for test raw data
178           //    break;               // comment this line for test raw data
179
180
181         case CALIBRATION_EVENT:
182           break;  // uncomment this line for test raw data
183         case PHYSICS_EVENT: // uncomment this line for test raw data
184           printf(" event number = %i \n",iev);
185           AliRawReader *rawReader = new AliRawReaderDate((void*)event);
186
187
188           Int_t evtyp=0;
189           while(rawReader->ReadHeader()){
190             const UInt_t *subev = rawReader->GetSubEventAttributes();
191             if(subev[0]==0 && subev[1]==0 && subev[2]==0) evtyp=1; 
192           }
193
194           rawReader->Reset();
195           for(Int_t iddl=0; iddl<kTotDDL;iddl++){
196             for(Int_t imod=0; imod<kModPerDDL;imod++){
197               for(Int_t isid=0;isid<kSides;isid++){
198                 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
199                 histo[index]->Reset();
200               }
201             }
202           }
203           AliITSRawStreamSDD s(rawReader);
204
205           while(s.Next()){
206             Int_t iDDL=rawReader->GetDDLID();
207             Int_t iCarlos=s.GetCarlosId();
208             if(s.IsCompletedModule()) continue;
209             if(s.IsCompletedDDL()) continue;
210             if(iDDL>=0 && iDDL<kTotDDL){ 
211               Int_t index=kSides*(kModPerDDL*iDDL+iCarlos)+s.GetChannel(); 
212               histo[index]->Fill(s.GetCoord2(),s.GetCoord1(),s.GetSignal());
213             }
214           }
215           delete rawReader;
216           for(Int_t iddl=0; iddl<kTotDDL;iddl++){
217             for(Int_t imod=0; imod<kModPerDDL;imod++){
218               for(Int_t isid=0;isid<kSides;isid++){
219                 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
220                 if(iStep==0) base[index]->AddEvent(histo[index]);
221                 if(iStep==1) corr[index]->AddEvent(histo[index]);
222               }
223             }
224           }
225
226           /* free resources */
227           free(event);
228         }
229       }
230     }
231
232     if(iStep==0){
233       for(Int_t iddl=0; iddl<kTotDDL;iddl++){
234         for(Int_t imod=0; imod<kModPerDDL;imod++){
235           for(Int_t isid=0;isid<kSides;isid++){
236             Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
237             base[index]->ValidateAnodes();
238             base[index]->WriteToASCII();
239           }
240         }
241       }  
242     }
243   }
244
245   /* write report */
246   printf("Run #%s, received %d calibration events\n",getenv("DATE_RUN_NUMBER"),iev);
247
248   /* report progress */
249   daqDA_progressReport(90);
250
251
252
253   TFile *fh=new TFile("SDDbaseHistos.root","RECREATE");
254   for(Int_t iddl=0; iddl<kTotDDL;iddl++){
255     for(Int_t imod=0; imod<kModPerDDL;imod++){
256       for(Int_t isid=0;isid<kSides;isid++){
257         Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
258         corr[index]->ValidateAnodes();
259         corr[index]->WriteToASCII();
260         corr[index]->WriteToROOT(fh);
261       }
262     }
263   }
264   fh->Close();
265
266
267   /* report progress */
268   daqDA_progressReport(100);
269
270
271
272   return status;
273 }