]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/ITSSDDGAINda.cxx
SDD detector algorithms udated due to a change in the equipment id
[u/mrichter/AliRoot.git] / ITS / ITSSDDGAINda.cxx
1 //////////////////////////////////////////////////////////////////////////////
2 // Detector Algorithm for analysis of SDD test pulse 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 // 5. Gain                                                                  //
11 // Files are written 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 #include <TSystem.h>
34
35 // AliRoot includes
36 #include "AliRawReaderDate.h"
37 #include "AliITSOnlineSDDTP.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   Int_t eqOffset = 256;
59   Int_t DDLrange = 24;
60   Int_t maxNEvents=10; // maximum number of events to be analyzed
61   const Int_t nSDDmodules=12;  // temp for test raw data
62   AliITSOnlineSDDTP **tpan=new AliITSOnlineSDDTP*[2*nSDDmodules];
63   TH2F **histo=new TH2F*[2*nSDDmodules];
64
65   Char_t hisnam[20];
66   for(Int_t imod=0; imod<nSDDmodules;imod++){
67     for(Int_t isid=0;isid<2;isid++){
68       Int_t index=2*imod+isid;
69       tpan[index]=new AliITSOnlineSDDTP(imod,isid,411.);
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=0;
78   /* read the data files */
79   int n;
80   for (n=1;n<argc;n++) {
81    
82     status=monitorSetDataSource( argv[n] );
83     if (status!=0) {
84       printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
85       return -1;
86     }
87
88     /* report progress */
89     /* in this example, indexed on the number of files */
90     // Progress report inside the event loop as well?
91     daqDA_progressReport(10+80*n/argc);
92
93     /* read the file */
94     for(;;) {
95       struct eventHeaderStruct *event;
96       eventTypeType eventT;
97
98       /* get next event */
99       status=monitorGetEventDynamic((void **)&event);
100       if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
101       if (status!=0) {
102         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
103         return -1;
104       }
105
106       /* retry if got no event */
107       if (event==NULL) {
108         break;
109       }
110
111       iev++; 
112       if(iev>maxNEvents) break;
113       
114       /* use event - here, just write event id to result file */
115       eventT=event->eventType;
116       switch (event->eventType){
117         
118         /* START OF RUN */
119       case START_OF_RUN:
120         break;
121         
122         /* END OF RUN */
123       case END_OF_RUN:
124         break;
125         
126         //      case PHYSICS_EVENT:  // comment this line for test raw data
127         //      break;               // comment this line for test raw data
128
129       case CALIBRATION_EVENT:
130         break;  // uncomment this line for test raw data
131       case PHYSICS_EVENT: // uncomment this line for test raw data
132         printf(" event number = %i \n",iev);
133         AliRawReader *rawReader = new AliRawReaderDate((void*)event);
134         rawReader->RequireHeader(kFALSE);
135         rawReader->SelectEquipment(17,eqOffset+1,eqOffset+DDLrange);
136
137         Int_t evtyp=0;
138         while(rawReader->ReadHeader()){
139           const UInt_t *subev = rawReader->GetSubEventAttributes();
140           if(subev[0]==0 && subev[1]==0 && subev[2]==0) evtyp=1; 
141         }
142
143         rawReader->Reset();
144         for(Int_t imod=0; imod<nSDDmodules;imod++){
145           for(Int_t isid=0; isid<2;isid++){
146             Int_t index=2*imod+isid;
147             histo[index]->Reset();
148           }
149         }
150         AliITSRawStreamSDD s(rawReader);
151         
152         while(s.Next()){
153           Int_t iddl=rawReader->GetDDLID();
154           iddl=0; // temporary for test raw data
155           Int_t isddmod=s.GetModuleNumber(iddl,s.GetCarlosId()); 
156           isddmod-=240;  // to have SDD modules from 0 to 259
157           isddmod=s.GetCarlosId(); // temporary for test raw data
158           if(isddmod<nSDDmodules){ 
159             Int_t index=2*isddmod+s.GetChannel(); 
160             histo[index]->Fill(s.GetCoord2(),s.GetCoord1(),s.GetSignal());
161           }
162         }
163         delete rawReader;
164         for(Int_t imod=0; imod<nSDDmodules;imod++){
165           for(Int_t isid=0; isid<2;isid++){
166             Int_t index=2*imod+isid;
167             tpan[index]->AddEvent(histo[index]);    
168           }
169         }
170         
171         /* free resources */
172         free(event);
173       }
174     }
175     
176   }
177     
178   TFile *fh=new TFile("SDDgainHistos.root","RECREATE");
179   Char_t filnam[100],command[120];
180   for(Int_t imod=0; imod<nSDDmodules;imod++){
181     for(Int_t isid=0; isid<2;isid++){
182       Int_t index=2*imod+isid;
183       tpan[index]->ValidateAnodes();
184       tpan[index]->WriteToASCII();
185       tpan[index]->WriteToROOT(fh);
186       sprintf(filnam,"SDDbase_mod%03d_sid%d.data",imod,isid);
187       sprintf(command,"tar -rf SDDbase_LDC1.tar %s",filnam);
188       gSystem->Exec(command);
189     }  
190   }
191   fh->Close();
192
193   /* write report */
194   printf("Run #%s, received %d calibration events\n",getenv("DATE_RUN_NUMBER"),iev);
195
196   /* report progress */
197   daqDA_progressReport(90);
198
199
200
201   fh->Close();
202
203
204
205   status=daqDA_FES_storeFile("./SDDbase_LDC1.tar","SDD_Calib");
206
207   /* report progress */
208   daqDA_progressReport(100);
209
210
211
212   return status;
213 }