2 TPC DA for online calibration
4 Contact: Haavard.Helstrup@cern.ch
8 Number of events needed: 100
10 Output Files: tpcPedestal.root, to be exported to the DAQ FXS
12 Trigger types used: CALIBRATION_EVENT
18 TPCda_pedestal.cxx - calibration algorithm for TPC pedestal runs
20 10/06/2007 sylvain.chapeland@cern.ch : first version - clean skeleton based on DAQ DA case1
21 19/10/2007 christian.lippmann@cern.ch : Possibility to write output to ASCII file
22 24/10/2007 christian.lippmann@cern.ch : Including pedestal calibration for time bins
23 23/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.
25 28/11/2007 christian.lippmann@cern.ch : TPC mapping file is read from DaqDetDB
26 18/09/2008 christian.lippmann@cern.ch : Noisy channels are output to ASCII file. Use max noise in ALTRO.
27 19/09/2008 J.Wiechula@gsi.de: Added export of the calibration data to the AMORE data base.
28 Added support for configuration files.
30 contact: marian.ivanov@cern.ch
32 This process reads RAW data from the files provided as command line arguments
33 and save results in a file (named from RESULT_FILE define - see below).
37 #define RESULT_FILE "tpcPedestal.root"
38 #define FILE_ID "pedestals"
39 #define MAPPING_FILE "tpcMapping.root"
40 #define CONFIG_FILE "TPCPEDESTALda.conf"
41 #define AliDebugLevel() -1
59 #include "TPluginManager.h"
62 #include "TObjString.h"
67 #include "AliRawReader.h"
68 #include "AliRawReaderDate.h"
69 #include "AliTPCmapper.h"
70 #include "AliTPCRawStream.h"
71 #include "AliTPCROC.h"
72 #include "AliTPCCalROC.h"
73 #include "AliTPCCalPad.h"
74 #include "AliMathBase.h"
75 #include "TTreeStream.h"
77 #include "AliTPCConfigDA.h"
85 // TPC calibration algorithm includes
87 #include "AliTPCCalibPedestal.h"
90 Main routine, TPC pedestal detector algorithm to be run on TPC LDC
91 Arguments: list of DATE raw data files
94 int main(int argc, char **argv) {
96 // Main for TPC pedestal detector algorithm
98 /* log start of process */
99 printf("TPC DA started - %s\n",__FILE__);
102 printf("Wrong number of arguments\n");
105 AliLog::SetClassDebugLevel("AliTPCRawStream",-5);
106 AliLog::SetClassDebugLevel("AliRawReaderDate",-5);
107 AliLog::SetClassDebugLevel("AliTPCAltroMapping",-5);
108 AliLog::SetModuleDebugLevel("RAW",-5);
111 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
117 /* declare monitoring program */
119 status=monitorDeclareMp( __FILE__ );
121 printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
126 AliTPCmapper *mapping = 0; // The TPC mapping
128 unsigned long32 runNb=0; // run number
129 // configuration options
130 Bool_t timeAnalysis = kTRUE;
131 Bool_t fastDecoding = kFALSE;
134 /* copy locally the mapping file from daq detector config db */
135 sprintf(localfile,"./%s",MAPPING_FILE);
136 status = daqDA_DB_getFile(MAPPING_FILE,localfile);
138 printf("Failed to get mapping file (%s) from DAQdetDB, status=%d\n", MAPPING_FILE, status);
142 /* open the mapping file and retrieve mapping object */
143 TFile *fileMapping = new TFile(MAPPING_FILE, "read");
144 mapping = (AliTPCmapper*) fileMapping->Get("tpcMapping");
149 printf("Failed to get mapping object from %s. ...\n", MAPPING_FILE);
152 printf("Got mapping object from %s\n", MAPPING_FILE);
156 // DA configuration from configuration file
158 // retrieve configuration file
159 sprintf(localfile,"./%s",CONFIG_FILE);
160 status = daqDA_DB_getFile(CONFIG_FILE,localfile);
162 printf("Failed to get configuration file (%s) from DAQdetDB, status=%d\n", CONFIG_FILE, status);
165 AliTPCConfigDA config(CONFIG_FILE);
166 // check configuration
167 if ( (Int_t)config.GetValue("NoTimeAnalysis") == 1 ) {
168 printf("WARNING: Time analysis was switched off in the configuration file!\n");
171 if ( (Int_t)config.GetValue("UseFastDecoder") == 1 ){
172 printf("Info: The fast decoder will be used for the processing.\n");
176 // create calibration object
177 AliTPCCalibPedestal calibPedestal(config.GetConfigurationMap()); // pedestal and noise calibration
178 calibPedestal.SetAltroMapping(mapping->GetAltroMapping()); // Use altro mapping we got from daqDetDb
180 //===========================//
181 // loop over RAW data files //
182 //==========================//
184 for ( i=1; i<argc; i++ ) {
186 /* define data source : this is argument i */
187 printf("Processing file %s\n", argv[i]);
188 status=monitorSetDataSource( argv[i] );
190 printf("monitorSetDataSource() failed. Error=%s. Exiting ...\n", monitorDecodeError(status));
196 struct eventHeaderStruct *event;
198 /* check shutdown condition */
199 if (daqDA_checkShutdown()) {break;}
201 /* get next event (blocking call until timeout) */
202 status=monitorGetEventDynamic((void **)&event);
203 if (status==MON_ERR_EOF) {
204 printf ("End of File %d (%s) detected\n", i, argv[i]);
205 break; /* end of monitoring file has been reached */
209 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
213 /* skip start/end of run events */
214 if ( (event->eventType != physicsEvent) && (event->eventType != calibrationEvent) )
217 /* retry if got no event */
222 // get the run number
223 runNb = event->eventRunNb;
224 // Pedestal calibration
225 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
226 if ( fastDecoding ) calibPedestal.ProcessEventFast(rawReader);
227 else calibPedestal.ProcessEvent(rawReader);
236 // Analyse pedestals and write them to rootfile
238 calibPedestal.Analyse();
239 calibPedestal.AnalyseTime(nevents);
240 printf ("%d physics/calibration events processed.\n",nevents);
242 TFile *fileTPC = new TFile(RESULT_FILE, "recreate");
243 calibPedestal.Write("tpcCalibPedestal");
245 printf("Wrote %s.\n",RESULT_FILE);
247 /* store the result file on FES */
248 status=daqDA_FES_storeFile(RESULT_FILE,FILE_ID);
253 //Send objects to the AMORE DB
255 printf ("AMORE part\n");
256 const char *amoreDANameorig=gSystem->Getenv("AMORE_DA_NAME");
257 //cheet a little -- temporary solution (hopefully)
259 //currently amoreDA uses the environment variable AMORE_DA_NAME to create the mysql
260 //table in which the calib objects are stored. This table is dropped each time AmoreDA
261 //is initialised. This of course makes a problem if we would like to store different
262 //calibration entries in the AMORE DB. Therefore in each DA which writes to the AMORE DB
263 //the AMORE_DA_NAME env variable is overwritten.
265 //find processed sector
268 for ( Int_t roc = 0; roc < 72; roc++ ) {
269 if ( !calibPedestal.GetCalRocPedestal(roc) ) continue;
270 if (mapping->GetSideFromRoc(roc)==1) sideName='C';
271 sector = mapping->GetSectorFromRoc(roc);
273 gSystem->Setenv("AMORE_DA_NAME",Form("TPC-%c%02d-%s",sideName,sector,FILE_ID));
278 TObjString info(Form("Run: %u; Date: %s",runNb,time.AsString()));
280 amore::da::AmoreDA amoreDA(amore::da::AmoreDA::kSender);
282 statusDA+=amoreDA.Send("Pedestals",calibPedestal.GetCalPadPedestal());
283 statusDA+=amoreDA.Send("Noise",calibPedestal.GetCalPadRMS());
284 statusDA+=amoreDA.Send("Info",&info);
286 printf("Waring: Failed to write one of the calib objects to the AMORE database\n");
288 if (amoreDANameorig) gSystem->Setenv("AMORE_DA_NAME",amoreDANameorig);
292 // Now prepare ASCII files for local ALTRO configuration through DDL.
297 ofstream noisychannelfile;
300 sprintf(filename,"tpcPedestals.data");
301 pedfile.open(filename);
302 sprintf(filename,"tpcNoise.data");
303 noisefile.open(filename);
304 sprintf(filename,"tpcPedestalMem.data");
305 pedmemfile.open(filename);
306 sprintf(filename,"tpcNoisyChannels.data");
307 noisychannelfile.open(filename);
309 TArrayF **timePed = calibPedestal.GetTimePedestals(); // pedestal values for each time bin
311 Int_t ctr_channel = 0;
313 Int_t ctr_pattern = 0;
316 pedfile << 10 << std::endl; // Mark file to contain PEDESTALS per channel
317 noisefile << 11 << std::endl; // Mark file to contain NOISE per altro
318 pedmemfile << 12 << std::endl; // Mark file to contain PEDESTALs per time bin
319 noisychannelfile << 14 << std::endl; // Mark file to contain NOISY CHANNELS
321 for ( Int_t roc = 0; roc < 72; roc++ ) {
322 if ( !calibPedestal.GetCalRocPedestal(roc) ) continue;
323 Int_t side = mapping->GetSideFromRoc(roc);
324 Int_t sector = mapping->GetSectorFromRoc(roc);
325 //printf("Analysing ROC %d (side %d, sector %d) ...\n", roc, side, sector);
326 Int_t nru = mapping->IsIROC(roc) ? 2 : 4;
327 for ( int rcu = 0; rcu < nru; rcu++ ) {
328 Int_t patch = mapping->IsIROC(roc) ? rcu : rcu+2;
329 for ( int branch = 0; branch < 2; branch++ ) {
330 for ( int fec = 0; fec < mapping->GetNfec(patch, branch); fec++ ) {
331 for ( int altro = 0; altro < 8; altro++ ) {
334 Float_t ctr_altrochannel = 0.;
335 for ( int channel = 0; channel < 16; channel++ ) {
336 Int_t hwadd = mapping->CodeHWAddress(branch, fec, altro, channel);
337 Int_t row = mapping->GetPadRow(patch, hwadd); // row in a ROC
338 Int_t globalrow = mapping->GetGlobalPadRow(patch, hwadd); // row in full sector
339 Int_t pad = mapping->GetPad(patch, hwadd);
340 Float_t ped = calibPedestal.GetCalRocPedestal(roc)->GetValue(row,pad);
342 pedfile << ctr_channel++ << "\t" << side << "\t" << sector << "\t" << patch << "\t"
343 << hwadd << "\t" << ped << std::endl;
345 if ( timePed && fabs(timePed[globalrow][pad].GetSum()) > 1e-10 ) {
346 pedmemfile << ctr_pattern++ << "\t" << side << "\t" << sector << "\t" << patch
348 for ( Int_t timebin = 0; timebin < 1024; timebin++ )
349 pedmemfile << "\t" << timePed[globalrow][pad].At(timebin);
350 pedmemfile << std::endl;
353 Float_t rms2 = calibPedestal.GetCalRocRMS(roc)->GetValue(row,pad);
354 if ( fabs(ped) < 1.e-10 ) { // dead channel
355 noisychannelfile << ctr_noisy++ << "\t" << side << "\t" << sector << "\t"
356 << patch << "\t" << hwadd << "\t" << rms2 << std::endl;
357 } else if ( (ped > 1.e-10) && (rms2 > 1.e-10) ) { // not dead
358 // Find noisy channels
359 if ( ((roc<36) && (rms2 > 2.0)) || // IROC
360 ((roc>35) && (row<65) && (rms2 > 2.0)) || // OROC, small pads
361 ((roc>35) && (row>64) && (rms2 > 3.0)) ) { // OROC, large pads (50% more signal)
362 noisychannelfile << ctr_noisy++ << "\t" << side << "\t" << sector << "\t"
363 << patch << "\t" << hwadd << "\t" << rms2 << std::endl;
365 // Not noisy. Get average and maximum noise in this ALTRO
367 ctr_altrochannel += 1.;
368 if (rms2 > maxrms) maxrms = rms2;
370 } // end if some signal
371 } // end channel for loop
372 Int_t hwadd = mapping->CodeHWAddress(branch, fec, altro, 0); // ALTRO address
373 // Noise data (rms) averaged over all channels in this ALTRO.
374 if ( ctr_altrochannel > 1.e-10 ) {
376 // average noise of this ALTRO (excluding high-noise channels)
377 noisefile << ctr_altro << "\t" << side << "\t" << sector << "\t" << patch << "\t"
378 << hwadd << "\t" << rms/ctr_altrochannel << std::endl;
380 // maximum noise of this ALTRO (excluding high-noise channels)
381 noisefile << ctr_altro << "\t" << side << "\t" << sector << "\t" << patch << "\t"
382 << hwadd << "\t" << maxrms << std::endl;
385 } // end altro for loop
386 } // end fec for loop
387 } // end branch for loop
388 } // end rcu for loop
394 noisychannelfile.close();
395 printf("Wrote ASCII files. Found %d noisy channels.\n", ctr_noisy);