CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PedsFullNoiseAlgorithm.cc
Go to the documentation of this file.
6 #include "TProfile.h"
7 #include "TH1.h"
8 #include "TH2.h"
9 #include "TF1.h"
10 #include <iostream>
11 #include <iomanip>
12 #include <cmath>
13 
14 using namespace sistrip;
15 
16 // ----------------------------------------------------------------------------
17 //
19  : CommissioningAlgorithm(anal),
20  hPeds_(0,""),
21  hNoise_(0,""),
22  hNoise1D_(0,""),
23  deadStripMax_(pset.getParameter<double>("DeadStripMax")),
24  noisyStripMin_(pset.getParameter<double>("NoisyStripMin")),
25  noiseDef_(pset.getParameter<std::string>("NoiseDefinition")),
26  ksProbCut_(pset.getParameter<double>("KsProbCut"))
27 {
28  //LogDebug(mlCommissioning_)
29  // << "[PedsFullNoiseAlgorithm::" << __func__ << "]"
30  // << " Set maximum noise deviation for dead strip determination to: " << deadStripMax_;
31  // LogDebug(mlCommissioning_)
32  // << "[PedsFullNoiseAlgorithm::" << __func__ << "]"
33  // << " Set minimal noise deviation for noise strip determination to: " << noisyStripMin_;
34 }
35 
36 // ----------------------------------------------------------------------------
37 //
38 
39 void PedsFullNoiseAlgorithm::extract( const std::vector<TH1*>& histos ) {
40 
41  if ( !anal() ) {
43  << "[PedsFullNoiseAlgorithm::" << __func__ << "]"
44  << " NULL pointer to Analysis object!";
45  return;
46  }
47 
48  // Check number of histograms
49  if ( histos.size() != 3 ) {
51  }
52 
53  // Extract FED key from histo title
54  if ( !histos.empty() ) {
55  anal()->fedKey( extractFedKey( histos.front() ) );
56  }
57 
58  // Extract 1D histograms
59  std::vector<TH1*>::const_iterator ihis = histos.begin();
60  for ( ; ihis != histos.end(); ihis++ ) {
61 
62  // Check for NULL pointer
63  if ( !(*ihis) ) { continue; }
64 
65 // TO BE UPDATED!!!
66  // Check run type
67  SiStripHistoTitle title( (*ihis)->GetName() );
68 /* if ( title.runType() != sistrip::PEDS_FULL_NOISE ) {
69  anal()->addErrorCode(sistrip::unexpectedTask_);
70  continue;
71  }
72 */
73  // Extract peds histos
74  if ( title.extraInfo().find(sistrip::extrainfo::roughPedestals_) != std::string::npos ) {
75  //@@ something here for rough peds?
76  } else if ( title.extraInfo().find(sistrip::extrainfo::pedestals_) != std::string::npos ) {
77  hPeds_.first = *ihis;
78  hPeds_.second = (*ihis)->GetName();
79  } else if ( title.extraInfo().find(sistrip::extrainfo::commonMode_) != std::string::npos ) {
80  //@@ something here for CM plots?
81  } else if ( title.extraInfo().find(sistrip::extrainfo::noiseProfile_) != std::string::npos ) {
82  //@@ something here for noise profile plot?
83  hNoise1D_.first = *ihis;
84  hNoise1D_.second = (*ihis)->GetName();
85  } else if ( title.extraInfo().find(sistrip::extrainfo::noise2D_) != std::string::npos ) {
86  hNoise_.first = *ihis;
87  hNoise_.second = (*ihis)->GetName();
88  } else {
90  }
91  }
92 
93 }
94 
95 // -----------------------------------------------------------------------------
96 //
98 
99  if ( !anal() ) {
101  << "[PedsFullNoiseAlgorithm::" << __func__ << "]"
102  << " NULL pointer to base Analysis object!";
103  return;
104  }
105 
107  PedsFullNoiseAnalysis* ana = dynamic_cast<PedsFullNoiseAnalysis*>( tmp );
108  if ( !ana ) {
110  << "[PedsFullNoiseAlgorithm::" << __func__ << "]"
111  << " NULL pointer to derived Analysis object!";
112  return;
113  }
114 
115  if ( !hPeds_.first ) {
117  return;
118  }
119 
120  if ( !hNoise_.first ) {
122  return;
123  }
124 
125  TProfile * peds_histo = dynamic_cast<TProfile *>(hPeds_.first);
126  TH2S * noise_histo = dynamic_cast<TH2S *>(hNoise_.first);
127  if ( !peds_histo ) {
129  return;
130  }
131 
132  if ( !noise_histo ) {
134  return;
135  }
136 
137  if ( peds_histo->GetNbinsX() != 256 ) {
139  return;
140  }
141 
142  if ( noise_histo->GetNbinsY() != 256 ) { // X range is configurable
144  return;
145  }
146 
147  // Iterate through APVs
148  for ( uint16_t iapv = 0; iapv < 2; iapv++ ) {
149 
150  // Used to calc mean and rms for peds and noise
151  float p_sum = 0., p_sum2 = 0., p_max = -1.*sistrip::invalid_, p_min = sistrip::invalid_;
152  float n_sum = 0., n_sum2 = 0., n_max = -1.*sistrip::invalid_, n_min = sistrip::invalid_;
153  float r_sum = 0., r_sum2 = 0., r_max = -1.*sistrip::invalid_, r_min = sistrip::invalid_;
154  // Iterate through strips of APV
155  for ( uint16_t istr = 0; istr < 128; istr++ ) {
156 
157  ana->ksProb_[iapv].push_back(0);
158  ana->noiseGaus_[iapv].push_back(0);
159  ana->noiseBin84_[iapv].push_back(0);
160  ana->noiseRMS_[iapv].push_back(0);
161  ana->noiseSignif_[iapv].push_back(0);
162 
163  // pedestals and raw noise
164  if ( peds_histo ) {
165  if ( peds_histo->GetBinEntries(iapv*128 + istr + 1) ) {
166  ana->peds_[iapv][istr] = peds_histo->GetBinContent(iapv*128 + istr + 1);
167  p_sum += ana->peds_[iapv][istr];
168  p_sum2 += (ana->peds_[iapv][istr] * ana->peds_[iapv][istr]);
169  if ( ana->peds_[iapv][istr] > p_max ) { p_max = ana->peds_[iapv][istr];}
170  if ( ana->peds_[iapv][istr] < p_min ) { p_min = ana->peds_[iapv][istr];}
171  ana->raw_[iapv][istr] = peds_histo->GetBinError(iapv*128 + istr + 1);
172  r_sum += ana->raw_[iapv][istr];
173  r_sum2 += (ana->raw_[iapv][istr] * ana->raw_[iapv][istr]);
174  if ( ana->raw_[iapv][istr] > r_max ) { r_max = ana->raw_[iapv][istr]; }
175  if ( ana->raw_[iapv][istr] < r_min ) { r_min = ana->raw_[iapv][istr]; }
176  }
177  }
178  // Noise from Full Distribution
179  if ( noise_histo ) {
180  // Fit the ADC Distribution from TH2S by projecting it out and fitting.
181 
182  TH1S * noisehist = new TH1S("noisehist","",noise_histo->GetNbinsX(),
183  -noise_histo->GetNbinsX()/2,noise_histo->GetNbinsX()/2);
184 
185 
186  for(int i=0;i<=noise_histo->GetNbinsX()+1;i++){
187  noisehist->SetBinContent(i,noise_histo->GetBinContent(i,iapv*128 + istr + 1));
188  }
189  // If the histogram is valid.
190  if(noisehist->Integral() > 0){
191  ana->noiseRMS_[iapv][istr] = noisehist->GetRMS();
192  noisehist->Fit("gaus","Q");
193  ana->noiseGaus_[iapv][istr] = noisehist->GetFunction("gaus")->GetParameter(2);
194 
195  // new Bin84 method
196  std::vector<float> integralFrac;
197  integralFrac.push_back(1.*noisehist->GetBinContent(0)/noisehist->Integral(0,noisehist->GetNbinsX()));
198  // Calculate the integral of distro as a function of bins.
199  for(int i = 1; i < noisehist->GetNbinsX();i++){
200  integralFrac.push_back(float(noisehist->GetBinContent(i))/
201  noisehist->Integral(0,noisehist->GetNbinsX())+integralFrac[i-1]);
202  //Take the two bins next to 84% and solve for x in 0.84 = mx+b
203  if (integralFrac[i] >= 0.84135 && integralFrac[i-1] < 0.84135) {
204  // my quadratic noise calculation
205  double w = noisehist->GetBinWidth(i);
206  double a = noisehist->GetBinContent(i-1);
207  double b = noisehist->GetBinContent(i);
208  double f = w*(0.84135 -integralFrac[i-1])/(integralFrac[i]-integralFrac[i-1]);
209  double x = 0;
210  if (a==b) {
211  x = f;
212  } else {
213  double aa = (b-a)/(2*w);
214  double bb = (b+a)/2;
215  double cc = -b*f;
216  double dd = bb*bb-4*aa*cc; //if (dd<0) dd=0;
217  x = (-bb+sqrt(dd))/(2*aa);
218  }
219  ana->noiseBin84_[iapv][istr] = noisehist->GetBinLowEdge(i) + x;
220  }
221  }
222  // Compare shape of ADC to a histogram made of Gaus Fit for KSProb, Chi2Prob, Etc...
223  TH1D * FitHisto = new TH1D("FitHisto","FitHisto",noisehist->GetNbinsX(),
224  -noisehist->GetNbinsX()/2,noisehist->GetNbinsX()/2);
225  FitHisto->Add(noisehist->GetFunction("gaus"));
226  FitHisto->Sumw2();
227  noisehist->Sumw2();
228 
229  if(FitHisto->Integral() > 0){
230  // This is stored as a float but will be plotted as an int at the summary histos.
231  // This forces the histo to draw 10000 bins instead of 1.
232  ana->ksProb_[iapv][istr] = noisehist->KolmogorovTest(FitHisto)*10000;
233  }
234  delete FitHisto;
235  }
236  delete noisehist;
237  }
238  // Assigning the actual noise values used for Upload!!!!!!!!!!!!!!!!!!!!
239  if (noiseDef_ == "Bin84") {
240  if (ana->noiseBin84_[iapv][istr] > 0) {
241  ana->noise_[iapv][istr] = ana->noiseBin84_[iapv][istr];
242  } else {
243  ana->noise_[iapv][istr] = ana->noiseRMS_[iapv][istr];
244  }
245  } else if (noiseDef_ == "RMS") {
246  ana->noise_[iapv][istr] = ana->noiseRMS_[iapv][istr];
247  } else edm::LogWarning(mlCommissioning_)<< "[PedsFullNoiseAlgorithm::"
248  << __func__ << "]"<< " Unknown noise definition!!!";
249 
250  // Setting Sum of RMS and RMS^2 for Dead/Noisy Strip calculations
251  n_sum += ana->noise_[iapv][istr];
252  n_sum2 += (ana->noise_[iapv][istr] * ana->noise_[iapv][istr]);
253  if ( ana->noise_[iapv][istr] > n_max ) { n_max = ana->noise_[iapv][istr]; }
254  if ( ana->noise_[iapv][istr] < n_min ) { n_min = ana->noise_[iapv][istr]; }
255 
256  } // strip loop
257 
258  // Calc mean and rms for peds
259  if ( !ana->peds_[iapv].empty() ) {
260  p_sum /= static_cast<float>( ana->peds_[iapv].size() );
261  p_sum2 /= static_cast<float>( ana->peds_[iapv].size() );
262  ana->pedsMean_[iapv] = p_sum;
263  ana->pedsSpread_[iapv] = sqrt( fabs(p_sum2 - p_sum*p_sum) );
264  }
265 
266  // Calc mean and rms for noise using noiseRMS.
267  if ( !ana->noise_[iapv].empty() ) {
268  n_sum /= static_cast<float>( ana->noise_[iapv].size() );
269  n_sum2 /= static_cast<float>( ana->noise_[iapv].size() );
270  ana->noiseMean_[iapv] = n_sum;
271  ana->noiseSpread_[iapv] = sqrt( fabs(n_sum2 - n_sum*n_sum) );
272  }
273 
274  // Calc mean and rms for raw noise
275  if ( !ana->raw_[iapv].empty() ) {
276  r_sum /= static_cast<float>( ana->raw_[iapv].size() );
277  r_sum2 /= static_cast<float>( ana->raw_[iapv].size() );
278  ana->rawMean_[iapv] = r_sum;
279  ana->rawSpread_[iapv] = sqrt( fabs(r_sum2 - r_sum*r_sum) );
280  }
281 
282  // Set max and min values for peds, noise and raw noise
283  if ( p_max > -1.*sistrip::maximum_ ) { ana->pedsMax_[iapv] = p_max; }
284  if ( p_min < 1.*sistrip::maximum_ ) { ana->pedsMin_[iapv] = p_min; }
285  if ( n_max > -1.*sistrip::maximum_ ) { ana->noiseMax_[iapv] = n_max; }
286  if ( n_min < 1.*sistrip::maximum_ ) { ana->noiseMin_[iapv] = n_min; }
287  if ( r_max > -1.*sistrip::maximum_ ) { ana->rawMax_[iapv] = r_max; }
288  if ( r_min < 1.*sistrip::maximum_ ) { ana->rawMin_[iapv] = r_min; }
289 
290  // Set dead and noisy strips
291  for ( uint16_t istr = 0; istr < 128; istr++ ) { // strip loop
292  // Set the significance of the noise of each strip also compared to apv avg.
293  ana->noiseSignif_[iapv][istr] = (ana->noise_[iapv][istr]-ana->noiseMean_[iapv])/ana->noiseSpread_[iapv];
294  if ( ana->noiseMin_[iapv] > sistrip::maximum_ || ana->noiseMax_[iapv] > sistrip::maximum_ ) {
295  continue;
296  }
297  // Strip Masking for Dead Strips
298  if(ana->noiseSignif_[iapv][istr] < -deadStripMax_){
299  ana->dead_[iapv].push_back(istr);
300  SiStripFecKey fec_key(ana->fecKey());
301  LogTrace(mlDqmClient_)<<"DeadSignif "<<ana->noiseSignif_[iapv][istr]
302  <<" "<<fec_key.fecCrate()
303  <<" "<<fec_key.fecSlot()
304  <<" "<<fec_key.fecRing()
305  <<" "<<fec_key.ccuAddr()
306  <<" "<<fec_key.ccuChan()
307  <<" "<<fec_key.lldChan()
308  <<" "<<iapv*128+istr<<std::endl;
309  } // Strip Masking for Dead Strips
310  // Laurent's Method for Flagging bad strips in TIB
311  else if((ana->noiseMax_[iapv]/ana->noiseMean_[iapv] > 3 || ana->noiseSpread_[iapv] > 3)
312  && ana->noiseSignif_[iapv][istr] > 1){
313  ana->noisy_[iapv].push_back(istr);
314  SiStripFecKey fec_key(ana->fecKey());
315  LogTrace(mlDqmClient_)<<"NoisyLM "<<ana->noiseMax_[iapv]/ana->noiseMean_[iapv]
316  <<" "<<fec_key.fecCrate()
317  <<" "<<fec_key.fecSlot()
318  <<" "<<fec_key.fecRing()
319  <<" "<<fec_key.ccuAddr()
320  <<" "<<fec_key.ccuChan()
321  <<" "<<fec_key.lldChan()
322  <<" "<<iapv*128+istr<<std::endl;
323  } // if NoisyLM
324  //Strip Masking for Non Gassian Strips which are also noisy.
325  else if(ana->ksProb_[iapv][istr] < ksProbCut_){
326  ana->noisy_[iapv].push_back(istr);
327  SiStripFecKey fec_key(ana->fecKey());
328  LogTrace(mlDqmClient_)<<"NoisyKS "<<ana->ksProb_[iapv][istr]
329  <<" "<<fec_key.fecCrate()
330  <<" "<<fec_key.fecSlot()
331  <<" "<<fec_key.fecRing()
332  <<" "<<fec_key.ccuAddr()
333  <<" "<<fec_key.ccuChan()
334  <<" "<<fec_key.lldChan()
335  <<" "<<iapv*128+istr<<std::endl;
336  } //Strip Masking for Non Gassian Strips which are also noisy.
337  else if(ana->noiseSignif_[iapv][istr] > noisyStripMin_){
338  ana->noisy_[iapv].push_back(istr);
339  SiStripFecKey fec_key(ana->fecKey());
340  LogTrace(mlDqmClient_)<<"NoisySignif "<<ana->noiseSignif_[iapv][istr]
341  <<" "<<fec_key.fecCrate()
342  <<" "<<fec_key.fecSlot()
343  <<" "<<fec_key.fecRing()
344  <<" "<<fec_key.ccuAddr()
345  <<" "<<fec_key.ccuChan()
346  <<" "<<fec_key.lldChan()
347  <<" "<<iapv*128+istr<<std::endl;
348  } // if Signif
349  }// strip loop to set dead or noisy strips
350  } // apv loop
351  //std::cout << std::endl;
352 }
int i
Definition: DBlmapReader.cc:9
const uint32_t & fedKey() const
const double w
Definition: UKUtility.cc:23
Utility class that holds histogram title.
static const char numberOfHistos_[]
void extract(const std::vector< TH1 * > &)
static const char mlDqmClient_[]
static const char unexpectedExtraInfo_[]
static const char numberOfBins_[]
Histogram-based analysis for pedestal run.
Utility class that identifies a position within the strip tracker control structure, down to the level of an APV25.
Definition: SiStripFecKey.h:45
static const char noise2D_[]
static const char mlCommissioning_[]
T sqrt(T t)
Definition: SSEVec.h:48
uint32_t extractFedKey(const TH1 *const )
static const uint16_t maximum_
Definition: Constants.h:20
static const char roughPedestals_[]
virtual void addErrorCode(const std::string &error)
const uint32_t & fecKey() const
double f[11][100]
static const char commonMode_[]
static const char noiseProfile_[]
#define LogTrace(id)
static const char pedestals_[]
static const uint16_t invalid_
Definition: Constants.h:16
double b
Definition: hdecay.h:120
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
double a
Definition: hdecay.h:121
Abstract base for derived classes that provide analysis of commissioning histograms.
Definition: DDAxes.h:10
if(conf.exists("allCellsPositionCalc"))
CommissioningAnalysis *const anal() const
static const char nullPtr_[]