]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/TPCPEDESTALda.cxx
Fixed warning.
[u/mrichter/AliRoot.git] / TPC / TPCPEDESTALda.cxx
1 /*
2 TPC DA for online calibration
3
4 Contact: Haavard.Helstrup@cern.ch
5 Link:
6 Run Type: PEDESTAL
7 DA Type: LDC
8 Number of events needed: 100
9 Input Files: 
10 Output Files: tpcPedestal.root, to be exported to the DAQ FXS
11 fileId:   pedestals    
12 Trigger types used: CALIBRATION_EVENT
13
14 */
15
16 /*
17
18 TPCda_pedestal.cxx - calibration algorithm for TPC pedestal runs
19
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.
29
30 contact: marian.ivanov@cern.ch
31
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).
34
35 */
36
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
42
43 extern "C" {
44 #include <daqDA.h>
45 }
46 #include "event.h"
47 #include "monitor.h"
48
49 #include "stdio.h"
50 #include "stdlib.h"
51 #include <fstream>
52
53 //
54 //Root includes
55 //
56 #include "TFile.h"
57 #include "TArrayF.h"
58 #include "TROOT.h"
59 #include "TPluginManager.h"
60 #include "TSystem.h"
61 #include "TString.h"
62 #include "TObjString.h"
63 #include "TDatime.h"
64 //
65 //AliRoot includes
66 //
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"
76 #include "AliLog.h"
77 #include "AliTPCConfigDA.h"
78
79 //
80 //AMORE
81 //
82 #include <AmoreDA.h>
83
84 //
85 // TPC calibration algorithm includes
86 //
87 #include "AliTPCCalibPedestal.h"
88
89 /*
90   Main routine, TPC pedestal detector algorithm to be run on TPC LDC
91   Arguments: list of DATE raw data files
92 */
93
94 int main(int argc, char **argv) {
95   //
96   // Main for TPC pedestal detector algorithm
97   //
98   /* log start of process */
99   printf("TPC DA started - %s\n",__FILE__);
100   
101   if (argc<2) {
102     printf("Wrong number of arguments\n");
103     return -1;
104   }
105   AliLog::SetClassDebugLevel("AliTPCRawStream",-5);
106   AliLog::SetClassDebugLevel("AliRawReaderDate",-5);
107   AliLog::SetClassDebugLevel("AliTPCAltroMapping",-5);
108   AliLog::SetModuleDebugLevel("RAW",-5);
109
110   /* magic line */
111   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
112           "*",
113           "TStreamerInfo",
114           "RIO",
115           "TStreamerInfo()");
116
117   /* declare monitoring program */
118   int i, status;
119   status=monitorDeclareMp( __FILE__ );
120   if (status!=0) {
121     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
122     return -1;
123   }
124   
125   // variables
126   AliTPCmapper *mapping = 0;   // The TPC mapping
127   char localfile[255];
128   unsigned long32 runNb=0;     // run number
129   // configuration options 
130   Bool_t timeAnalysis = kTRUE;
131   Bool_t fastDecoding = kFALSE;
132
133   if (!mapping){
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);
137     if (status) {
138       printf("Failed to get mapping file (%s) from DAQdetDB, status=%d\n", MAPPING_FILE, status);
139       return -1;
140     }
141
142     /* open the mapping file and retrieve mapping object */
143     TFile *fileMapping = new TFile(MAPPING_FILE, "read");
144     mapping = (AliTPCmapper*) fileMapping->Get("tpcMapping");
145     delete fileMapping;
146   }
147
148   if (mapping == 0) {
149     printf("Failed to get mapping object from %s.  ...\n", MAPPING_FILE);
150     return -1;
151   } else {
152     printf("Got mapping object from %s\n", MAPPING_FILE);
153   }
154
155   //
156   // DA configuration from configuration file
157   //
158   // retrieve configuration file
159   sprintf(localfile,"./%s",CONFIG_FILE);
160   status = daqDA_DB_getFile(CONFIG_FILE,localfile);
161   if (status) {
162     printf("Failed to get configuration file (%s) from DAQdetDB, status=%d\n", CONFIG_FILE, status);
163     return -1;
164   }
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");
169     timeAnalysis=kFALSE;
170   }
171   if ( (Int_t)config.GetValue("UseFastDecoder") == 1 ){
172     printf("Info: The fast decoder will be used for the processing.\n");
173     fastDecoding=kTRUE;
174   }
175   
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
179   calibPedestal.SetTimeAnalysis(timeAnalysis);               // pedestal(t) calibration 
180   
181   //===========================//
182   // loop over RAW data files //
183   //==========================//
184   int nevents=0;
185   for ( i=1; i<argc; i++ ) {
186     
187     /* define data source : this is argument i */
188     printf("Processing file %s\n", argv[i]);
189     status=monitorSetDataSource( argv[i] );
190     if (status!=0) {
191       printf("monitorSetDataSource() failed. Error=%s. Exiting ...\n", monitorDecodeError(status));
192       return -1;
193     }
194     
195     /* read until EOF */
196     for ( ; ; ) {
197       struct eventHeaderStruct *event;
198       
199       /* check shutdown condition */
200       if (daqDA_checkShutdown()) {break;}
201       
202       /* get next event (blocking call until timeout) */
203       status=monitorGetEventDynamic((void **)&event);
204       if (status==MON_ERR_EOF) {
205         printf ("End of File %d (%s) detected\n", i, argv[i]);
206         break; /* end of monitoring file has been reached */
207       }
208       
209       if (status!=0) {
210         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
211         break;
212       }
213
214       /* retry if got no event */
215       if (event==NULL)
216         continue;
217       
218       /* skip start/end of run events */
219       if ( (event->eventType != physicsEvent) && (event->eventType != calibrationEvent) )
220         continue;
221
222       nevents++;
223       // get the run number
224       runNb = event->eventRunNb;
225       //  Pedestal calibration
226       calibPedestal.ProcessEvent(event);
227
228       /* free resources */
229       free(event);
230     }
231   }
232
233   //
234   // Analyse pedestals and write them to rootfile
235   //
236   calibPedestal.Analyse();
237   calibPedestal.AnalyseTime(nevents);
238   printf ("%d physics/calibration events processed.\n",nevents);
239
240   TFile *fileTPC = new TFile(RESULT_FILE, "recreate");
241   calibPedestal.Write("tpcCalibPedestal");
242   delete fileTPC;
243   printf("Wrote %s.\n",RESULT_FILE);
244
245   /* store the result file on FES */
246   status=daqDA_FES_storeFile(RESULT_FILE,FILE_ID);
247   if (status) {
248     status = -2;
249   }
250   //
251   //Send objects to the AMORE DB
252   //
253   printf ("AMORE part\n");
254   const char *amoreDANameorig=gSystem->Getenv("AMORE_DA_NAME");
255   //cheet a little -- temporary solution (hopefully)
256   //
257   //currently amoreDA uses the environment variable AMORE_DA_NAME to create the mysql
258   //table in which the calib objects are stored. This table is dropped each time AmoreDA
259   //is initialised. This of course makes a problem if we would like to store different
260   //calibration entries in the AMORE DB. Therefore in each DA which writes to the AMORE DB
261   //the AMORE_DA_NAME env variable is overwritten.  
262
263   //find processed sector
264   Char_t sideName='A';
265   Int_t sector = -1;
266   for ( Int_t roc = 0; roc < 72; roc++ ) {
267     if ( !calibPedestal.GetCalRocPedestal(roc) ) continue;
268     if (mapping->GetSideFromRoc(roc)==1) sideName='C';
269     sector = mapping->GetSectorFromRoc(roc);
270   }
271   gSystem->Setenv("AMORE_DA_NAME",Form("TPC-%c%02d-%s",sideName,sector,FILE_ID));
272   
273   //
274   // end cheet
275   if (sector>-1){
276     TDatime time;
277     TObjString info(Form("Run: %u; Date: %s",runNb,time.AsSQLString()));
278     
279     amore::da::AmoreDA amoreDA(amore::da::AmoreDA::kSender);
280     Int_t statusDA=0;
281     statusDA+=amoreDA.Send("Pedestals",calibPedestal.GetCalPadPedestal());
282     statusDA+=amoreDA.Send("Noise",calibPedestal.GetCalPadRMS());
283     statusDA+=amoreDA.Send("Info",&info);
284     if ( statusDA )
285       printf("Waring: Failed to write one of the calib objects to the AMORE database\n");
286   }  else {
287     printf("Waring: No data found!\n");
288   }
289   // reset env var
290   if (amoreDANameorig) gSystem->Setenv("AMORE_DA_NAME",amoreDANameorig);
291   
292   //
293   // Now prepare ASCII files for local ALTRO configuration through DDL.
294   //
295   ofstream pedfile;
296   ofstream noisefile;
297   ofstream pedmemfile;
298   ofstream noisychannelfile;
299
300   char filename[255];
301   sprintf(filename,"tpcPedestals.data");
302   pedfile.open(filename);
303   sprintf(filename,"tpcNoise.data");
304   noisefile.open(filename);
305   sprintf(filename,"tpcPedestalMem.data");
306   pedmemfile.open(filename);
307   sprintf(filename,"tpcNoisyChannels.data");
308   noisychannelfile.open(filename);
309
310   TArrayF **timePed = calibPedestal.GetTimePedestals();  // pedestal values for each time bin
311
312   Int_t ctr_channel = 0;
313   Int_t ctr_altro   = 0;
314   Int_t ctr_pattern = 0;
315   Int_t ctr_noisy   = 0;
316
317   pedfile          << 10 << std::endl; // Mark file to contain PEDESTALS per channel
318   noisefile        << 11 << std::endl; // Mark file to contain NOISE per altro
319   pedmemfile       << 12 << std::endl; // Mark file to contain PEDESTALs per time bin
320   noisychannelfile << 14 << std::endl; // Mark file to contain NOISY CHANNELS
321
322   for ( Int_t roc = 0; roc < 72; roc++ ) {
323     if ( !calibPedestal.GetCalRocPedestal(roc) ) continue;
324     Int_t side   = mapping->GetSideFromRoc(roc);
325     Int_t sector = mapping->GetSectorFromRoc(roc);
326     //printf("Analysing ROC %d (side %d, sector %d) ...\n", roc, side, sector);
327     Int_t nru = mapping->IsIROC(roc) ? 2 : 4;
328     for ( int rcu = 0; rcu < nru; rcu++ ) {
329       Int_t patch = mapping->IsIROC(roc) ? rcu : rcu+2;
330       for ( int branch = 0; branch < 2; branch++ ) {
331         for ( int fec = 0; fec < mapping->GetNfec(patch, branch); fec++ ) {
332           for ( int altro = 0; altro < 8; altro++ ) {
333             Float_t rms = 0.;
334             Float_t maxrms = 0.;
335             Float_t ctr_altrochannel = 0.;
336             for ( int channel = 0; channel < 16; channel++ ) {
337               Int_t hwadd     = mapping->CodeHWAddress(branch, fec, altro, channel);
338               Int_t row       = mapping->GetPadRow(patch, hwadd);        // row in a ROC
339               Int_t globalrow = mapping->GetGlobalPadRow(patch, hwadd);  // row in full sector
340               Int_t pad       = mapping->GetPad(patch, hwadd);
341               Float_t ped     = calibPedestal.GetCalRocPedestal(roc)->GetValue(row,pad);
342               // fixed pedestal
343               pedfile << ctr_channel++ << "\t" << side << "\t" << sector << "\t" << patch << "\t"
344                   << hwadd << "\t" << ped << std::endl;
345               // pedestal(t)
346               if ( timePed && fabs(timePed[globalrow][pad].GetSum()) > 1e-10 ) {
347                 pedmemfile << ctr_pattern++ << "\t" << side << "\t" << sector << "\t" << patch
348                     << "\t" << hwadd;
349                 for ( Int_t timebin = 0; timebin < 1024; timebin++ )
350                   pedmemfile << "\t" << timePed[globalrow][pad].At(timebin);
351                 pedmemfile << std::endl;
352               }
353               // rms=noise
354               Float_t rms2 = calibPedestal.GetCalRocRMS(roc)->GetValue(row,pad);
355               if ( fabs(ped) < 1.e-10 ) {                        // dead channel
356                 noisychannelfile << ctr_noisy++ << "\t" << side << "\t" << sector << "\t"
357                     << patch << "\t" << hwadd << "\t" << rms2 << std::endl;
358               } else if ( (ped > 1.e-10) && (rms2 > 1.e-10) ) {  // not dead
359               // Find noisy channels
360                 if ( ((roc<36)             && (rms2 > 2.0))  ||  // IROC
361                        ((roc>35) && (row<65) && (rms2 > 2.0))  ||  // OROC, small pads
362                        ((roc>35) && (row>64) && (rms2 > 3.0)) ) {  // OROC, large pads (50% more signal)
363                   noisychannelfile << ctr_noisy++ << "\t" << side << "\t" << sector << "\t"
364                       << patch << "\t" << hwadd << "\t" << rms2 << std::endl;
365                        } else {
366                        // Not noisy. Get average and maximum noise in this ALTRO
367                          rms += rms2;
368                          ctr_altrochannel += 1.;
369                          if (rms2 > maxrms) maxrms = rms2;
370                        } // end if noisy
371               } // end if some signal
372             } // end channel for loop
373             Int_t hwadd = mapping->CodeHWAddress(branch, fec, altro, 0);   // ALTRO address
374       // Noise data (rms) averaged over all channels in this ALTRO.
375             if ( ctr_altrochannel > 1.e-10 ) {
376         /*
377         // average noise of this ALTRO (excluding high-noise channels)
378               noisefile << ctr_altro << "\t" << side << "\t" << sector << "\t" << patch << "\t"
379               << hwadd << "\t" << rms/ctr_altrochannel << std::endl;
380         */
381         // maximum noise of this ALTRO (excluding high-noise channels)
382               noisefile << ctr_altro << "\t" << side << "\t" << sector << "\t" << patch << "\t"
383                   << hwadd << "\t" << maxrms << std::endl;
384               ctr_altro++;
385             }
386           } // end altro for loop
387         } // end fec for loop
388       } // end branch for loop
389     } // end rcu for loop
390   } // end roc loop
391   
392   pedfile.close();
393   noisefile.close();
394   pedmemfile.close();
395   noisychannelfile.close();
396   printf("Wrote ASCII files. Found %d noisy channels.\n", ctr_noisy);
397
398   return status;
399
400 }