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