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