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