]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/ITSSDDBASda.cxx
dq/dt and dq/dl reference plots
[u/mrichter/AliRoot.git] / ITS / ITSSDDBASda.cxx
CommitLineData
00252d8f 1/*
2- Contact: - prino@to.infn.it
f2c71038 3- Link: - http://www.to.infn.it/~prino/alice/RawData/run11161.date
4- Run Type: - PEDESTAL_RUN
00252d8f 5- DA Type: - LDC
6- Number of events needed: 100
7- Input Files: -
979b5a5f 8- Output Files: - SDDbase_step1_ddl*c*_sid*.data SDDbase_step2_ddl*c*_sid*.data
00252d8f 9- Trigger types used:
10*/
11
12
e44f571c 13//////////////////////////////////////////////////////////////////////////////
14// Detector Algorithm for analysis of SDD baseline runs. //
15// //
16// Produces ASCII and ROOT output files with: //
17// 1. anode quality bit //
18// 1. Baseline values //
19// 2. Raw noise //
20// 3. Common mode coefficients //
21// 4. Noise corrected for common mode //
22// Files are stored locally. //
23// The next DAQ-DA step on Test Pulse run writes files to FXS //
24// //
25// Author: F. Prino (prino@to.infn.it) //
26// //
27//////////////////////////////////////////////////////////////////////////////
28
29
30
31extern "C" {
32#include "daqDA.h"
33}
34
35#include "event.h"
36#include "monitor.h"
37
38#include <stdio.h>
39#include <stdlib.h>
40
41// ROOT includes
42#include <TFile.h>
42e868a7 43#include <TSystem.h>
e44f571c 44#include <TH1F.h>
45#include <TH2F.h>
00252d8f 46#include <TROOT.h>
47#include <TPluginManager.h>
e44f571c 48
49// AliRoot includes
50#include "AliRawReaderDate.h"
51#include "AliITSOnlineSDDBase.h"
52#include "AliITSOnlineSDDCMN.h"
53#include "AliITSRawStreamSDD.h"
a8cd2ab2 54#include "AliITSRawStreamSDDCompressed.h"
e44f571c 55/* Main routine
56 Arguments: list of DATE raw data files
57*/
58int main(int argc, char **argv) {
59
60 int status = 0;
61
00252d8f 62 // line added to solve IO problems
63 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
64 "*",
65 "TStreamerInfo",
66 "RIO",
67 "TStreamerInfo()");
e44f571c 68
69 /* log start of process */
979b5a5f 70 printf("ITS SDD BASELINE+NOISE algorithm program started\n");
e44f571c 71
72
73 /* check that we got some arguments = list of files */
74 if (argc<2) {
75 printf("Wrong number of arguments\n");
76 return -1;
77 }
78
79
5065707f 80 Int_t maxNEvents=10; // maximum number of events to be analyzed
979b5a5f 81 const Int_t kTotDDL=24;
82 const Int_t kModPerDDL=12;
83 const Int_t kSides=2;
42e868a7 84 Int_t adcSamplFreq=40;
5065707f 85 Bool_t readfeeconf=kFALSE;
a8cd2ab2 86 Int_t dataformat=1;
5065707f 87 gSystem->Exec("rm -f SDDbase_*.data");
42e868a7 88 if(gSystem->Getenv("DAQ_DETDB_LOCAL")!=NULL){
89 const char* dir=gSystem->Getenv("DAQ_DETDB_LOCAL");
90 TString filnam=Form("%s/fee.conf",dir);
91 FILE* feefil=fopen(filnam.Data(),"r");
5065707f 92 if(feefil){
93 fscanf(feefil,"%d \n",&adcSamplFreq);
a8cd2ab2 94 fscanf(feefil,"%d \n",&dataformat);
5065707f 95 fclose(feefil);
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");
42e868a7 101
979b5a5f 102 AliITSOnlineSDDBase **base=new AliITSOnlineSDDBase*[kTotDDL*kModPerDDL*kSides];
103 AliITSOnlineSDDCMN **corr=new AliITSOnlineSDDCMN*[kTotDDL*kModPerDDL*kSides];
104 TH2F **histo=new TH2F*[kTotDDL*kModPerDDL*kSides];
e44f571c 105
106 Char_t hisnam[20];
979b5a5f 107 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
108 for(Int_t imod=0; imod<kModPerDDL;imod++){
109 for(Int_t isid=0;isid<kSides;isid++){
110 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
111 base[index]=new AliITSOnlineSDDBase(iddl,imod,isid);
42e868a7 112 if(adcSamplFreq==20) base[index]->SetLastGoodTB(126);
113 else base[index]->SetLastGoodTB(254);
979b5a5f 114 sprintf(hisnam,"h%02dc%02ds%d",iddl,imod,isid);
115 histo[index]=new TH2F(hisnam,"",256,-0.5,255.5,256,-0.5,255.5);
116 }
e44f571c 117 }
118 }
119
120 /* report progress */
3083967f 121 daqDA_progressReport(8);
e44f571c 122 Int_t iev;
5065707f 123 Int_t ievPed;
e44f571c 124 for(Int_t iStep=0;iStep<2;iStep++){
e44f571c 125 printf("Start Analysis Step %d\n",iStep);
126 iev=0;
5065707f 127 ievPed=0;
e44f571c 128 if(iStep==1){
979b5a5f 129 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
130 for(Int_t imod=0; imod<kModPerDDL;imod++){
131 for(Int_t isid=0;isid<kSides;isid++){
132 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
133 corr[index]=new AliITSOnlineSDDCMN(iddl,imod,isid);
42e868a7 134 if(adcSamplFreq==20) corr[index]->SetLastGoodTB(126);
135 else corr[index]->SetLastGoodTB(254);
979b5a5f 136 }
e44f571c 137 }
138 }
139 }
140
141 /* read the data files */
142 int n;
143 for (n=1;n<argc;n++) {
144
145 status=monitorSetDataSource( argv[n] );
146 if (status!=0) {
147 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
148 return -1;
149 }
150
151 /* report progress */
152 /* in this example, indexed on the number of files */
153 // Progress report inside the event loop as well?
3083967f 154 daqDA_progressReport(10+40*iStep*n/argc);
e44f571c 155
156 /* read the file */
157 for(;;) {
158 struct eventHeaderStruct *event;
159 eventTypeType eventT;
160
161 /* get next event */
162 status=monitorGetEventDynamic((void **)&event);
163 if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
164 if (status!=0) {
165 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
166 return -1;
167 }
168
169 /* retry if got no event */
170 if (event==NULL) {
171 break;
172 }
173
174 iev++;
e44f571c 175
176 /* use event - here, just write event id to result file */
177 eventT=event->eventType;
178 switch (event->eventType){
179
180 /* START OF RUN */
181 case START_OF_RUN:
182 break;
183
184 /* END OF RUN */
185 case END_OF_RUN:
186 break;
187
3083967f 188 // case PHYSICS_EVENT: // comment this line for test raw data
189 // break; // comment this line for test raw data
190
e44f571c 191
192 case CALIBRATION_EVENT:
3083967f 193 break; // uncomment this line for test raw data
194 case PHYSICS_EVENT: // uncomment this line for test raw data
e44f571c 195 printf(" event number = %i \n",iev);
5065707f 196 ievPed++;
e44f571c 197 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
e44f571c 198
e44f571c 199
200 Int_t evtyp=0;
201 while(rawReader->ReadHeader()){
202 const UInt_t *subev = rawReader->GetSubEventAttributes();
203 if(subev[0]==0 && subev[1]==0 && subev[2]==0) evtyp=1;
204 }
205
206 rawReader->Reset();
979b5a5f 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 histo[index]->Reset();
212 }
e44f571c 213 }
214 }
a8cd2ab2 215 AliITSRawStream* s;
216 if(dataformat==0){
217 s=new AliITSRawStreamSDD(rawReader);
218 }else{
219 s=new AliITSRawStreamSDDCompressed(rawReader);
220 if(dataformat==1) s->SetADCEncoded(kTRUE);
221 }
e44f571c 222
a8cd2ab2 223 while(s->Next()){
979b5a5f 224 Int_t iDDL=rawReader->GetDDLID();
a8cd2ab2 225 Int_t iCarlos=s->GetCarlosId();
226 if(s->IsCompletedModule()) continue;
227 if(s->IsCompletedDDL()) continue;
de075dae 228 if(iDDL>=0 && iDDL<kTotDDL){
a8cd2ab2 229 Int_t index=kSides*(kModPerDDL*iDDL+iCarlos)+s->GetChannel();
230 histo[index]->Fill(s->GetCoord2(),s->GetCoord1(),s->GetSignal());
e44f571c 231 }
232 }
a8cd2ab2 233 delete s;
e44f571c 234 delete rawReader;
979b5a5f 235 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
236 for(Int_t imod=0; imod<kModPerDDL;imod++){
237 for(Int_t isid=0;isid<kSides;isid++){
238 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
239 if(iStep==0) base[index]->AddEvent(histo[index]);
240 if(iStep==1) corr[index]->AddEvent(histo[index]);
241 }
e44f571c 242 }
243 }
244
245 /* free resources */
246 free(event);
247 }
5065707f 248 if(ievPed>=maxNEvents) break;
e44f571c 249 }
e44f571c 250 }
979b5a5f 251
252 if(iStep==0){
253 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
254 for(Int_t imod=0; imod<kModPerDDL;imod++){
255 for(Int_t isid=0;isid<kSides;isid++){
256 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
257 base[index]->ValidateAnodes();
258 base[index]->WriteToASCII();
259 }
e44f571c 260 }
5065707f 261 }
979b5a5f 262 }
e44f571c 263 }
264
265 /* write report */
5065707f 266 printf("Run #%s, received %d calibration events\n",getenv("DATE_RUN_NUMBER"),ievPed);
e44f571c 267
268 /* report progress */
269 daqDA_progressReport(90);
270
271
272
3083967f 273 TFile *fh=new TFile("SDDbaseHistos.root","RECREATE");
979b5a5f 274 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
275 for(Int_t imod=0; imod<kModPerDDL;imod++){
276 for(Int_t isid=0;isid<kSides;isid++){
277 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
278 corr[index]->ValidateAnodes();
279 corr[index]->WriteToASCII();
280 corr[index]->WriteToROOT(fh);
281 }
e44f571c 282 }
283 }
284 fh->Close();
285
286
e44f571c 287 /* report progress */
288 daqDA_progressReport(100);
289
290
291
292 return status;
293}