2 TPC DA for online calibration
4 Contact: Haavard.Helstrup@cern.ch
8 Number of events needed: 100
10 Output Files: tpcPedestal.root, to be exported to the DAQ FXS
11 Trigger types used: CALIBRATION_EVENT
17 TPCda_pedestal.cxx - calibration algorithm for TPC pedestal runs
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
26 contact: marian.ivanov@cern.ch
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).
33 #define RESULT_FILE "tpcPedestal.root"
34 #define MAPPING_FILE "tpcMapping.root"
53 #include "TPluginManager.h"
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"
69 // TPC calibration algorithm includes
71 #include "AliTPCCalibPedestal.h"
74 Main routine, TPC pedestal detector algorithm to be run on TPC LDC
75 Arguments: list of DATE raw data files
78 int main(int argc, char **argv) {
80 // Main for TPC pedestal detector algorithm
83 Bool_t timeAnalysis = kTRUE;
86 printf("Wrong number of arguments\n");
91 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
98 /* log start of process */
99 printf("TPC DA started - %s\n",__FILE__);
101 /* declare monitoring program */
102 status=monitorDeclareMp( __FILE__ );
104 printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
108 /* copy locally the mapping file from daq detector config db */
109 status = daqDA_DB_getFile(MAPPING_FILE,"./tpcMapping.root");
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 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
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");
121 printf("Failed to get mapping object from %s. Exiting ...\n", MAPPING_FILE);
125 printf("Got mapping object from %s\n", MAPPING_FILE);
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
133 /* loop over RAW data files */
135 for ( i=1; i<argc; i++ ) {
137 /* define data source : this is argument i */
138 printf("Processing file %s\n", argv[i]);
139 status=monitorSetDataSource( argv[i] );
141 printf("monitorSetDataSource() failed. Error=%s. Exiting ...\n", monitorDecodeError(status));
148 struct eventHeaderStruct *event;
150 /* check shutdown condition */
151 if (daqDA_checkShutdown()) {break;}
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 */
160 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
164 /* skip start/end of run events */
165 if ( (event->eventType != physicsEvent) && (event->eventType != calibrationEvent) )
168 /* retry if got no event */
174 // Pedestal calibration
175 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
176 //calibPedestal.ProcessEvent(rawReader);
177 calibPedestal.ProcessEventFast(rawReader); // fast data reader
186 // Analyse pedestals and write them to rootfile
189 calibPedestal.Analyse();
190 calibPedestal.AnalyseTime(nevents);
191 printf ("%d physics/calibration events processed.\n",nevents);
193 TFile *fileTPC = new TFile(RESULT_FILE, "recreate");
194 calibPedestal.Write("calibPedestal");
196 printf("Wrote %s.\n",RESULT_FILE);
198 /* store the result file on FES */
200 status=daqDA_FES_storeFile(RESULT_FILE,RESULT_FILE);
207 // Now prepare ASCII files for local ALTRO configuration through DDL.
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);
221 TArrayF **timePed = calibPedestal.GetTimePedestals(); // pedestal values for each time bin
223 Int_t ctr_channel = 0;
225 Int_t ctr_pattern = 0;
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
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++ ) {
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);
251 if ( ped > 1.e-10 ) {
252 pedfile << ctr_channel << "\t" << side << "\t" << sector << "\t" << patch << "\t"
253 << hwadd << "\t" << ped << std::endl;
257 if ( timePed && fabs(timePed[globalrow][pad].GetSum()) > 1e-10 ) {
258 pedmemfile << ctr_pattern << "\t" << side << "\t" << sector << "\t" << patch
260 for ( Int_t timebin = 0; timebin < 1024; timebin++ )
261 pedmemfile << "\t" << timePed[globalrow][pad].At(timebin);
262 pedmemfile << std::endl;
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;
276 } // end altro for loop
277 } // end fec for loop
278 } // end branch for loop
279 } // end rcu for loop
285 printf("Wrote ASCII files.\n");