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 
156  //Need our own copy for thread safety
157  TF1 mygaus("mygaus","gaus");
158  for ( uint16_t istr = 0; istr < 128; istr++ ) {
159 
160  ana->ksProb_[iapv].push_back(0);
161  ana->noiseGaus_[iapv].push_back(0);
162  ana->noiseBin84_[iapv].push_back(0);
163  ana->noiseRMS_[iapv].push_back(0);
164  ana->noiseSignif_[iapv].push_back(0);
165 
166  // pedestals and raw noise
167  if ( peds_histo ) {
168  if ( peds_histo->GetBinEntries(iapv*128 + istr + 1) ) {
169  ana->peds_[iapv][istr] = peds_histo->GetBinContent(iapv*128 + istr + 1);
170  p_sum += ana->peds_[iapv][istr];
171  p_sum2 += (ana->peds_[iapv][istr] * ana->peds_[iapv][istr]);
172  if ( ana->peds_[iapv][istr] > p_max ) { p_max = ana->peds_[iapv][istr];}
173  if ( ana->peds_[iapv][istr] < p_min ) { p_min = ana->peds_[iapv][istr];}
174  ana->raw_[iapv][istr] = peds_histo->GetBinError(iapv*128 + istr + 1);
175  r_sum += ana->raw_[iapv][istr];
176  r_sum2 += (ana->raw_[iapv][istr] * ana->raw_[iapv][istr]);
177  if ( ana->raw_[iapv][istr] > r_max ) { r_max = ana->raw_[iapv][istr]; }
178  if ( ana->raw_[iapv][istr] < r_min ) { r_min = ana->raw_[iapv][istr]; }
179  }
180  }
181  // Noise from Full Distribution
182  if ( noise_histo ) {
183  // Fit the ADC Distribution from TH2S by projecting it out and fitting.
184 
185  TH1S * noisehist = new TH1S("noisehist","",noise_histo->GetNbinsX(),
186  -noise_histo->GetNbinsX()/2,noise_histo->GetNbinsX()/2);
187 
188 
189  for(int i=0;i<=noise_histo->GetNbinsX()+1;i++){
190  noisehist->SetBinContent(i,noise_histo->GetBinContent(i,iapv*128 + istr + 1));
191  }
192  // If the histogram is valid.
193  if(noisehist->Integral() > 0){
194  ana->noiseRMS_[iapv][istr] = noisehist->GetRMS();
195  noisehist->Fit(&mygaus,"Q");
196  ana->noiseGaus_[iapv][istr] = mygaus.GetParameter(2);
197 
198  // new Bin84 method
199  std::vector<float> integralFrac;
200  integralFrac.push_back(1.*noisehist->GetBinContent(0)/noisehist->Integral(0,noisehist->GetNbinsX()));
201  // Calculate the integral of distro as a function of bins.
202  for(int i = 1; i < noisehist->GetNbinsX();i++){
203  integralFrac.push_back(float(noisehist->GetBinContent(i))/
204  noisehist->Integral(0,noisehist->GetNbinsX())+integralFrac[i-1]);
205  //Take the two bins next to 84% and solve for x in 0.84 = mx+b
206  if (integralFrac[i] >= 0.84135 && integralFrac[i-1] < 0.84135) {
207  // my quadratic noise calculation
208  double w = noisehist->GetBinWidth(i);
209  double a = noisehist->GetBinContent(i-1);
210  double b = noisehist->GetBinContent(i);
211  double f = w*(0.84135 -integralFrac[i-1])/(integralFrac[i]-integralFrac[i-1]);
212  double x = 0;
213  if (a==b) {
214  x = f;
215  } else {
216  double aa = (b-a)/(2*w);
217  double bb = (b+a)/2;
218  double cc = -b*f;
219  double dd = bb*bb-4*aa*cc; //if (dd<0) dd=0;
220  x = (-bb+sqrt(dd))/(2*aa);
221  }
222  ana->noiseBin84_[iapv][istr] = noisehist->GetBinLowEdge(i) + x;
223  }
224  }
225  // Compare shape of ADC to a histogram made of Gaus Fit for KSProb, Chi2Prob, Etc...
226  TH1D * FitHisto = new TH1D("FitHisto","FitHisto",noisehist->GetNbinsX(),
227  -noisehist->GetNbinsX()/2,noisehist->GetNbinsX()/2);
228  FitHisto->Add(&mygaus);
229  FitHisto->Sumw2();
230  noisehist->Sumw2();
231 
232  if(FitHisto->Integral() > 0){
233  // This is stored as a float but will be plotted as an int at the summary histos.
234  // This forces the histo to draw 10000 bins instead of 1.
235  ana->ksProb_[iapv][istr] = noisehist->KolmogorovTest(FitHisto)*10000;
236  }
237  delete FitHisto;
238  }
239  delete noisehist;
240  }
241  // Assigning the actual noise values used for Upload!!!!!!!!!!!!!!!!!!!!
242  if (noiseDef_ == "Bin84") {
243  if (ana->noiseBin84_[iapv][istr] > 0) {
244  ana->noise_[iapv][istr] = ana->noiseBin84_[iapv][istr];
245  } else {
246  ana->noise_[iapv][istr] = ana->noiseRMS_[iapv][istr];
247  }
248  } else if (noiseDef_ == "RMS") {
249  ana->noise_[iapv][istr] = ana->noiseRMS_[iapv][istr];
250  } else edm::LogWarning(mlCommissioning_)<< "[PedsFullNoiseAlgorithm::"
251  << __func__ << "]"<< " Unknown noise definition!!!";
252 
253  // Setting Sum of RMS and RMS^2 for Dead/Noisy Strip calculations
254  n_sum += ana->noise_[iapv][istr];
255  n_sum2 += (ana->noise_[iapv][istr] * ana->noise_[iapv][istr]);
256  if ( ana->noise_[iapv][istr] > n_max ) { n_max = ana->noise_[iapv][istr]; }
257  if ( ana->noise_[iapv][istr] < n_min ) { n_min = ana->noise_[iapv][istr]; }
258 
259  } // strip loop
260 
261  // Calc mean and rms for peds
262  if ( !ana->peds_[iapv].empty() ) {
263  p_sum /= static_cast<float>( ana->peds_[iapv].size() );
264  p_sum2 /= static_cast<float>( ana->peds_[iapv].size() );
265  ana->pedsMean_[iapv] = p_sum;
266  ana->pedsSpread_[iapv] = sqrt( fabs(p_sum2 - p_sum*p_sum) );
267  }
268 
269  // Calc mean and rms for noise using noiseRMS.
270  if ( !ana->noise_[iapv].empty() ) {
271  n_sum /= static_cast<float>( ana->noise_[iapv].size() );
272  n_sum2 /= static_cast<float>( ana->noise_[iapv].size() );
273  ana->noiseMean_[iapv] = n_sum;
274  ana->noiseSpread_[iapv] = sqrt( fabs(n_sum2 - n_sum*n_sum) );
275  }
276 
277  // Calc mean and rms for raw noise
278  if ( !ana->raw_[iapv].empty() ) {
279  r_sum /= static_cast<float>( ana->raw_[iapv].size() );
280  r_sum2 /= static_cast<float>( ana->raw_[iapv].size() );
281  ana->rawMean_[iapv] = r_sum;
282  ana->rawSpread_[iapv] = sqrt( fabs(r_sum2 - r_sum*r_sum) );
283  }
284 
285  // Set max and min values for peds, noise and raw noise
286  if ( p_max > -1.*sistrip::maximum_ ) { ana->pedsMax_[iapv] = p_max; }
287  if ( p_min < 1.*sistrip::maximum_ ) { ana->pedsMin_[iapv] = p_min; }
288  if ( n_max > -1.*sistrip::maximum_ ) { ana->noiseMax_[iapv] = n_max; }
289  if ( n_min < 1.*sistrip::maximum_ ) { ana->noiseMin_[iapv] = n_min; }
290  if ( r_max > -1.*sistrip::maximum_ ) { ana->rawMax_[iapv] = r_max; }
291  if ( r_min < 1.*sistrip::maximum_ ) { ana->rawMin_[iapv] = r_min; }
292 
293  // Set dead and noisy strips
294  for ( uint16_t istr = 0; istr < 128; istr++ ) { // strip loop
295  // Set the significance of the noise of each strip also compared to apv avg.
296  ana->noiseSignif_[iapv][istr] = (ana->noise_[iapv][istr]-ana->noiseMean_[iapv])/ana->noiseSpread_[iapv];
297  if ( ana->noiseMin_[iapv] > sistrip::maximum_ || ana->noiseMax_[iapv] > sistrip::maximum_ ) {
298  continue;
299  }
300  // Strip Masking for Dead Strips
301  if(ana->noiseSignif_[iapv][istr] < -deadStripMax_){
302  ana->dead_[iapv].push_back(istr);
303  SiStripFecKey fec_key(ana->fecKey());
304  LogTrace(mlDqmClient_)<<"DeadSignif "<<ana->noiseSignif_[iapv][istr]
305  <<" "<<fec_key.fecCrate()
306  <<" "<<fec_key.fecSlot()
307  <<" "<<fec_key.fecRing()
308  <<" "<<fec_key.ccuAddr()
309  <<" "<<fec_key.ccuChan()
310  <<" "<<fec_key.lldChan()
311  <<" "<<iapv*128+istr<<std::endl;
312  } // Strip Masking for Dead Strips
313  // Laurent's Method for Flagging bad strips in TIB
314  else if((ana->noiseMax_[iapv]/ana->noiseMean_[iapv] > 3 || ana->noiseSpread_[iapv] > 3)
315  && ana->noiseSignif_[iapv][istr] > 1){
316  ana->noisy_[iapv].push_back(istr);
317  SiStripFecKey fec_key(ana->fecKey());
318  LogTrace(mlDqmClient_)<<"NoisyLM "<<ana->noiseMax_[iapv]/ana->noiseMean_[iapv]
319  <<" "<<fec_key.fecCrate()
320  <<" "<<fec_key.fecSlot()
321  <<" "<<fec_key.fecRing()
322  <<" "<<fec_key.ccuAddr()
323  <<" "<<fec_key.ccuChan()
324  <<" "<<fec_key.lldChan()
325  <<" "<<iapv*128+istr<<std::endl;
326  } // if NoisyLM
327  //Strip Masking for Non Gassian Strips which are also noisy.
328  else if(ana->ksProb_[iapv][istr] < ksProbCut_){
329  ana->noisy_[iapv].push_back(istr);
330  SiStripFecKey fec_key(ana->fecKey());
331  LogTrace(mlDqmClient_)<<"NoisyKS "<<ana->ksProb_[iapv][istr]
332  <<" "<<fec_key.fecCrate()
333  <<" "<<fec_key.fecSlot()
334  <<" "<<fec_key.fecRing()
335  <<" "<<fec_key.ccuAddr()
336  <<" "<<fec_key.ccuChan()
337  <<" "<<fec_key.lldChan()
338  <<" "<<iapv*128+istr<<std::endl;
339  } //Strip Masking for Non Gassian Strips which are also noisy.
340  else if(ana->noiseSignif_[iapv][istr] > noisyStripMin_){
341  ana->noisy_[iapv].push_back(istr);
342  SiStripFecKey fec_key(ana->fecKey());
343  LogTrace(mlDqmClient_)<<"NoisySignif "<<ana->noiseSignif_[iapv][istr]
344  <<" "<<fec_key.fecCrate()
345  <<" "<<fec_key.fecSlot()
346  <<" "<<fec_key.fecRing()
347  <<" "<<fec_key.ccuAddr()
348  <<" "<<fec_key.ccuChan()
349  <<" "<<fec_key.lldChan()
350  <<" "<<iapv*128+istr<<std::endl;
351  } // if Signif
352  }// strip loop to set dead or noisy strips
353  } // apv loop
354  //std::cout << std::endl;
355 }
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.
if(conf.exists("allCellsPositionCalc"))
CommissioningAnalysis *const anal() const
static const char nullPtr_[]