]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/ITSSDDGAINda.cxx
Updated comment for compliance with the parsing algorithm
[u/mrichter/AliRoot.git] / ITS / ITSSDDGAINda.cxx
1 /*
2 - Contact: - prino@to.infn.it
3 - Link: - http://www.to.infn.it/~prino/alice/RawData/run11173.date
4 - Run Type: - PULSER_RUN
5 - DA Type: - LDC
6 - Number of events needed: 100
7 - Input Files: - SDDbase_step1_ddl*c*_sid*.data
8 - Output Files: - SDDbase_ddl*c*_sid*.data
9 - Trigger types used: 
10 */
11
12 //////////////////////////////////////////////////////////////////////////////
13 // Detector Algorithm for analysis of SDD test pulse runs.                  //
14 //                                                                          //
15 // Produces ASCII and ROOT output files with:                               //
16 // 1. anode quality bit                                                     //
17 // 1. Baseline values                                                       //
18 // 2. Raw noise                                                             //
19 // 3. Common mode coefficients                                              //
20 // 4. Noise corrected for common mode                                       //
21 // 5. Gain                                                                  //
22 // Files are written to FXS                                                 //
23 //                                                                          //
24 // Author: F. Prino (prino@to.infn.it)                                      //
25 //                                                                          //
26 //////////////////////////////////////////////////////////////////////////////
27
28
29
30 extern "C" {
31 #include "daqDA.h"
32 }
33
34 #include "event.h"
35 #include "monitor.h"
36
37 #include <stdio.h>
38 #include <stdlib.h>
39
40 // ROOT includes
41 #include <TFile.h>
42 #include <TH1F.h>
43 #include <TH2F.h>
44 #include <TSystem.h>
45 #include <TROOT.h>
46 #include <TPluginManager.h>
47
48 // AliRoot includes
49 #include "AliRawReaderDate.h"
50 #include "AliITSOnlineSDDTP.h"
51 #include "AliITSRawStreamSDD.h"
52 /* Main routine
53       Arguments: list of DATE raw data files
54 */
55 int main(int argc, char **argv) {
56
57   int status = 0;
58
59   // line added to solve IO problems
60   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
61                                         "*",
62                                         "TStreamerInfo",
63                                         "RIO",
64                                         "TStreamerInfo()");
65
66   /* log start of process */
67   printf("ITS SDD TEST-PULSE algorithm program started\n");  
68
69
70   /* check that we got some arguments = list of files */
71   if (argc<2) {
72     printf("Wrong number of arguments\n");
73     return -1;
74   }
75
76
77   Int_t eqOffset = 256;
78   Int_t DDLrange = 24;
79   Int_t maxNEvents=15; // 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
84   AliITSOnlineSDDTP **tpan=new AliITSOnlineSDDTP*[kTotDDL*kModPerDDL*kSides];
85   TH2F **histo=new TH2F*[kTotDDL*kModPerDDL*kSides];
86   Bool_t isFilled[kTotDDL*kModPerDDL*kSides];
87   Char_t hisnam[20];
88   for(Int_t iddl=0; iddl<kTotDDL;iddl++){
89     for(Int_t imod=0; imod<kModPerDDL;imod++){
90       for(Int_t isid=0;isid<kSides;isid++){
91         Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
92         tpan[index]=new AliITSOnlineSDDTP(iddl,imod,isid,100.);
93         sprintf(hisnam,"h%02dc%02ds%d",iddl,imod,isid);
94         histo[index]=new TH2F(hisnam,"",256,-0.5,255.5,256,-0.5,255.5);
95         isFilled[index]=0;
96       }
97     }
98   }
99   
100   /* report progress */
101   daqDA_progressReport(10);
102   Int_t iev=0,iAnalyzedEv=0;
103   /* read the data files */
104   int n;
105   for (n=1;n<argc;n++) {
106    
107     status=monitorSetDataSource( argv[n] );
108     if (status!=0) {
109       printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
110       return -1;
111     }
112
113     /* report progress */
114     /* in this example, indexed on the number of files */
115     // Progress report inside the event loop as well?
116     daqDA_progressReport(10+80*n/argc);
117
118     /* read the file */
119     for(;;) {
120       struct eventHeaderStruct *event;
121       eventTypeType eventT;
122
123       /* get next event */
124       status=monitorGetEventDynamic((void **)&event);
125       if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
126       if (status!=0) {
127         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
128         return -1;
129       }
130
131       /* retry if got no event */
132       if (event==NULL) {
133         break;
134       }
135
136       if(iAnalyzedEv>=maxNEvents) break;
137       iev++; 
138       
139       /* use event - here, just write event id to result file */
140       eventT=event->eventType;
141       switch (event->eventType){
142         
143         /* START OF RUN */
144       case START_OF_RUN:
145         break;
146         
147         /* END OF RUN */
148       case END_OF_RUN:
149         break;
150         
151         //      case PHYSICS_EVENT:  // comment this line for test raw data
152         //      break;               // comment this line for test raw data
153
154       case CALIBRATION_EVENT:
155         break;  // uncomment this line for test raw data
156       case PHYSICS_EVENT: // uncomment this line for test raw data
157         printf(" event number = %i \n",iev);
158         AliRawReader *rawReader = new AliRawReaderDate((void*)event);
159         rawReader->SelectEquipment(17,eqOffset+1,eqOffset+DDLrange);
160
161         Int_t evtyp=0;
162         while(rawReader->ReadHeader()){
163           const UInt_t *subev = rawReader->GetSubEventAttributes();
164           if(subev[0]==0 && subev[1]==0 && subev[2]==0) evtyp=1; 
165         }
166
167         rawReader->Reset();
168         for(Int_t iddl=0; iddl<kTotDDL;iddl++){
169           for(Int_t imod=0; imod<kModPerDDL;imod++){
170             for(Int_t isid=0;isid<kSides;isid++){
171               Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
172               histo[index]->Reset();
173             }
174           }
175         }
176         AliITSRawStreamSDD s(rawReader);
177         
178         while(s.Next()){
179           Int_t iDDL=rawReader->GetDDLID();
180           Int_t iCarlos=s.GetCarlosId();
181           if(iDDL>=0 && iDDL<kTotDDL && s.IsCompletedModule()==kFALSE){ 
182             Int_t index=kSides*(kModPerDDL*iDDL+iCarlos)+s.GetChannel(); 
183             histo[index]->Fill(s.GetCoord2(),s.GetCoord1(),s.GetSignal());
184             isFilled[index]=1;
185           }
186         }
187         delete rawReader;
188         for(Int_t iddl=0; iddl<kTotDDL;iddl++){
189           for(Int_t imod=0; imod<kModPerDDL;imod++){
190             for(Int_t isid=0;isid<kSides;isid++){
191               Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
192               if(isFilled[index]) tpan[index]->AddEvent(histo[index]);    
193             }
194           }
195         }
196         
197         /* free resources */
198         iAnalyzedEv++;
199         free(event);
200       }
201     }
202     
203   }
204     
205   TFile *fh=new TFile("SDDgainHistos.root","RECREATE");
206   Char_t filnam[100],command[120];
207   for(Int_t iddl=0; iddl<kTotDDL;iddl++){
208     for(Int_t imod=0; imod<kModPerDDL;imod++){
209       for(Int_t isid=0;isid<kSides;isid++){
210         Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
211         if(isFilled[index]){
212           tpan[index]->ValidateAnodes();
213           tpan[index]->WriteToASCII();
214           tpan[index]->WriteToROOT(fh);
215           sprintf(filnam,"SDDbase_ddl%02dc%02d_sid%d.data",iddl,imod,isid);
216           sprintf(command,"tar -rf SDDbase_LDC.tar %s",filnam);
217           gSystem->Exec(command);
218         }
219       }
220     }  
221   }
222   fh->Close();
223
224   /* write report */
225   printf("Run #%s, received %d calibration events\n",getenv("DATE_RUN_NUMBER"),iAnalyzedEv);
226
227   /* report progress */
228   daqDA_progressReport(90);
229
230
231
232   fh->Close();
233
234
235
236   status=daqDA_FES_storeFile("./SDDbase_LDC.tar","SDD_Calib");
237
238   /* report progress */
239   daqDA_progressReport(100);
240
241
242
243   return status;
244 }