CMS 3D CMS Logo

OptoScanAlgorithm.cc
Go to the documentation of this file.
7 #include "TProfile.h"
8 #include "TH1.h"
9 #include <iostream>
10 #include <sstream>
11 #include <iomanip>
12 
13 using namespace sistrip;
14 
15 // ----------------------------------------------------------------------------
16 //
18  : CommissioningAlgorithm(anal),
19  histos_( 4, std::vector<Histo>( 3, Histo(0,"") ) ),
20  targetGain_(pset.getParameter<double>("TargetGain"))
21 {
23  << "[PedestalsAlgorithm::" << __func__ << "]"
24  << " Set target gain to: " << targetGain_;
25 }
26 
27 // ----------------------------------------------------------------------------
28 //
29 void OptoScanAlgorithm::extract( const std::vector<TH1*>& histos ) {
30 
31  if ( !anal() ) {
33  << "[OptoScanAlgorithm::" << __func__ << "]"
34  << " NULL pointer to Analysis object!";
35  return;
36  }
37 
38  // Check number of histograms
39  if ( histos.size() != 12 ) {
41  }
42 
43  // Extract FED key from histo title
44  if ( !histos.empty() ) anal()->fedKey( extractFedKey( histos.front() ) );
45 
46  // Extract histograms
47  std::vector<TH1*>::const_iterator ihis = histos.begin();
48  for ( ; ihis != histos.end(); ihis++ ) {
49 
50  // Check for NULL pointer
51  if ( !(*ihis) ) { continue; }
52 
53  // Check name
54  SiStripHistoTitle title( (*ihis)->GetName() );
55  if ( title.runType() != sistrip::OPTO_SCAN ) {
57  continue;
58  }
59 
60  // Extract gain setting and digital high/low info
61  uint16_t gain = sistrip::invalid_;
62  if ( title.extraInfo().find(sistrip::extrainfo::gain_) != std::string::npos ) {
63  std::stringstream ss;
64  ss << title.extraInfo().substr( title.extraInfo().find(sistrip::extrainfo::gain_) + (sizeof(sistrip::extrainfo::gain_) - 1), 1 );
65  ss >> std::dec >> gain;
66  }
67  uint16_t digital = sistrip::invalid_;
68  if ( title.extraInfo().find(sistrip::extrainfo::digital_) != std::string::npos ) {
69  std::stringstream ss;
70  ss << title.extraInfo().substr( title.extraInfo().find(sistrip::extrainfo::digital_) + (sizeof(sistrip::extrainfo::digital_) - 1), 1 );
71  ss >> std::dec >> digital;
72  }
73  bool baseline_rms = false;
74  if ( title.extraInfo().find(sistrip::extrainfo::baselineRms_) != std::string::npos ) {
75  baseline_rms = true;
76  }
77 
78  if ( gain <= 3 ) {
79  if ( digital <= 1 ) {
80  histos_[gain][digital].first = *ihis;
81  histos_[gain][digital].second = (*ihis)->GetName();
82  } else if ( baseline_rms ) {
83  histos_[gain][2].first = *ihis;
84  histos_[gain][2].second = (*ihis)->GetName();
85  } else {
87  }
88  } else {
90  }
91 
92  }
93 
94 }
95 
96 // ----------------------------------------------------------------------------
97 //
99 
100  if ( !anal() ) {
102  << "[OptoScanAlgorithm::" << __func__ << "]"
103  << " NULL pointer to base Analysis object!";
104  return;
105  }
106 
108  OptoScanAnalysis* anal = dynamic_cast<OptoScanAnalysis*>( tmp );
109  if ( !anal ) {
111  << "[OptoScanAlgorithm::" << __func__ << "]"
112  << " NULL pointer to derived Analysis object!";
113  return;
114  }
115 
116  // Iterate through four gain settings
117  for ( uint16_t igain = 0; igain < 4; igain++ ) {
118 
119  // Select histos appropriate for gain setting
120  TH1* base_his = histos_[igain][0].first;
121  TH1* peak_his = histos_[igain][1].first;
122  TH1* noise_his = histos_[igain][2].first;
123 
124  if ( !base_his ) {
126  return;
127  }
128 
129  if ( !peak_his ) {
131  return;
132  }
133 
134  if ( !noise_his ) {
136  return;
137  }
138 
139  TProfile* base_histo = dynamic_cast<TProfile*>(base_his);
140  if ( !base_histo ) {
142  return;
143  }
144 
145  TProfile* peak_histo = dynamic_cast<TProfile*>(peak_his);
146  if ( !peak_histo ) {
148  return;
149  }
150 
151  TProfile* noise_histo = dynamic_cast<TProfile*>(noise_his);
152  if ( !noise_histo ) {
154  return;
155  }
156 
157  // Check histogram binning
158  uint16_t nbins = static_cast<uint16_t>( peak_histo->GetNbinsX() );
159  if ( static_cast<uint16_t>( base_histo->GetNbinsX() ) != nbins ) {
161  if ( static_cast<uint16_t>( base_histo->GetNbinsX() ) < nbins ) {
162  nbins = static_cast<uint16_t>( base_histo->GetNbinsX() );
163  }
164  }
165 
166  // Some containers
167  std::vector<float> peak_contents(0);
168  std::vector<float> peak_errors(0);
169  std::vector<float> peak_entries(0);
170  std::vector<float> base_contents(0);
171  std::vector<float> base_errors(0);
172  std::vector<float> base_entries(0);
173  std::vector<float> noise_contents(0);
174  std::vector<float> noise_errors(0);
175  std::vector<float> noise_entries(0);
176  float peak_max = -1.*sistrip::invalid_;
177  float peak_min = 1.*sistrip::invalid_;
178  float base_max = -1.*sistrip::invalid_;
179  float base_min = 1.*sistrip::invalid_;
180  float noise_max = -1.*sistrip::invalid_;
181  float noise_min = 1.*sistrip::invalid_;
182 
183  // Transfer histogram contents/errors/stats to containers
184  for ( uint16_t ibin = 0; ibin < nbins; ibin++ ) {
185 
186  // Peak histogram
187  peak_contents.push_back( peak_histo->GetBinContent(ibin+1) );
188  peak_errors.push_back( peak_histo->GetBinError(ibin+1) );
189  peak_entries.push_back( peak_histo->GetBinEntries(ibin+1) );
190  if ( peak_entries[ibin] ) {
191  if ( peak_contents[ibin] > peak_max ) { peak_max = peak_contents[ibin]; }
192  if ( peak_contents[ibin] < peak_min && ibin ) { peak_min = peak_contents[ibin]; }
193  }
194 
195  // Base histogram
196  base_contents.push_back( base_histo->GetBinContent(ibin+1) );
197  base_errors.push_back( base_histo->GetBinError(ibin+1) );
198  base_entries.push_back( base_histo->GetBinEntries(ibin+1) );
199  if ( base_entries[ibin] ) {
200  if ( base_contents[ibin] > base_max ) { base_max = base_contents[ibin]; }
201  if ( base_contents[ibin] < base_min && ibin ) { base_min = base_contents[ibin]; }
202  }
203 
204  // Noise histogram
205  noise_contents.push_back( noise_histo->GetBinContent(ibin+1) );
206  noise_errors.push_back( noise_histo->GetBinError(ibin+1) );
207  noise_entries.push_back( noise_histo->GetBinEntries(ibin+1) );
208  if ( noise_entries[ibin] ) {
209  if ( noise_contents[ibin] > noise_max ) { noise_max = noise_contents[ibin]; }
210  if ( noise_contents[ibin] < noise_min && ibin ) { noise_min = noise_contents[ibin]; }
211  }
212 
213  }
214 
215  // Find "zero light" level and error
216  //@@ record bias setting used for zero light level
217  //@@ zero light error changes wrt gain setting ???
218  float zero_light_level = sistrip::invalid_;
219  float zero_light_error = sistrip::invalid_;
220  for ( uint16_t ibin = 0; ibin < nbins; ibin++ ) {
221  if ( base_entries[ibin] ) {
222  zero_light_level = base_contents[ibin];
223  zero_light_error = base_errors[ibin];
224  break;
225  }
226  }
227 
228  float zero_light_thres = sistrip::invalid_;
229  if ( zero_light_level <= sistrip::maximum_ &&
230  zero_light_error <= sistrip::maximum_ ) {
231  zero_light_thres = zero_light_level + 5. * zero_light_error;
232  } else {
233  std::stringstream ss;
234  ss << sistrip::invalidZeroLightLevel_ << "ForGain" << igain;
235  anal->addErrorCode( ss.str() );
236  continue;
237  }
238 
239  // Find range of base histogram
240  float base_range = base_max - base_min;
241 
242  // Find overlapping max/min region that constrains range of linear fit
243  float max = peak_max < base_max ? peak_max : base_max;
244  float min = peak_min > base_min ? peak_min : base_min;
245  float range = max - min;
246 
247  // Container identifying whether samples from 'base' histo are above "zero light"
248  std::vector<bool> above_zero_light;
249  above_zero_light.resize(3,true);
250 
251  // Linear fits to top of peak and base curves and one to bottom of base curve
252  sistrip::LinearFit peak_high;
253  sistrip::LinearFit base_high;
254  sistrip::LinearFit base_low;
255 
256  // Iterate through histogram bins
257  uint16_t peak_bin = 0;
258  uint16_t base_bin = 0;
259  uint16_t low_bin = 0;
260  for ( uint16_t ibin = 0; ibin < nbins; ibin++ ) {
261 
262  // Record whether consecutive samples from 'base' histo are above the "zero light" level
263  if ( base_entries[ibin] ) {
264  above_zero_light.erase( above_zero_light.begin() );
265  if ( base_contents[ibin] > zero_light_thres ) { above_zero_light.push_back( true ); }
266  else { above_zero_light.push_back( false ); }
267  if ( above_zero_light.size() != 3 ) { above_zero_light.resize(3,false); }
268  }
269 
270  // Linear fit to peak histogram
271  if ( peak_entries[ibin] &&
272  peak_contents[ibin] > ( min + 0.2*range ) &&
273  peak_contents[ibin] < ( min + 0.8*range ) ) {
274  if ( !peak_bin ) { peak_bin = ibin; }
275  if ( ( ibin - peak_bin ) < 10 ) {
276  peak_high.add( ibin, peak_contents[ibin] ); //@@ should weight using bin error or bin contents (sqrt(N)/N)
277  }
278  }
279  // Linear fit to base histogram
280  if ( base_entries[ibin] &&
281  base_contents[ibin] > ( min + 0.2*range ) &&
282  base_contents[ibin] < ( min + 0.8*range ) ) {
283  if ( !base_bin ) { base_bin = ibin; }
284  if ( ( ibin - base_bin ) < 10 ) {
285  base_high.add( ibin, base_contents[ibin] ); //@@ should weight using bin error or bin contents (sqrt(N)/N)
286  }
287  }
288  // Low linear fit to base histogram
289  if ( base_entries[ibin] &&
290  //@@ above_zero_light[0] && above_zero_light[1] && above_zero_light[2] &&
291  base_contents[ibin] > ( base_min + 0.2*base_range ) &&
292  base_contents[ibin] < ( base_min + 0.6*base_range ) ) {
293  if ( !low_bin ) { low_bin = ibin; }
294  if ( ( ibin - low_bin ) < 10 ) {
295  base_low.add( ibin, base_contents[ibin] ); //@@ should weight using bin error or bin contents (sqrt(N)/N)
296  }
297  }
298 
299  }
300 
301  // Extract width between two curves at midpoint within range
302  float mid = min + 0.5*range;
303  sistrip::LinearFit::Params peak_params;
304  sistrip::LinearFit::Params base_params;
305  peak_high.fit( peak_params );
306  base_high.fit( base_params );
307 
308  float peak_pos = sistrip::invalid_;
309  float base_pos = sistrip::invalid_;
310  float width = sistrip::invalid_;
311  if ( peak_params.b_ > 0. ) {
312  peak_pos = ( mid - peak_params.a_ ) / peak_params.b_;
313  }
314  if ( base_params.b_ > 0. ) {
315  base_pos = ( mid - base_params.a_ ) / base_params.b_;
316  }
317  if ( base_pos < sistrip::valid_ &&
318  peak_pos < sistrip::valid_ ) {
319  width = base_pos - peak_pos;
320  }
321 
322  // Extrapolate to zero light to find "lift off"
323  sistrip::LinearFit::Params low_params;
324  base_low.fit( low_params );
325  float lift_off = sistrip::invalid_;
326  if ( low_params.b_ > 0. ) {
327  lift_off = ( zero_light_level - low_params.a_ ) / low_params.b_;
328  }
329 
330  // ---------- Set all parameters ----------
331 
332  // Slope of baseline
333  if ( low_params.b_ > 0. ) {
334  anal->baseSlope_[igain] = low_params.b_;
335  }
336 
337  // Check "lift off" value and set bias setting accordingly
338  if ( lift_off <= sistrip::maximum_ ) {
339  anal->bias_[igain] = static_cast<uint16_t>( lift_off ) + 2;
340  } else { anal->bias_[igain] = OptoScanAnalysis::defaultBiasSetting_; }
341 
342  // Calculate "lift off" (in mA)
343  if ( lift_off <= sistrip::maximum_ ) {
344  anal->liftOff_[igain] = 0.45 * lift_off;
345  }
346 
347  // Calculate laser threshold (in mA)
348  if ( width < sistrip::invalid_ ) {
349  anal->threshold_[igain] = 0.45 * ( lift_off - width/2. );
350  }
351 
352  // Set "zero light" level
353  anal->zeroLight_[igain] = zero_light_level;
354 
355  // Set link noise
356  uint16_t bin_number = sistrip::invalid_;
357  if ( anal->threshold_[igain] < sistrip::valid_ ) {
358  // Old calculation, used in commissioning in 2008
359  // always leads to zero link noise
360  // bin_number = static_cast<uint16_t>( anal->threshold_[igain] / 0.45 );
361  // New calculation asked by Karl et al, for commissioning in 2009
362  bin_number = (uint16_t) (lift_off + width / 3.);
363  }
364  if ( bin_number < sistrip::valid_ ) {
365  if ( bin_number < noise_contents.size() ) {
366  anal->linkNoise_[igain] = noise_contents[bin_number];
368  }
369 
370  // Calculate tick mark height
371  if ( low_params.b_ <= sistrip::maximum_ &&
372  width <= sistrip::maximum_ ) {
373  anal->tickHeight_[igain] = width * low_params.b_;
374  }
375 
376  // Set measured gain
377  if ( anal->tickHeight_[igain] < sistrip::invalid_-1. ) {
378  anal->measGain_[igain] = anal->tickHeight_[igain] * OptoScanAnalysis::fedAdcGain_ / 0.800;
379  } else { anal->measGain_[igain] = sistrip::invalid_; }
380 
381  } // gain loop
382 
383  // Iterate through four gain settings and identify optimum gain setting
384  const float target_gain = targetGain_; //0.863; // changed from 0.8 to avoid choice of low tickheights (Xtof, SL, 15/6/2009)
385 
386  float diff_in_gain = sistrip::invalid_;
387  for ( uint16_t igain = 0; igain < 4; igain++ ) {
388 
389  // Check for sensible gain value
390  if ( anal->measGain_[igain] > sistrip::maximum_ ) { continue; }
391 
392  // Find optimum gain setting
393  if ( fabs( anal->measGain_[igain] - target_gain ) < diff_in_gain ) {
394  anal->gain_ = igain;
395  diff_in_gain = fabs( anal->measGain_[igain] - target_gain );
396  }
397 
398  }
399 
400  // Check optimum gain setting
402 
403 }
404 
405 // ----------------------------------------------------------------------------
406 //
408  const uint16_t& digital_level ) const {
409  if ( gain <= 3 && digital_level <= 1 ) { return histos_[gain][digital_level]; }
410  else { return Histo(0,""); }
411 }
Histogram-based analysis for opto bias/gain scan.
std::vector< std::vector< Histo > > histos_
static const char unexpectedTask_[]
const uint32_t & fedKey() const
Histo histo(const uint16_t &gain, const uint16_t &digital_level) const
Utility class that holds histogram title.
std::pair< TH1 *, std::string > Histo
static const char numberOfHistos_[]
static const char unexpectedExtraInfo_[]
static const uint16_t valid_
Definition: Constants.h:17
void extract(const std::vector< TH1 * > &)
static const char numberOfBins_[]
sistrip classes
static const char baselineRms_[]
static const float fedAdcGain_
static const char mlCommissioning_[]
uint32_t extractFedKey(const TH1 *const )
static const uint16_t maximum_
Definition: Constants.h:20
virtual void addErrorCode(const std::string &error)
T min(T a, T b)
Definition: MathUtil.h:58
static const char digital_[]
void add(const float &value_x, const float &value_y)
Definition: Utility.cc:17
static const uint16_t defaultBiasSetting_
static const uint16_t invalid_
Definition: Constants.h:16
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
static const char gain_[]
Abstract base for derived classes that provide analysis of commissioning histograms.
static const uint16_t defaultGainSetting_
void fit(Params &fit_params)
Definition: Utility.cc:47
static const char unexpectedBinNumber_[]
CommissioningAnalysis *const anal() const
static const char nullPtr_[]
static const char invalidZeroLightLevel_[]