]> git.uio.no Git - u/mrichter/AliRoot.git/blame - FMD/AliFMDPedestalDA.cxx
Compilation warnings
[u/mrichter/AliRoot.git] / FMD / AliFMDPedestalDA.cxx
CommitLineData
3bd993ba 1/**************************************************************************
2 * Copyright(c) 2004, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16/** @file AliFMDPedestalDA.cxx
17 @author Hans Hjersing Dalsgaard <canute@nbi.dk>
18 @date Mon Mar 10 09:46:05 2008
19 @brief Derived class for the pedestal detector algorithm.
20*/
21//
22// This class implements the virtual functions of the AliFMDBaseDA class.
23// The most important of these functions, FillChannels(..) and Analyse(...) collect
24// and analyse the data of each channel. The resulting pedestal and noise values are
25// written to a comma separated values (csv) file on the go. The csv files produced
26// in this way are the basic input to the AliFMDPreprocessor.
27//
28
29#include "AliFMDPedestalDA.h"
30#include "iostream"
31#include "fstream"
32#include "AliLog.h"
33#include "TF1.h"
34
35//_____________________________________________________________________
36ClassImp(AliFMDPedestalDA)
37
38//_____________________________________________________________________
39AliFMDPedestalDA::AliFMDPedestalDA() : AliFMDBaseDA()
40{
41 fOutputFile.open("peds.csv");
42
43}
44
45//_____________________________________________________________________
46AliFMDPedestalDA::AliFMDPedestalDA(const AliFMDPedestalDA & pedDA) :
47 AliFMDBaseDA(pedDA)
48{
49
50}
51
52//_____________________________________________________________________
53AliFMDPedestalDA::~AliFMDPedestalDA() {
54
55}
56
57//_____________________________________________________________________
58void AliFMDPedestalDA::Init() {
59
60 SetRequiredEvents(1000);
61}
62
63//_____________________________________________________________________
64void AliFMDPedestalDA::AddChannelContainer(TObjArray* sectorArray,
65 UShort_t det,
66 Char_t ring,
67 UShort_t sec,
68 UShort_t strip) {
69
70 TH1S* hChannel = new TH1S(Form("hFMD%d%c_%d_%d",det,ring,sec,strip),
71 Form("hFMD%d%c_%d_%d",det,ring,sec,strip),
72 1024,0,1023);
73
74 hChannel->SetDirectory(0);
75 sectorArray->AddAtAndExpand(hChannel,strip);
76}
77
78//_____________________________________________________________________
79void AliFMDPedestalDA::FillChannels(AliFMDDigit* digit) {
80 UShort_t det = digit->Detector();
81 Char_t ring = digit->Ring();
82 UShort_t sec = digit->Sector();
83 UShort_t strip = digit->Strip();
84
85 TH1S* hChannel = GetChannel(det, ring, sec, strip);
86 hChannel->Fill(digit->Counts());
87
88}
89
90//_____________________________________________________________________
91void AliFMDPedestalDA::Analyse(UShort_t det,
92 Char_t ring,
93 UShort_t sec,
94 UShort_t strip) {
95
96 TH1S* hChannel = GetChannel(det, ring, sec, strip);
80fdb9f3 97 if(hChannel->GetEntries() == 0) {
98 // AliWarning(Form("No entries for FMD%d%c, sector %d, strip %d",det,ring,sec,strip));
3bd993ba 99 return;
100 }
80fdb9f3 101
102 //std::cout<<"Fitting Channel #"<<Form("FMD%d%c_%d_%d",det,ring,sec,strip)<<" with "<<hChannel->GetEntries()<<" entries \r"<<std::flush;
3bd993ba 103
104 Float_t mean = hChannel->GetMean();
105 Float_t rms = hChannel->GetRMS();
106
80fdb9f3 107 hChannel->GetXaxis()->SetRangeUser(mean-4*rms,mean+4*rms);
3bd993ba 108
109 mean = hChannel->GetMean();
110 rms = hChannel->GetRMS();
111
112 hChannel->Fit("gaus","QO","QO",mean-5*rms,mean+5*rms);
113 TF1* fitFunc = hChannel->GetFunction("gaus");
114
115 UInt_t ddl, board, chip, channel;
116
117 UShort_t relStrip = strip%128;
118 AliFMDParameters* pars = AliFMDParameters::Instance();
119 pars->Detector2Hardware(det,ring,sec,strip,ddl,board,chip,channel);
120 Float_t chi2ndf = 0;
121 if(fitFunc->GetNDF())
122 chi2ndf = fitFunc->GetChisquare() / fitFunc->GetNDF();
123 ddl = ddl + kBaseDDL;
124
125 fOutputFile << ddl << ','
126 << board << ','
127 << chip << ','
128 << channel << ','
129 << relStrip << ','
130 << 0 << ','
131 << 0 << ','
132 << mean << ','
133 << rms << ','
134 << fitFunc->GetParameter(1) << ','
135 << fitFunc->GetParameter(2) << ','
136 << chi2ndf <<"\n";
137
138 if(fSaveHistograms) {
139 gDirectory->cd(Form("%s:FMD%d%c/sector_%d/strip_%d",fDiagnosticsFilename,det,ring,sec,strip));
80fdb9f3 140 hChannel->GetXaxis()->SetRange(0,1023);
3bd993ba 141 hChannel->Write();
142
143 }
144
145}
146
147//_____________________________________________________________________
148void AliFMDPedestalDA::WriteHeaderToFile() {
80fdb9f3 149 AliFMDParameters* pars = AliFMDParameters::Instance();
150 fOutputFile.write(Form("# %s \n",pars->GetPedestalShuttleID()),13);
3bd993ba 151 fOutputFile.write("# Rcu, Board, Chip, Channel, Strip, Sample, TimeBin, Pedestal, Noise, Mu, Sigma, Chi2/NDF \n",91);
152
153}
154
155//_____________________________________________________________________
156TH1S* AliFMDPedestalDA::GetChannel(UShort_t det, Char_t ring, UShort_t sec, UShort_t strip) {
157
158 UShort_t Ring = 1;
159 if(ring == 'O')
160 Ring = 0;
161
162
163 TObjArray* detArray = static_cast<TObjArray*>(fDetectorArray.At(det));
164 TObjArray* ringArray = static_cast<TObjArray*>(detArray->At(Ring));
165 TObjArray* secArray = static_cast<TObjArray*>(ringArray->At(sec));
166 TH1S* hChannel = static_cast<TH1S*>(secArray->At(strip));
167
168 return hChannel;
169}
170//_____________________________________________________________________
171//
172//EOF
173//