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