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