Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226
#include "DQM/SiStripCommissioningAnalysis/interface/SamplingAlgorithm.h"
#include "CondFormats/SiStripObjects/interface/SamplingAnalysis.h"
#include "DataFormats/SiStripCommon/interface/SiStripHistoTitle.h"
#include "DataFormats/SiStripCommon/interface/SiStripEnumsAndStrings.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "DQM/SiStripCommissioningAnalysis/interface/SiStripPulseShape.h"
#include "TProfile.h"
#include "TF1.h"
#include <iostream>
#include <sstream>
#include <iomanip>
#include <cmath>

using namespace sistrip;

// ----------------------------------------------------------------------------
//
SamplingAlgorithm::SamplingAlgorithm(const edm::ParameterSet& pset, SamplingAnalysis* const anal, uint32_t latencyCode)
    : CommissioningAlgorithm(anal),
      histo_(nullptr, ""),
      deconv_fitter_(nullptr),
      peak_fitterA_(nullptr),
      peak_fitterB_(nullptr),
      latencyCode_(latencyCode),
      samp_(nullptr) {
  peak_fitterA_ = new TF1("peak_fitterA", fpeak_convoluted, -4800, 0, 5);
  peak_fitterA_->SetNpx(2000);
  peak_fitterA_->FixParameter(0, 0);
  peak_fitterA_->SetParLimits(1, 0, 4800);
  peak_fitterA_->SetParLimits(2, 0, 20);
  peak_fitterA_->FixParameter(3, 50);
  peak_fitterA_->SetParLimits(4, 0, 75);
  peak_fitterA_->SetParameters(0., 1250, 10, 50, 10);

  peak_fitterB_ = new TF1("peak_fitterB", fpeak_convoluted, -100, 100, 5);
  peak_fitterB_->SetNpx(200);
  peak_fitterB_->FixParameter(0, 0);
  peak_fitterB_->SetParLimits(1, -100, 100);
  peak_fitterB_->SetParLimits(2, 0, 20);
  peak_fitterB_->FixParameter(3, 50);
  peak_fitterB_->SetParLimits(4, 0, 75);
  peak_fitterB_->SetParameters(0., -50, 10, 50, 10);

  deconv_fitter_ = new TF1("deconv_fitter", fdeconv_convoluted, -50, 50, 5);
  deconv_fitter_->SetNpx(1000);
  deconv_fitter_->FixParameter(0, 0);
  deconv_fitter_->SetParLimits(1, -50, 50);
  deconv_fitter_->SetParLimits(2, 0, 200);
  deconv_fitter_->SetParLimits(3, 5, 100);
  deconv_fitter_->FixParameter(3, 50);
  deconv_fitter_->SetParLimits(4, 0, 75);
  deconv_fitter_->SetParameters(0., -2.82, 0.96, 50, 20);
}

// ----------------------------------------------------------------------------
//
void SamplingAlgorithm::extract(const std::vector<TH1*>& histos) {
  if (!anal()) {
    edm::LogWarning(mlCommissioning_) << "[SamplingAlgorithm::" << __func__ << "]"
                                      << " NULL pointer to Analysis object!";
    return;
  }

  CommissioningAnalysis* tmp = const_cast<CommissioningAnalysis*>(anal());
  samp_ = dynamic_cast<SamplingAnalysis*>(tmp);
  if (!samp_) {
    edm::LogWarning(mlCommissioning_) << "[SamplingAlgorithm::" << __func__ << "]"
                                      << " NULL pointer to derived Analysis object!";
    return;
  }

  // Check
  if (histos.size() != 1 && histos.size() != 2) {
    samp_->addErrorCode(sistrip::numberOfHistos_);
  }

  // Extract FED key from histo title
  if (!histos.empty()) {
    samp_->fedKey(extractFedKey(histos.front()));
  }

  // Extract
  std::vector<TH1*>::const_iterator ihis = histos.begin();
  for (; ihis != histos.end(); ihis++) {
    // Check pointer
    if (!(*ihis)) {
      edm::LogWarning(mlCommissioning_) << " NULL pointer to histogram!";
      continue;
    }

    // Check name
    SiStripHistoTitle title((*ihis)->GetName());
    if (title.runType() != sistrip::APV_LATENCY && title.runType() != sistrip::FINE_DELAY) {
      samp_->addErrorCode(sistrip::unexpectedTask_);
      continue;
    }
    // Set the mode for later fits
    samp_->runType_ = title.runType();

    // Set the granularity
    samp_->granularity_ = title.granularity();

    // Extract timing histo
    if (title.extraInfo().find(sistrip::extrainfo::clusterCharge_) != std::string::npos) {
      histo_.first = *ihis;
      histo_.second = (*ihis)->GetName();
    }
  }
}

// ----------------------------------------------------------------------------
//
void SamplingAlgorithm::analyse() {
  if (!samp_) {
    edm::LogWarning(mlCommissioning_) << "[SamplingAlgorithm::" << __func__ << "]"
                                      << " NULL pointer to derived Analysis object!";
    return;
  }

  TProfile* prof = (TProfile*)(histo_.first);
  if (!prof) {
    edm::LogWarning(mlCommissioning_) << " NULL pointer to histogram!";
    return;
  }

  // set the right error mode: rms
  prof->SetErrorOption(" ");

  //that should not be needed, but it seems histos are stored with error option " " and errors "s" in all cases.
  //it MUST be removed if the DQM (?) bug is solved
  for (int i = 0; i < prof->GetNbinsX(); ++i) {
    if (prof->GetBinEntries(i) > 0)
      prof->SetBinError(i, prof->GetBinError(i) / sqrt(prof->GetBinEntries(i)));
  }

  // prune the profile
  pruneProfile(prof);

  // correct for the binning
  correctBinning(prof);

  // correct for clustering effects
  correctProfile(prof, samp_->sOnCut_);

  // fit depending on the mode
  if (samp_->runType_ == sistrip::APV_LATENCY) {
    // initialize  the fit (overal latency)
    float max = prof->GetBinCenter(prof->GetMaximumBin());
    float ampl = prof->GetMaximum();
    peak_fitterA_->SetParameters(0., 50 - max, ampl / 20., 50, 10);

    // fit
    if (prof->Fit(peak_fitterA_, "Q") == 0)
      prof->Fit(peak_fitterA_, "QEM");

    // Set monitorables
    samp_->max_ = peak_fitterA_->GetMaximumX();
    samp_->error_ = peak_fitterA_->GetParError(1);

  } else {  // sistrip::FINE_DELAY

    // initialize  the fit (overal latency)
    float max = prof->GetBinCenter(prof->GetMaximumBin());
    float ampl = prof->GetMaximum();
    deconv_fitter_->SetParameters(0., -max, ampl / 10., 50, 20);
    peak_fitterB_->SetParameters(0., 50 - max, ampl / 20., 50, 10);
    if (latencyCode_ & 0x80) {  // deconv mode
      // fit
      if (prof->Fit(deconv_fitter_, "Q") == 0)
        prof->Fit(deconv_fitter_, "QEM");
      // Set monitorables
      samp_->max_ = deconv_fitter_->GetMaximumX();
      samp_->error_ = deconv_fitter_->GetParError(1);
    } else {  // peak mode
      // fit
      if (prof->Fit(peak_fitterB_, "Q") == 0)
        prof->Fit(peak_fitterB_, "QEM");
      // Set monitorables
      samp_->max_ = peak_fitterB_->GetMaximumX();
      samp_->error_ = peak_fitterB_->GetParError(1);
    }
  }
}

// ----------------------------------------------------------------------------
//
void SamplingAlgorithm::pruneProfile(TProfile* profile) const {
  // loop over bins to find the max stat
  uint32_t nbins = profile->GetNbinsX();
  uint32_t max = 0;
  for (uint32_t bin = 1; bin <= nbins; ++bin)
    max = max < profile->GetBinEntries(bin) ? uint32_t(profile->GetBinEntries(bin)) : max;
  // loop over bins to clean
  uint32_t min = max / 10;
  for (uint32_t bin = 1; bin <= nbins; ++bin)
    if (profile->GetBinEntries(bin) < min) {
      profile->SetBinContent(bin, 0.);
      profile->SetBinError(bin, 0.);
    }
}

// ----------------------------------------------------------------------------
//
void SamplingAlgorithm::correctBinning(TProfile* prof) const {
  prof->GetXaxis()->SetLimits(prof->GetXaxis()->GetXmin() - prof->GetBinWidth(1) / 2.,
                              prof->GetXaxis()->GetXmax() - prof->GetBinWidth(1) / 2.);
}

// ----------------------------------------------------------------------------
//
void SamplingAlgorithm::correctProfile(TProfile* profile, float SoNcut) const {
  if (!samp_) {
    return;
  }
  uint32_t nbins = profile->GetNbinsX();
  float min = samp_->limit(SoNcut);
  for (uint32_t bin = 1; bin <= nbins; ++bin)
    if (profile->GetBinContent(bin) < min) {
      profile->SetBinContent(bin, 0.);
      profile->SetBinError(bin, 0.);
      profile->SetBinEntries(bin, 0);
    } else {
      profile->SetBinContent(
          bin, profile->GetBinEntries(bin) * samp_->correctMeasurement(profile->GetBinContent(bin), SoNcut));
    }
}