3 TPCda_pedestal.cxx - calibration algorithm for TPC pedestal runs
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
12 contact: marian.ivanov@cern.ch
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).
19 #define RESULT_FILE "tpcPedestal.root"
20 #define MAPPING_FILE "tpcMapping.root"
39 #include "TPluginManager.h"
44 #include "AliRawReader.h"
45 #include "AliRawReaderDate.h"
46 #include "AliTPCmapper.h"
47 #include "AliTPCRawStream.h"
48 #include "AliTPCROC.h"
49 #include "AliTPCCalROC.h"
50 #include "AliTPCCalPad.h"
51 #include "AliMathBase.h"
52 #include "TTreeStream.h"
55 // TPC calibration algorithm includes
57 #include "AliTPCCalibPedestal.h"
60 Main routine, TPC pedestal detector algorithm to be run on TPC LDC
61 Arguments: list of DATE raw data files
64 int main(int argc, char **argv) {
66 // Main for TPC pedestal detector algorithm
69 Bool_t timeAnalysis = kTRUE;
72 printf("Wrong number of arguments\n");
77 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
84 /* log start of process */
85 printf("TPC DA started - %s\n",__FILE__);
87 /* declare monitoring program */
88 status=monitorDeclareMp( __FILE__ );
90 printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
94 /* copy locally the mapping file from daq detector config db */
95 status = daqDA_DB_getFile(MAPPING_FILE,"./tpcMapping.root");
97 printf("Failed to get mapping file (%s) from DAQdetDB, status=%d\n", MAPPING_FILE, status);
98 printf("Continue anyway ... maybe it works?\n"); // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
99 //return -1; // temporarily uncommented for testing on pcald47 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
102 /* open the mapping file and retrieve mapping object */
103 AliTPCmapper *mapping = 0; // The TPC mapping
104 TFile *fileMapping = new TFile(MAPPING_FILE, "read");
105 mapping = (AliTPCmapper*) fileMapping->Get("tpcMapping");
107 printf("Failed to get mapping object from %s. Exiting ...\n", MAPPING_FILE);
111 printf("Got mapping object from %s\n", MAPPING_FILE);
114 AliTPCCalibPedestal calibPedestal; // pedestal and noise calibration
115 calibPedestal.SetTimeAnalysis(timeAnalysis); // pedestal(t) calibration
116 calibPedestal.SetAltroMapping(mapping->GetAltroMapping()); // Use altro mapping we got from daqDetDb
118 /* loop over RAW data files */
120 for ( i=1; i<argc; i++ ) {
122 /* define data source : this is argument i */
123 printf("Processing file %s\n", argv[i]);
124 status=monitorSetDataSource( argv[i] );
126 printf("monitorSetDataSource() failed. Error=%s. Exiting ...\n", monitorDecodeError(status));
133 struct eventHeaderStruct *event;
135 /* check shutdown condition */
136 if (daqDA_checkShutdown()) {break;}
138 /* get next event (blocking call until timeout) */
139 status=monitorGetEventDynamic((void **)&event);
140 if (status==MON_ERR_EOF) {
141 printf ("End of File %d (%s) detected\n", i, argv[i]);
142 break; /* end of monitoring file has been reached */
145 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
149 /* skip start/end of run events */
150 if ( (event->eventType != physicsEvent) && (event->eventType != calibrationEvent) )
153 /* retry if got no event */
159 // Pedestal calibration
160 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
161 //calibPedestal.ProcessEvent(rawReader);
162 calibPedestal.ProcessEventFast(rawReader); // fast data reader
171 // Analyse pedestals and write them to rootfile
174 calibPedestal.Analyse();
175 calibPedestal.AnalyseTime(nevents);
176 printf ("%d physics/calibration events processed.\n",nevents);
178 TFile *fileTPC = new TFile(RESULT_FILE, "recreate");
179 calibPedestal.Write("calibPedestal");
181 printf("Wrote %s.\n",RESULT_FILE);
184 // Now prepare ASCII files for local ALTRO configuration through DDL.
191 sprintf(filename,"tpcPedestals.data");
192 pedfile.open(filename);
193 sprintf(filename,"tpcNoise.data");
194 noisefile.open(filename);
195 sprintf(filename,"tpcPedestalMem.data");
196 pedmemfile.open(filename);
198 TArrayF **timePed = calibPedestal.GetTimePedestals(); // pedestal values for each time bin
200 Int_t ctr_channel = 0;
202 Int_t ctr_pattern = 0;
204 pedfile << 10 << std::endl; // mark file to contain PEDESTALS per channel
205 noisefile << 11 << std::endl; // mark file to contain NOISE per altro
206 pedmemfile << 12 << std::endl; // mark file to contain PEDESTALs per time bin
208 for ( Int_t roc = 0; roc < 72; roc++ ) {
209 if ( !calibPedestal.GetCalRocPedestal(roc) ) continue;
210 Int_t side = mapping->GetSideFromRoc(roc);
211 Int_t sector = mapping->GetSectorFromRoc(roc);
212 //printf("Analysing ROC %d (side %d, sector %d) ...\n", roc, side, sector);
213 Int_t nru = mapping->IsIROC(roc) ? 2 : 4;
214 for ( int rcu = 0; rcu < nru; rcu++ ) {
215 Int_t patch = mapping->IsIROC(roc) ? rcu : rcu+2;
216 for ( int branch = 0; branch < 2; branch++ ) {
217 for ( int fec = 0; fec < mapping->GetNfec(patch, branch); fec++ ) {
218 for ( int altro = 0; altro < 8; altro++ ) {
221 for ( int channel = 0; channel < 16; channel++ ) {
222 Int_t hwadd = mapping->CodeHWAddress(branch, fec, altro, channel);
223 Int_t row = mapping->GetPadRow(patch, hwadd); // row in a ROC
224 Int_t globalrow = mapping->GetGlobalPadRow(patch, hwadd); // row in full sector
225 Int_t pad = mapping->GetPad(patch, hwadd);
226 Float_t ped = calibPedestal.GetCalRocPedestal(roc)->GetValue(row,pad);
228 if ( ped > 1.e-10 ) {
229 pedfile << ctr_channel << "\t" << side << "\t" << sector << "\t" << patch << "\t"
230 << hwadd << "\t" << ped << std::endl;
234 if ( timePed && fabs(timePed[globalrow][pad].GetSum()) > 1e-10 ) {
235 pedmemfile << ctr_pattern << "\t" << side << "\t" << sector << "\t" << patch
237 for ( Int_t timebin = 0; timebin < 1024; timebin++ )
238 pedmemfile << "\t" << timePed[globalrow][pad].At(timebin);
239 pedmemfile << std::endl;
243 Float_t rms2 = calibPedestal.GetCalRocRMS(roc)->GetValue(row,pad);
244 if ( rms2 > 1.e-10 ) { rms += rms2; ctr += 1.; }
245 } // end channel for loop
246 // noise data (rms) averaged over all channels in this ALTRO.
247 Int_t hwadd = mapping->CodeHWAddress(branch, fec, altro, 0);
248 if ( ctr > 1.e-10 ) {
249 noisefile << ctr_altro << "\t" << side << "\t" << sector << "\t" << patch << "\t"
250 << hwadd << "\t" << rms/ctr << std::endl;
253 } // end altro for loop
254 } // end fec for loop
255 } // end branch for loop
256 } // end rcu for loop
262 printf("Wrote ASCII files.\n");