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