]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/ITSSDDGAINda.cxx
DA for SDD pedestal runs updated with info for DAQ team (link to test data) and usage...
[u/mrichter/AliRoot.git] / ITS / ITSSDDGAINda.cxx
CommitLineData
00252d8f 1/*
2- Contact: - prino@to.infn.it
3- Link: -
4- Run Type: -
5- DA Type: - LDC
6- Number of events needed: 100
7- Input Files: - SDDbase_step1_mod*_sid*.data
8- Output Files: - SDDbase_mod*_sid*.data
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"
52/* Main routine
53 Arguments: list of DATE raw data files
54*/
55int main(int argc, char **argv) {
56
57 int status = 0;
58
00252d8f 59 // line added to solve IO problems
60 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
61 "*",
62 "TStreamerInfo",
63 "RIO",
64 "TStreamerInfo()");
3083967f 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
767934e8 77 Int_t eqOffset = 256;
78 Int_t DDLrange = 24;
3083967f 79 Int_t maxNEvents=10; // maximum number of events to be analyzed
80 const Int_t nSDDmodules=12; // temp for test raw data
81 AliITSOnlineSDDTP **tpan=new AliITSOnlineSDDTP*[2*nSDDmodules];
82 TH2F **histo=new TH2F*[2*nSDDmodules];
83
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,411.);
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 }
92 }
93
94 /* report progress */
95 daqDA_progressReport(10);
96 Int_t iev=0;
97 /* read the data files */
98 int n;
99 for (n=1;n<argc;n++) {
100
101 status=monitorSetDataSource( argv[n] );
102 if (status!=0) {
103 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
104 return -1;
105 }
106
107 /* report progress */
108 /* in this example, indexed on the number of files */
109 // Progress report inside the event loop as well?
110 daqDA_progressReport(10+80*n/argc);
111
112 /* read the file */
113 for(;;) {
114 struct eventHeaderStruct *event;
115 eventTypeType eventT;
116
117 /* get next event */
118 status=monitorGetEventDynamic((void **)&event);
119 if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
120 if (status!=0) {
121 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
122 return -1;
123 }
124
125 /* retry if got no event */
126 if (event==NULL) {
127 break;
128 }
129
130 iev++;
131 if(iev>maxNEvents) break;
132
133 /* use event - here, just write event id to result file */
134 eventT=event->eventType;
135 switch (event->eventType){
136
137 /* START OF RUN */
138 case START_OF_RUN:
139 break;
140
141 /* END OF RUN */
142 case END_OF_RUN:
143 break;
144
145 // case PHYSICS_EVENT: // comment this line for test raw data
146 // break; // comment this line for test raw data
147
148 case CALIBRATION_EVENT:
149 break; // uncomment this line for test raw data
150 case PHYSICS_EVENT: // uncomment this line for test raw data
151 printf(" event number = %i \n",iev);
152 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
767934e8 153 rawReader->SelectEquipment(17,eqOffset+1,eqOffset+DDLrange);
3083967f 154
155 Int_t evtyp=0;
156 while(rawReader->ReadHeader()){
157 const UInt_t *subev = rawReader->GetSubEventAttributes();
158 if(subev[0]==0 && subev[1]==0 && subev[2]==0) evtyp=1;
159 }
160
161 rawReader->Reset();
162 for(Int_t imod=0; imod<nSDDmodules;imod++){
163 for(Int_t isid=0; isid<2;isid++){
164 Int_t index=2*imod+isid;
165 histo[index]->Reset();
166 }
167 }
168 AliITSRawStreamSDD s(rawReader);
3083967f 169
170 while(s.Next()){
171 Int_t iddl=rawReader->GetDDLID();
172 iddl=0; // temporary for test raw data
173 Int_t isddmod=s.GetModuleNumber(iddl,s.GetCarlosId());
174 isddmod-=240; // to have SDD modules from 0 to 259
175 isddmod=s.GetCarlosId(); // temporary for test raw data
5dfa68c5 176 if(isddmod<nSDDmodules&&s.IsCompletedModule()==kFALSE){
3083967f 177 Int_t index=2*isddmod+s.GetChannel();
178 histo[index]->Fill(s.GetCoord2(),s.GetCoord1(),s.GetSignal());
179 }
180 }
181 delete rawReader;
182 for(Int_t imod=0; imod<nSDDmodules;imod++){
183 for(Int_t isid=0; isid<2;isid++){
184 Int_t index=2*imod+isid;
185 tpan[index]->AddEvent(histo[index]);
186 }
187 }
188
189 /* free resources */
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 tpan[index]->ValidateAnodes();
202 tpan[index]->WriteToASCII();
203 tpan[index]->WriteToROOT(fh);
204 sprintf(filnam,"SDDbase_mod%03d_sid%d.data",imod,isid);
205 sprintf(command,"tar -rf SDDbase_LDC1.tar %s",filnam);
206 gSystem->Exec(command);
207 }
208 }
209 fh->Close();
210
211 /* write report */
212 printf("Run #%s, received %d calibration events\n",getenv("DATE_RUN_NUMBER"),iev);
213
214 /* report progress */
215 daqDA_progressReport(90);
216
217
218
219 fh->Close();
220
221
222
223 status=daqDA_FES_storeFile("./SDDbase_LDC1.tar","SDD_Calib");
224
225 /* report progress */
226 daqDA_progressReport(100);
227
228
229
230 return status;
231}