]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/TPCPEDESTALda.cxx
47159a14bc9f9a20f2f153f5a2c3d867e88eef62
[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 #define AliDebugLevel() -1
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 #include "AliLog.h"
68 #include "TSystem.h"
69
70 //
71 // TPC calibration algorithm includes
72 //
73 #include "AliTPCCalibPedestal.h"
74
75 /*
76   Main routine, TPC pedestal detector algorithm to be run on TPC LDC
77   Arguments: list of DATE raw data files
78 */
79
80 int main(int argc, char **argv) {
81   //
82   // Main for TPC pedestal detector algorithm
83   //
84
85   AliLog::SetClassDebugLevel("AliTPCRawStream",AliLog::kFatal);
86   AliLog::SetClassDebugLevel("AliRawReaderDate",AliLog::kFatal);
87   AliLog::SetClassDebugLevel("AliTPCAltroMapping",AliLog::kFatal);
88   AliLog::SetModuleDebugLevel("TPC",AliLog::kFatal);
89   AliLog::SetModuleDebugLevel("RAW",AliLog::kFatal);
90
91   Bool_t timeAnalysis = kTRUE;
92
93   if (argc<2) {
94     printf("Wrong number of arguments\n");
95     return -1;
96   }
97
98   /* magic line */
99   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
100                                         "*",
101                                         "TStreamerInfo",
102                                         "RIO",
103                                         "TStreamerInfo()"); 
104   int i, status;
105
106   /* log start of process */
107   printf("TPC DA started - %s\n",__FILE__);
108
109   /* declare monitoring program */
110   status=monitorDeclareMp( __FILE__ );
111   if (status!=0) {
112     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
113     return -1;
114   }
115
116   AliTPCmapper *mapping = 0;   // The TPC mapping
117    
118   if (!mapping){
119     /* copy locally the mapping file from daq detector config db */
120     status = daqDA_DB_getFile(MAPPING_FILE,"./tpcMapping.root");
121     if (status) {
122       printf("Failed to get mapping file (%s) from DAQdetDB, status=%d\n", MAPPING_FILE, status);
123       printf("Continue anyway ... maybe it works?\n");              // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
124       return -1;   // temporarily uncommented for testing on pcald47 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
125     }
126
127     /* open the mapping file and retrieve mapping object */
128     TFile *fileMapping = new TFile(MAPPING_FILE, "read");
129     mapping = (AliTPCmapper*) fileMapping->Get("tpcMapping");
130     delete fileMapping;
131   }
132
133   if (mapping == 0) {
134     printf("Failed to get mapping object from %s.  ...\n", MAPPING_FILE);
135     //return -1;
136   } else {
137     printf("Got mapping object from %s\n", MAPPING_FILE);
138   }
139
140   AliTPCCalibPedestal calibPedestal;                         // pedestal and noise calibration
141   calibPedestal.SetRangeTime(60,940);                        // set time bin range
142   calibPedestal.SetTimeAnalysis(timeAnalysis);               // pedestal(t) calibration
143   if (mapping){
144     calibPedestal.SetAltroMapping(mapping->GetAltroMapping()); // Use altro mapping we got from daqDetDb
145   }
146   /* loop over RAW data files */
147   int nevents=0;
148   for ( i=1; i<argc; i++ ) {
149
150     /* define data source : this is argument i */
151     printf("Processing file %s\n", argv[i]);
152     status=monitorSetDataSource( argv[i] );
153     if (status!=0) {
154       printf("monitorSetDataSource() failed. Error=%s. Exiting ...\n", monitorDecodeError(status));
155       return -1;
156     }
157
158     /* read until EOF */
159     for ( ; ; ) {
160       struct eventHeaderStruct *event;
161
162       /* check shutdown condition */
163       if (daqDA_checkShutdown()) {break;}
164
165       /* get next event (blocking call until timeout) */
166       status=monitorGetEventDynamic((void **)&event);
167       if (status==MON_ERR_EOF) {
168         printf ("End of File %d (%s) detected\n", i, argv[i]);
169         break; /* end of monitoring file has been reached */
170       }
171       if (status!=0) {
172         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
173         break;
174       }
175
176       /* skip start/end of run events */
177       if ( (event->eventType != physicsEvent) && (event->eventType != calibrationEvent) )
178         continue;
179
180       /* retry if got no event */
181       if (event==NULL)
182         continue;
183
184       nevents++;
185
186       //  Pedestal calibration
187       AliRawReader *rawReader = new AliRawReaderDate((void*)event);
188       calibPedestal.ProcessEvent(rawReader);
189       //calibPedestal.ProcessEventFast(rawReader);   // fast data reader
190       delete rawReader;
191
192       /* free resources */
193       free(event);
194     }
195   }
196
197   //
198   // Analyse pedestals and write them to rootfile
199   //
200
201   calibPedestal.Analyse();
202   calibPedestal.AnalyseTime(nevents);
203   printf ("%d physics/calibration events processed.\n",nevents);
204
205   TFile *fileTPC = new TFile(RESULT_FILE, "recreate");
206   calibPedestal.Write("calibPedestal");
207   delete fileTPC;
208   printf("Wrote %s.\n",RESULT_FILE);
209
210  /* store the result file on FES */
211  
212    status=daqDA_FES_storeFile(RESULT_FILE,RESULT_FILE);
213    if (status) {
214      status = -2;
215    }
216
217
218   //
219   // Now prepare ASCII files for local ALTRO configuration through DDL.
220   //
221
222   ofstream pedfile;
223   ofstream noisefile;
224   ofstream pedmemfile;
225   char filename[255];
226   sprintf(filename,"tpcPedestals.data");
227   pedfile.open(filename);
228   sprintf(filename,"tpcNoise.data");
229   noisefile.open(filename);
230   sprintf(filename,"tpcPedestalMem.data");
231   pedmemfile.open(filename);
232
233   TArrayF **timePed = calibPedestal.GetTimePedestals();  // pedestal values for each time bin
234
235   Int_t ctr_channel = 0;
236   Int_t ctr_altro = 0;
237   Int_t ctr_pattern = 0;
238
239   pedfile    << 10 << std::endl; // mark file to contain PEDESTALS per channel
240   noisefile  << 11 << std::endl; // mark file to contain NOISE per altro
241   pedmemfile << 12 << std::endl; // mark file to contain PEDESTALs per time bin
242
243   for ( Int_t roc = 0; roc < 72; roc++ ) {
244     if ( !calibPedestal.GetCalRocPedestal(roc) ) continue;
245     Int_t side   = mapping->GetSideFromRoc(roc);
246     Int_t sector = mapping->GetSectorFromRoc(roc);
247     //printf("Analysing ROC %d (side %d, sector %d) ...\n", roc, side, sector);
248     Int_t nru = mapping->IsIROC(roc) ? 2 : 4;
249     for ( int rcu = 0; rcu < nru; rcu++ ) {
250       Int_t patch = mapping->IsIROC(roc) ? rcu : rcu+2;
251       for ( int branch = 0; branch < 2; branch++ ) {
252         for ( int fec = 0; fec < mapping->GetNfec(patch, branch); fec++ ) {
253           for ( int altro = 0; altro < 8; altro++ ) {
254             Float_t rms = 0.;
255             Float_t ctr = 0.;
256             for ( int channel = 0; channel < 16; channel++ ) {
257               Int_t hwadd     = mapping->CodeHWAddress(branch, fec, altro, channel);
258               Int_t row       = mapping->GetPadRow(patch, hwadd);        // row in a ROC
259               Int_t globalrow = mapping->GetGlobalPadRow(patch, hwadd);  // row in full sector
260               Int_t pad       = mapping->GetPad(patch, hwadd);
261               Float_t ped     = calibPedestal.GetCalRocPedestal(roc)->GetValue(row,pad);
262               // fixed pedestal
263               if ( ped > 1.e-10 ) {
264                 pedfile << ctr_channel << "\t" << side << "\t" << sector << "\t" << patch << "\t"
265                         << hwadd << "\t" << ped << std::endl;
266                 ctr_channel++;
267               }
268               // pedestal(t)
269               if ( timePed && fabs(timePed[globalrow][pad].GetSum()) > 1e-10 ) {
270                 pedmemfile << ctr_pattern << "\t" << side << "\t" << sector << "\t" << patch
271                            << "\t" << hwadd;
272                 for ( Int_t timebin = 0; timebin < 1024; timebin++ )
273                   pedmemfile << "\t" << timePed[globalrow][pad].At(timebin);
274                 pedmemfile << std::endl;
275                 ctr_pattern++;
276               }
277               // rms=noise
278               Float_t rms2 = calibPedestal.GetCalRocRMS(roc)->GetValue(row,pad);
279               if ( rms2 > 1.e-10 ) { rms += rms2; ctr += 1.; }
280             } // end channel for loop
281             // noise data (rms) averaged over all channels in this ALTRO.
282             Int_t hwadd = mapping->CodeHWAddress(branch, fec, altro, 0);
283             if ( ctr > 1.e-10 ) {
284               noisefile << ctr_altro << "\t" << side << "\t" << sector << "\t" << patch << "\t"
285                         << hwadd << "\t" << rms/ctr << std::endl;
286               ctr_altro++;
287             }
288           } // end altro for loop
289         } // end fec for loop
290       } // end branch for loop
291     } // end rcu for loop
292   } // end roc loop
293
294   pedfile.close();
295   noisefile.close();
296   pedmemfile.close();
297   printf("Wrote ASCII files.\n");
298
299
300   return status;
301 }