]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/ITSSDDGAINda.cxx
Adding new selection criteria for filtering of fake
[u/mrichter/AliRoot.git] / ITS / ITSSDDGAINda.cxx
CommitLineData
00252d8f 1/*
2- Contact: - prino@to.infn.it
1ccdd571 3- Link: - http://www.to.infn.it/~prino/alice/RawData/run11173.date
4- Run Type: - PULSER_RUN
00252d8f 5- DA Type: - LDC
6- Number of events needed: 100
979b5a5f 7- Input Files: - SDDbase_step1_ddl*c*_sid*.data
8- Output Files: - SDDbase_ddl*c*_sid*.data
00252d8f 9- Trigger types used:
10*/
11
3083967f 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
30extern "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>
00252d8f 45#include <TROOT.h>
46#include <TPluginManager.h>
3083967f 47
48// AliRoot includes
49#include "AliRawReaderDate.h"
50#include "AliITSOnlineSDDTP.h"
51#include "AliITSRawStreamSDD.h"
a8cd2ab2 52#include "AliITSRawStreamSDDCompressed.h"
3083967f 53/* Main routine
54 Arguments: list of DATE raw data files
55*/
56int main(int argc, char **argv) {
57
58 int status = 0;
59
00252d8f 60 // line added to solve IO problems
61 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
62 "*",
63 "TStreamerInfo",
64 "RIO",
65 "TStreamerInfo()");
3083967f 66
67 /* log start of process */
979b5a5f 68 printf("ITS SDD TEST-PULSE algorithm program started\n");
3083967f 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
1ccdd571 78 Int_t maxNEvents=15; // maximum number of events to be analyzed
979b5a5f 79 const Int_t kTotDDL=24;
80 const Int_t kModPerDDL=12;
81 const Int_t kSides=2;
42e868a7 82 Int_t adcSamplFreq=40;
5065707f 83 Bool_t readfeeconf=kFALSE;
a8cd2ab2 84 Int_t dataformat=1;
5065707f 85 gSystem->Exec("rm -f SDDbase_LDC.tar");
86 if(gSystem->Getenv("DAQ_DETDB_LOCAL")!=NULL){
42e868a7 87 const char* dir=gSystem->Getenv("DAQ_DETDB_LOCAL");
88 TString filnam=Form("%s/fee.conf",dir);
89 FILE* feefil=fopen(filnam.Data(),"r");
5065707f 90 if(feefil){
91 fscanf(feefil,"%d \n",&adcSamplFreq);
a8cd2ab2 92 fscanf(feefil,"%d \n",&dataformat);
93 fclose(feefil);
5065707f 94 TString shcomm=Form("tar -rf SDDbase_LDC.tar -C %s fee.conf",dir);
95 gSystem->Exec(shcomm.Data());
5065707f 96 readfeeconf=kTRUE;
a8cd2ab2 97 printf("ADC sampling frequency = %d MHz dataformat code =%d \n",adcSamplFreq,dataformat);
5065707f 98 }
42e868a7 99 }
a8cd2ab2 100 if(!readfeeconf) printf("File fee.conf not found, sampling frequency set to 40 MHz and data format to 1\n");
979b5a5f 101
102 AliITSOnlineSDDTP **tpan=new AliITSOnlineSDDTP*[kTotDDL*kModPerDDL*kSides];
103 TH2F **histo=new TH2F*[kTotDDL*kModPerDDL*kSides];
104 Bool_t isFilled[kTotDDL*kModPerDDL*kSides];
3083967f 105 Char_t hisnam[20];
979b5a5f 106 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
107 for(Int_t imod=0; imod<kModPerDDL;imod++){
108 for(Int_t isid=0;isid<kSides;isid++){
109 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
110 tpan[index]=new AliITSOnlineSDDTP(iddl,imod,isid,100.);
42e868a7 111 if(adcSamplFreq==20) tpan[index]->SetLastGoodTB(126);
112 else tpan[index]->SetLastGoodTB(254);
979b5a5f 113 sprintf(hisnam,"h%02dc%02ds%d",iddl,imod,isid);
114 histo[index]=new TH2F(hisnam,"",256,-0.5,255.5,256,-0.5,255.5);
115 isFilled[index]=0;
116 }
3083967f 117 }
118 }
119
120 /* report progress */
121 daqDA_progressReport(10);
1ccdd571 122 Int_t iev=0,iAnalyzedEv=0;
3083967f 123 /* read the data files */
124 int n;
125 for (n=1;n<argc;n++) {
126
127 status=monitorSetDataSource( argv[n] );
128 if (status!=0) {
129 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
130 return -1;
131 }
132
133 /* report progress */
134 /* in this example, indexed on the number of files */
135 // Progress report inside the event loop as well?
136 daqDA_progressReport(10+80*n/argc);
137
138 /* read the file */
139 for(;;) {
140 struct eventHeaderStruct *event;
141 eventTypeType eventT;
142
143 /* get next event */
144 status=monitorGetEventDynamic((void **)&event);
145 if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
146 if (status!=0) {
147 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
148 return -1;
149 }
150
151 /* retry if got no event */
152 if (event==NULL) {
153 break;
154 }
155
1ccdd571 156 if(iAnalyzedEv>=maxNEvents) break;
3083967f 157 iev++;
3083967f 158
159 /* use event - here, just write event id to result file */
160 eventT=event->eventType;
161 switch (event->eventType){
162
163 /* START OF RUN */
164 case START_OF_RUN:
165 break;
166
167 /* END OF RUN */
168 case END_OF_RUN:
169 break;
170
171 // case PHYSICS_EVENT: // comment this line for test raw data
172 // break; // comment this line for test raw data
173
174 case CALIBRATION_EVENT:
175 break; // uncomment this line for test raw data
176 case PHYSICS_EVENT: // uncomment this line for test raw data
177 printf(" event number = %i \n",iev);
178 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
3083967f 179
180 Int_t evtyp=0;
181 while(rawReader->ReadHeader()){
182 const UInt_t *subev = rawReader->GetSubEventAttributes();
183 if(subev[0]==0 && subev[1]==0 && subev[2]==0) evtyp=1;
184 }
185
186 rawReader->Reset();
979b5a5f 187 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
188 for(Int_t imod=0; imod<kModPerDDL;imod++){
189 for(Int_t isid=0;isid<kSides;isid++){
190 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
191 histo[index]->Reset();
192 }
3083967f 193 }
194 }
a8cd2ab2 195 AliITSRawStream* s;
196 if(dataformat==0){
197 s=new AliITSRawStreamSDD(rawReader);
198 }else{
199 s=new AliITSRawStreamSDDCompressed(rawReader);
200 if(dataformat==1) s->SetADCEncoded(kTRUE);
201 }
3083967f 202
a8cd2ab2 203 while(s->Next()){
979b5a5f 204 Int_t iDDL=rawReader->GetDDLID();
a8cd2ab2 205 Int_t iCarlos=s->GetCarlosId();
206 if(s->IsCompletedModule()) continue;
207 if(s->IsCompletedDDL()) continue;
de075dae 208 if(iDDL>=0 && iDDL<kTotDDL){
a8cd2ab2 209 Int_t index=kSides*(kModPerDDL*iDDL+iCarlos)+s->GetChannel();
210 histo[index]->Fill(s->GetCoord2(),s->GetCoord1(),s->GetSignal());
1ccdd571 211 isFilled[index]=1;
3083967f 212 }
213 }
214 delete rawReader;
979b5a5f 215 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
216 for(Int_t imod=0; imod<kModPerDDL;imod++){
217 for(Int_t isid=0;isid<kSides;isid++){
218 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
219 if(isFilled[index]) tpan[index]->AddEvent(histo[index]);
220 }
3083967f 221 }
222 }
223
224 /* free resources */
1ccdd571 225 iAnalyzedEv++;
3083967f 226 free(event);
227 }
228 }
229
230 }
231
232 TFile *fh=new TFile("SDDgainHistos.root","RECREATE");
233 Char_t filnam[100],command[120];
979b5a5f 234 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
235 for(Int_t imod=0; imod<kModPerDDL;imod++){
236 for(Int_t isid=0;isid<kSides;isid++){
237 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
238 if(isFilled[index]){
239 tpan[index]->ValidateAnodes();
240 tpan[index]->WriteToASCII();
241 tpan[index]->WriteToROOT(fh);
242 sprintf(filnam,"SDDbase_ddl%02dc%02d_sid%d.data",iddl,imod,isid);
243 sprintf(command,"tar -rf SDDbase_LDC.tar %s",filnam);
244 gSystem->Exec(command);
245 }
1ccdd571 246 }
3083967f 247 }
248 }
249 fh->Close();
250
251 /* write report */
1ccdd571 252 printf("Run #%s, received %d calibration events\n",getenv("DATE_RUN_NUMBER"),iAnalyzedEv);
3083967f 253
254 /* report progress */
255 daqDA_progressReport(90);
256
257
258
259 fh->Close();
260
261
262
1ccdd571 263 status=daqDA_FES_storeFile("./SDDbase_LDC.tar","SDD_Calib");
3083967f 264
265 /* report progress */
266 daqDA_progressReport(100);
267
268
269
270 return status;
271}