]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/TPCPEDESTALda.cxx
Adding simple example to load default debug streamer
[u/mrichter/AliRoot.git] / TPC / TPCPEDESTALda.cxx
CommitLineData
9c754ce7 1/*
2TPC DA for online calibration
3
4Contact: Haavard.Helstrup@cern.ch
5Link:
1b37fc95 6Run Type: PEDESTAL
9c754ce7 7DA Type: LDC
8Number of events needed: 100
9Input Files:
10Output Files: tpcPedestal.root, to be exported to the DAQ FXS
bd955ed2 11fileId: pedestals
9c754ce7 12Trigger types used: CALIBRATION_EVENT
13
14*/
15
c12208b8 16/*
17
18TPCda_pedestal.cxx - calibration algorithm for TPC pedestal runs
19
2010/06/2007 sylvain.chapeland@cern.ch : first version - clean skeleton based on DAQ DA case1
49efac78 2119/10/2007 christian.lippmann@cern.ch : Possibility to write output to ASCII file
2224/10/2007 christian.lippmann@cern.ch : Including pedestal calibration for time bins
2323/11/2007 christian.lippmann@cern.ch : Fix in order to avoid streamer problems in case of
24 invalid ROOTSTYS. The famous magic line provided by Rene.
2528/11/2007 christian.lippmann@cern.ch : TPC mapping file is read from DaqDetDB
c12208b8 26
27contact: marian.ivanov@cern.ch
28
c12208b8 29This process reads RAW data from the files provided as command line arguments
30and save results in a file (named from RESULT_FILE define - see below).
31
32*/
33
49efac78 34#define RESULT_FILE "tpcPedestal.root"
bd955ed2 35#define FILE_ID "pedestals"
49efac78 36#define MAPPING_FILE "tpcMapping.root"
b401648b 37#define AliDebugLevel() -1
c12208b8 38
39extern "C" {
40#include <daqDA.h>
41}
42#include "event.h"
43#include "monitor.h"
49efac78 44
45#include "stdio.h"
46#include "stdlib.h"
47#include <fstream>
c12208b8 48
49//
50//Root includes
51//
49efac78 52#include "TFile.h"
53#include "TArrayF.h"
54#include "TROOT.h"
55#include "TPluginManager.h"
c12208b8 56
57//
58//AliRoot includes
59//
60#include "AliRawReader.h"
61#include "AliRawReaderDate.h"
87d57858 62#include "AliTPCmapper.h"
c12208b8 63#include "AliTPCRawStream.h"
64#include "AliTPCROC.h"
65#include "AliTPCCalROC.h"
66#include "AliTPCCalPad.h"
67#include "AliMathBase.h"
68#include "TTreeStream.h"
b401648b 69#include "AliLog.h"
70#include "TSystem.h"
c12208b8 71
72//
73// TPC calibration algorithm includes
74//
75#include "AliTPCCalibPedestal.h"
76
bdf99a93 77/*
78 Main routine, TPC pedestal detector algorithm to be run on TPC LDC
79 Arguments: list of DATE raw data files
c12208b8 80*/
bdf99a93 81
c12208b8 82int main(int argc, char **argv) {
bdf99a93 83 //
84 // Main for TPC pedestal detector algorithm
85 //
c12208b8 86
bd955ed2 87 AliLog::SetClassDebugLevel("AliTPCRawStream",-5);
88 AliLog::SetClassDebugLevel("AliRawReaderDate",-5);
89 AliLog::SetClassDebugLevel("AliTPCAltroMapping",-5);
90 AliLog::SetModuleDebugLevel("RAW",-5);
b401648b 91
49efac78 92 Bool_t timeAnalysis = kTRUE;
9c754ce7 93
c12208b8 94 if (argc<2) {
95 printf("Wrong number of arguments\n");
96 return -1;
97 }
98
49efac78 99 /* magic line */
100 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
101 "*",
102 "TStreamerInfo",
103 "RIO",
104 "TStreamerInfo()");
105 int i, status;
106
c12208b8 107 /* log start of process */
108 printf("TPC DA started - %s\n",__FILE__);
109
c12208b8 110 /* declare monitoring program */
111 status=monitorDeclareMp( __FILE__ );
112 if (status!=0) {
113 printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
114 return -1;
115 }
116
b401648b 117 AliTPCmapper *mapping = 0; // The TPC mapping
5d694e7e 118
b401648b 119 if (!mapping){
120 /* copy locally the mapping file from daq detector config db */
121 status = daqDA_DB_getFile(MAPPING_FILE,"./tpcMapping.root");
122 if (status) {
123 printf("Failed to get mapping file (%s) from DAQdetDB, status=%d\n", MAPPING_FILE, status);
124 printf("Continue anyway ... maybe it works?\n"); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
125 return -1; // temporarily uncommented for testing on pcald47 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
126 }
49efac78 127
b401648b 128 /* open the mapping file and retrieve mapping object */
129 TFile *fileMapping = new TFile(MAPPING_FILE, "read");
130 mapping = (AliTPCmapper*) fileMapping->Get("tpcMapping");
49efac78 131 delete fileMapping;
b401648b 132 }
133
134 if (mapping == 0) {
135 printf("Failed to get mapping object from %s. ...\n", MAPPING_FILE);
136 //return -1;
49efac78 137 } else {
138 printf("Got mapping object from %s\n", MAPPING_FILE);
139 }
140
141 AliTPCCalibPedestal calibPedestal; // pedestal and noise calibration
142 calibPedestal.SetRangeTime(60,940); // set time bin range
143 calibPedestal.SetTimeAnalysis(timeAnalysis); // pedestal(t) calibration
b401648b 144 if (mapping){
145 calibPedestal.SetAltroMapping(mapping->GetAltroMapping()); // Use altro mapping we got from daqDetDb
146 }
c12208b8 147 /* loop over RAW data files */
148 int nevents=0;
bdf99a93 149 for ( i=1; i<argc; i++ ) {
c12208b8 150
151 /* define data source : this is argument i */
152 printf("Processing file %s\n", argv[i]);
153 status=monitorSetDataSource( argv[i] );
154 if (status!=0) {
49efac78 155 printf("monitorSetDataSource() failed. Error=%s. Exiting ...\n", monitorDecodeError(status));
c12208b8 156 return -1;
157 }
158
159 /* read until EOF */
bdf99a93 160 for ( ; ; ) {
c12208b8 161 struct eventHeaderStruct *event;
162
163 /* check shutdown condition */
164 if (daqDA_checkShutdown()) {break;}
165
166 /* get next event (blocking call until timeout) */
167 status=monitorGetEventDynamic((void **)&event);
168 if (status==MON_ERR_EOF) {
bdf99a93 169 printf ("End of File %d (%s) detected\n", i, argv[i]);
170 break; /* end of monitoring file has been reached */
c12208b8 171 }
c12208b8 172 if (status!=0) {
bdf99a93 173 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
174 break;
c12208b8 175 }
176
bdf99a93 177 /* skip start/end of run events */
178 if ( (event->eventType != physicsEvent) && (event->eventType != calibrationEvent) )
179 continue;
180
c12208b8 181 /* retry if got no event */
49efac78 182 if (event==NULL)
bdf99a93 183 continue;
bdf99a93 184
c12208b8 185 nevents++;
186
187 // Pedestal calibration
188 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
b401648b 189 calibPedestal.ProcessEvent(rawReader);
190 //calibPedestal.ProcessEventFast(rawReader); // fast data reader
c12208b8 191 delete rawReader;
192
193 /* free resources */
194 free(event);
195 }
196 }
197
49efac78 198 //
199 // Analyse pedestals and write them to rootfile
200 //
201
202 calibPedestal.Analyse();
203 calibPedestal.AnalyseTime(nevents);
204 printf ("%d physics/calibration events processed.\n",nevents);
c12208b8 205
bdf99a93 206 TFile *fileTPC = new TFile(RESULT_FILE, "recreate");
e73181c9 207 calibPedestal.Write("tpcCalibPedestal");
c12208b8 208 delete fileTPC;
49efac78 209 printf("Wrote %s.\n",RESULT_FILE);
f2c72763 210
49efac78 211 /* store the result file on FES */
212
bd955ed2 213 status=daqDA_FES_storeFile(RESULT_FILE,FILE_ID);
49efac78 214 if (status) {
215 status = -2;
216 }
f2c72763 217
218
bdf99a93 219 //
49efac78 220 // Now prepare ASCII files for local ALTRO configuration through DDL.
bdf99a93 221 //
87d57858 222
bdf99a93 223 ofstream pedfile;
224 ofstream noisefile;
225 ofstream pedmemfile;
87d57858 226 char filename[255];
bdf99a93 227 sprintf(filename,"tpcPedestals.data");
228 pedfile.open(filename);
229 sprintf(filename,"tpcNoise.data");
230 noisefile.open(filename);
231 sprintf(filename,"tpcPedestalMem.data");
232 pedmemfile.open(filename);
233
49efac78 234 TArrayF **timePed = calibPedestal.GetTimePedestals(); // pedestal values for each time bin
87d57858 235
bdf99a93 236 Int_t ctr_channel = 0;
237 Int_t ctr_altro = 0;
238 Int_t ctr_pattern = 0;
87d57858 239
49efac78 240 pedfile << 10 << std::endl; // mark file to contain PEDESTALS per channel
241 noisefile << 11 << std::endl; // mark file to contain NOISE per altro
242 pedmemfile << 12 << std::endl; // mark file to contain PEDESTALs per time bin
87d57858 243
bdf99a93 244 for ( Int_t roc = 0; roc < 72; roc++ ) {
245 if ( !calibPedestal.GetCalRocPedestal(roc) ) continue;
49efac78 246 Int_t side = mapping->GetSideFromRoc(roc);
247 Int_t sector = mapping->GetSectorFromRoc(roc);
bdf99a93 248 //printf("Analysing ROC %d (side %d, sector %d) ...\n", roc, side, sector);
49efac78 249 Int_t nru = mapping->IsIROC(roc) ? 2 : 4;
bdf99a93 250 for ( int rcu = 0; rcu < nru; rcu++ ) {
49efac78 251 Int_t patch = mapping->IsIROC(roc) ? rcu : rcu+2;
bdf99a93 252 for ( int branch = 0; branch < 2; branch++ ) {
49efac78 253 for ( int fec = 0; fec < mapping->GetNfec(patch, branch); fec++ ) {
bdf99a93 254 for ( int altro = 0; altro < 8; altro++ ) {
255 Float_t rms = 0.;
256 Float_t ctr = 0.;
257 for ( int channel = 0; channel < 16; channel++ ) {
49efac78 258 Int_t hwadd = mapping->CodeHWAddress(branch, fec, altro, channel);
259 Int_t row = mapping->GetPadRow(patch, hwadd); // row in a ROC
260 Int_t globalrow = mapping->GetGlobalPadRow(patch, hwadd); // row in full sector
261 Int_t pad = mapping->GetPad(patch, hwadd);
262 Float_t ped = calibPedestal.GetCalRocPedestal(roc)->GetValue(row,pad);
bdf99a93 263 // fixed pedestal
264 if ( ped > 1.e-10 ) {
49efac78 265 pedfile << ctr_channel << "\t" << side << "\t" << sector << "\t" << patch << "\t"
266 << hwadd << "\t" << ped << std::endl;
bdf99a93 267 ctr_channel++;
268 }
269 // pedestal(t)
49efac78 270 if ( timePed && fabs(timePed[globalrow][pad].GetSum()) > 1e-10 ) {
271 pedmemfile << ctr_pattern << "\t" << side << "\t" << sector << "\t" << patch
272 << "\t" << hwadd;
bdf99a93 273 for ( Int_t timebin = 0; timebin < 1024; timebin++ )
49efac78 274 pedmemfile << "\t" << timePed[globalrow][pad].At(timebin);
bdf99a93 275 pedmemfile << std::endl;
276 ctr_pattern++;
277 }
278 // rms=noise
279 Float_t rms2 = calibPedestal.GetCalRocRMS(roc)->GetValue(row,pad);
280 if ( rms2 > 1.e-10 ) { rms += rms2; ctr += 1.; }
281 } // end channel for loop
282 // noise data (rms) averaged over all channels in this ALTRO.
49efac78 283 Int_t hwadd = mapping->CodeHWAddress(branch, fec, altro, 0);
bdf99a93 284 if ( ctr > 1.e-10 ) {
49efac78 285 noisefile << ctr_altro << "\t" << side << "\t" << sector << "\t" << patch << "\t"
286 << hwadd << "\t" << rms/ctr << std::endl;
bdf99a93 287 ctr_altro++;
288 }
289 } // end altro for loop
290 } // end fec for loop
291 } // end branch for loop
292 } // end rcu for loop
293 } // end roc loop
294
295 pedfile.close();
296 noisefile.close();
297 pedmemfile.close();
49efac78 298 printf("Wrote ASCII files.\n");
bdf99a93 299
87d57858 300
c12208b8 301 return status;
302}