time bin limits introduced
[u/mrichter/AliRoot.git] / TPC / TPCPEDESTALda.cxx
CommitLineData
c12208b8 1/*
9c754ce7 2TPC DA for online calibration
3
4Contact: Haavard.Helstrup@cern.ch
5Link:
6Run Type: PEDESTAL_RUN
7DA Type: LDC
8Number of events needed: 100
9Input Files:
10Output Files: tpcPedestal.root, to be exported to the DAQ FXS
11Trigger types used: CALIBRATION_EVENT
12
13*/
14
15
16/*
c12208b8 17
18TPCda_pedestal.cxx - calibration algorithm for TPC pedestal runs
19
2010/06/2007 sylvain.chapeland@cern.ch : first version - clean skeleton based on DAQ DA case1
21
22contact: marian.ivanov@cern.ch
23
9c754ce7 24
c12208b8 25This process reads RAW data from the files provided as command line arguments
26and save results in a file (named from RESULT_FILE define - see below).
27
28*/
29
9c754ce7 30#define RESULT_FILE "tpcPedestal.root"
c12208b8 31
32
33extern "C" {
34#include <daqDA.h>
35}
36#include "event.h"
37#include "monitor.h"
9c754ce7 38#include <stdio.h>
39#include <stdlib.h>
40#include <fstream.h>
c12208b8 41
42//
43//Root includes
44//
9c754ce7 45#include <TFile.h>
46#include <TArrayF.h>
c12208b8 47
48//
49//AliRoot includes
50//
51#include "AliRawReader.h"
52#include "AliRawReaderDate.h"
87d57858 53#include "AliTPCmapper.h"
c12208b8 54#include "AliTPCRawStream.h"
55#include "AliTPCROC.h"
56#include "AliTPCCalROC.h"
57#include "AliTPCCalPad.h"
58#include "AliMathBase.h"
59#include "TTreeStream.h"
60
61//
62// TPC calibration algorithm includes
63//
64#include "AliTPCCalibPedestal.h"
65
bdf99a93 66/*
67 Main routine, TPC pedestal detector algorithm to be run on TPC LDC
68 Arguments: list of DATE raw data files
c12208b8 69*/
bdf99a93 70
c12208b8 71int main(int argc, char **argv) {
bdf99a93 72 //
73 // Main for TPC pedestal detector algorithm
74 //
75 Bool_t timeAnalysis = kTRUE;
c12208b8 76
9c754ce7 77 int i,status;
78 AliTPCCalibPedestal calibPedestal; // pedestal and noise calibration
ca273e23 79 calibPedestal.SetRangeTime(60,940); // set time bin range
9c754ce7 80 calibPedestal.SetTimeAnalysis(timeAnalysis); // pedestal(t) calibration
81
c12208b8 82 if (argc<2) {
83 printf("Wrong number of arguments\n");
84 return -1;
85 }
86
c12208b8 87 /* log start of process */
88 printf("TPC DA started - %s\n",__FILE__);
89
c12208b8 90 /* declare monitoring program */
91 status=monitorDeclareMp( __FILE__ );
92 if (status!=0) {
93 printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
94 return -1;
95 }
96
c12208b8 97 /* loop over RAW data files */
98 int nevents=0;
bdf99a93 99 for ( i=1; i<argc; i++ ) {
c12208b8 100
101 /* define data source : this is argument i */
102 printf("Processing file %s\n", argv[i]);
103 status=monitorSetDataSource( argv[i] );
104 if (status!=0) {
9c754ce7 105 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
c12208b8 106 return -1;
107 }
108
109 /* read until EOF */
bdf99a93 110 for ( ; ; ) {
c12208b8 111 struct eventHeaderStruct *event;
112
113 /* check shutdown condition */
114 if (daqDA_checkShutdown()) {break;}
115
116 /* get next event (blocking call until timeout) */
117 status=monitorGetEventDynamic((void **)&event);
118 if (status==MON_ERR_EOF) {
bdf99a93 119 printf ("End of File %d (%s) detected\n", i, argv[i]);
120 break; /* end of monitoring file has been reached */
c12208b8 121 }
c12208b8 122 if (status!=0) {
bdf99a93 123 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
124 break;
c12208b8 125 }
126
bdf99a93 127 /* skip start/end of run events */
128 if ( (event->eventType != physicsEvent) && (event->eventType != calibrationEvent) )
129 continue;
130
c12208b8 131 /* retry if got no event */
9c754ce7 132 if (event==NULL) {
bdf99a93 133 continue;
9c754ce7 134 }
bdf99a93 135
c12208b8 136 nevents++;
137
138 // Pedestal calibration
139 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
9c754ce7 140 calibPedestal.ProcessEvent(rawReader);
c12208b8 141 delete rawReader;
142
143 /* free resources */
144 free(event);
145 }
146 }
147
9c754ce7 148 calibPedestal.Analyse();
149 calibPedestal.AnalyseTime(nevents);
150 printf ("%d physics/calibration events processed\n",nevents);
c12208b8 151
bdf99a93 152 TFile *fileTPC = new TFile(RESULT_FILE, "recreate");
c12208b8 153 calibPedestal.Write("calibPedestal");
154 delete fileTPC;
9c754ce7 155 printf("Wrote %s\n",RESULT_FILE);
c12208b8 156
bdf99a93 157 //
9c754ce7 158 // Now prepare ASCII files for local ALTRO configuration through DDL
bdf99a93 159 //
87d57858 160
bdf99a93 161 ofstream pedfile;
162 ofstream noisefile;
163 ofstream pedmemfile;
87d57858 164 char filename[255];
bdf99a93 165 sprintf(filename,"tpcPedestals.data");
166 pedfile.open(filename);
167 sprintf(filename,"tpcNoise.data");
168 noisefile.open(filename);
169 sprintf(filename,"tpcPedestalMem.data");
170 pedmemfile.open(filename);
171
9c754ce7 172 TArrayF **timePed = calibPedestal.GetTimePedestals();
173 AliTPCmapper mapping;
87d57858 174
bdf99a93 175 Int_t ctr_channel = 0;
176 Int_t ctr_altro = 0;
177 Int_t ctr_pattern = 0;
87d57858 178
9c754ce7 179 pedfile << 10 << std::endl; // PEDESTALS per channel
180 noisefile << 11 << std::endl; // NOISE per altro
181 pedmemfile << 12 << std::endl; // PEDESTALs per time bin
87d57858 182
bdf99a93 183 for ( Int_t roc = 0; roc < 72; roc++ ) {
184 if ( !calibPedestal.GetCalRocPedestal(roc) ) continue;
9c754ce7 185 Int_t side = mapping.GetSideFromRoc(roc);
186 Int_t sector = mapping.GetSectorFromRoc(roc);
bdf99a93 187 //printf("Analysing ROC %d (side %d, sector %d) ...\n", roc, side, sector);
9c754ce7 188 Int_t nru = mapping.IsIROC(roc) ? 2 : 4;
bdf99a93 189 for ( int rcu = 0; rcu < nru; rcu++ ) {
9c754ce7 190 Int_t patch = mapping.IsIROC(roc) ? rcu : rcu+2;
bdf99a93 191 for ( int branch = 0; branch < 2; branch++ ) {
9c754ce7 192 for ( int fec = 0; fec < mapping.GetNfec(patch, branch); fec++ ) {
bdf99a93 193 for ( int altro = 0; altro < 8; altro++ ) {
194 Float_t rms = 0.;
195 Float_t ctr = 0.;
196 for ( int channel = 0; channel < 16; channel++ ) {
9c754ce7 197 Int_t hwadd = mapping.CodeHWAddress(branch, fec, altro, channel);
198 Int_t row = mapping.GetPadRow(patch, hwadd);
199 Int_t pad = mapping.GetPad(patch, hwadd);
200 Float_t ped = calibPedestal.GetCalRocPedestal(roc)->GetValue(row,pad);
bdf99a93 201 // fixed pedestal
202 if ( ped > 1.e-10 ) {
9c754ce7 203 pedfile << ctr_channel << "\t" << side << "\t" << sector << "\t" << patch << "\t" << hwadd << "\t" << ped << std::endl;
bdf99a93 204 ctr_channel++;
205 }
206 // pedestal(t)
9c754ce7 207 if ( timePed && fabs(timePed[row][pad].GetSum()) > 1e-10 ) {
208 pedmemfile << ctr_pattern << "\t" << side << "\t" << sector << "\t" << patch << "\t" << hwadd;
bdf99a93 209 for ( Int_t timebin = 0; timebin < 1024; timebin++ )
9c754ce7 210 pedmemfile << "\t" << timePed[row][pad].At(timebin);
bdf99a93 211 pedmemfile << std::endl;
212 ctr_pattern++;
213 }
214 // rms=noise
215 Float_t rms2 = calibPedestal.GetCalRocRMS(roc)->GetValue(row,pad);
216 if ( rms2 > 1.e-10 ) { rms += rms2; ctr += 1.; }
217 } // end channel for loop
218 // noise data (rms) averaged over all channels in this ALTRO.
9c754ce7 219 Int_t hwadd = mapping.CodeHWAddress(branch, fec, altro, 0);
bdf99a93 220 if ( ctr > 1.e-10 ) {
9c754ce7 221 noisefile << ctr_altro << "\t" << side << "\t" << sector << "\t" << patch << "\t" << hwadd << "\t" << rms/ctr << std::endl;
bdf99a93 222 ctr_altro++;
223 }
224 } // end altro for loop
225 } // end fec for loop
226 } // end branch for loop
227 } // end rcu for loop
228 } // end roc loop
229
230 pedfile.close();
231 noisefile.close();
232 pedmemfile.close();
233
9c754ce7 234 printf("Wrote ASCII files\n");
87d57858 235
c12208b8 236 return status;
237}