Macros

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 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570
//
//  The implementation of the PixelResolutinHistograms.cc class.  Please
//  look at PixelResolutionHistograms.h header file for the interface.
//
//------------------------------------------------------------------------------

// The switch, undefined in CMSSW release, and defined by standalone compilation script:

#ifdef SI_PIXEL_TEMPLATE_STANDALONE
//
//--- Stand-alone: Include a the header file from the local directory, as well as
//    dummy implementations of SimpleHistogramGenerator, LogInfo, LogError and LogDebug...
//
#include "PixelResolutionHistograms.h"
//
class TH1F;
class TH2F;
class SimpleHistogramGenerator {
public:
  SimpleHistogramGenerator(TH1F* hist) : hist_(hist) {}

private:
  TH1F* hist_;  // we don't own it
};
#define LOGDEBUG std::cout
#define LOGERROR std::cout
#define LOGINFO std::cout
//
#else
//--- We're inside a CMSSW release: Include the real thing.
//
#include "FastSimulation/TrackingRecHitProducer/interface/PixelResolutionHistograms.h"
#include "FastSimulation/Utilities/interface/RandomEngineAndDistribution.h"
#include "FastSimulation/Utilities/interface/SimpleHistogramGenerator.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
//
#define LOGDEBUG LogDebug("")
#define LOGERROR edm::LogError("Error")
#define LOGINFO edm::LogInfo("Info")
//
#endif

// Generic C stuff
#include <cmath>
#include <iostream>
#include <string>

// ROOT
#include <TFile.h>
#include <TH1F.h>
#include <TH2F.h>

// Global definitions
const float cmtomicron = 10000.0;

//------------------------------------------------------------------------------
//  Constructor: Books the FastSim Histograms, given the input parameters
//  which are provided as arguments. These variables are then const inside
//  the class. (That is, once we make the histograms, we can't change the
//  definition of the binning.)
//------------------------------------------------------------------------------
PixelResolutionHistograms::PixelResolutionHistograms(std::string filename,   // ROOT file for histograms
                                                     std::string rootdir,    // Subdirectory in the file, "" if none
                                                     std::string descTitle,  // Descriptive title
                                                     unsigned int detType,   // Where we are... (&&& do we need this?)
                                                     double cotbetaBinWidth,
                                                     double cotbetaLowEdge,
                                                     int cotbetaBins,
                                                     double cotalphaBinWidth,
                                                     double cotalphaLowEdge,
                                                     int cotalphaBins)  //,
                                                                        //int	qbinWidth,
                                                                        //int	qbins )
    : weOwnHistograms_(true),                                           // we'll be making some histos
      detType_(detType),
      cotbetaBinWidth_(cotbetaBinWidth),
      cotbetaLowEdge_(cotbetaLowEdge),
      cotbetaBins_(cotbetaBins),
      cotalphaBinWidth_(cotalphaBinWidth),
      cotalphaLowEdge_(cotalphaLowEdge),
      cotalphaBins_(cotalphaBins),
      qbinWidth_(1),
      qbins_(4),
      binningHisto_(nullptr),
      resMultiPixelXHist_(),
      resSinglePixelXHist_(),  // all to nullptr
      resMultiPixelYHist_(),
      resSinglePixelYHist_(),  // all to nullptr
      qbinHist_(),             // all to nullptr
      file_(nullptr),
      status_(0),
      resMultiPixelXGen_(),
      resSinglePixelXGen_(),
      resMultiPixelYGen_(),
      resSinglePixelYGen_(),
      qbinGen_() {
  file_ = std::make_unique<TFile>(filename.c_str(), "RECREATE");
  //Resolution binning
  // const double 	cotbetaBinWidth = 1.0;
  // const double 	cotbetaLowEdge	= -11.5 ;
  // const int	cotbetaBins	= 23;
  // const double	cotalphaBinWidth	= 0.08 ;
  // const double	cotalphaLowEdge	= -0.36 ;
  // const int	cotalphaBins	= 9;
  // const int	qbinWidth	= 1;
  // const int	qbins		= 4;

  // Dummy 2D histogram to store binning:
  binningHisto_ = new TH2F("ResHistoBinning",
                           descTitle.c_str(),
                           cotbetaBins,
                           cotbetaLowEdge,
                           cotbetaLowEdge + cotbetaBins * cotbetaBinWidth,
                           cotalphaBins,
                           cotalphaLowEdge,
                           cotalphaLowEdge + cotalphaBins * cotalphaBinWidth);

  // Store detType in the underflow bin
  binningHisto_->SetBinContent(0, 0, detType_);

  // All other histograms:
  Char_t histo[200];
  Char_t title[200];
  //
  //--- Histograms for clusters with multiple pixels hit in a given direction.
  //
  for (int ii = 0; ii < cotbetaBins_; ii++) {
    for (int jj = 0; jj < cotalphaBins_; jj++) {
      for (int kk = 0; kk < qbins_; kk++) {
        //
        sprintf(histo, "hx%d1%02d%d%d", detType_, ii + 1, jj + 1, kk + 1);  //information of bits of histogram names
        //--- First bit 1/0 barrel/forward, second 1/0 multi/single, cotbeta, cotalpha, qbins
        sprintf(title,
                "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin %d npixel>1 X",
                cotbetaLowEdge_ + ii * cotbetaBinWidth_,
                cotbetaLowEdge_ + (ii + 1) * cotbetaBinWidth_,
                cotalphaLowEdge_ + jj * cotalphaBinWidth_,
                cotalphaLowEdge_ + (jj + 1) * cotalphaBinWidth_,
                kk + 1);
        //
        resMultiPixelXHist_[ii][jj][kk] = new TH1F(histo, title, 1000, -0.05, 0.05);

        sprintf(histo, "hy%d1%02d%d%d", detType_, ii + 1, jj + 1, kk + 1);
        sprintf(title,
                "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin %d npixel>1 Y",
                cotbetaLowEdge_ + ii * cotbetaBinWidth_,
                cotbetaLowEdge_ + (ii + 1) * cotbetaBinWidth_,
                cotalphaLowEdge_ + jj * cotalphaBinWidth_,
                cotalphaLowEdge_ + (jj + 1) * cotalphaBinWidth_,
                kk + 1);
        //
        resMultiPixelYHist_[ii][jj][kk] = new TH1F(histo, title, 1000, -0.05, 0.05);
      }
    }
  }

  //
  //--- Histograms for clusters where only a single pixel was hit in a given direction.
  //
  for (int ii = 0; ii < cotbetaBins_; ii++) {
    for (int jj = 0; jj < cotalphaBins_; jj++) {
      sprintf(histo, "hx%d0%02d%d", detType_, ii + 1, jj + 1);  //information of bits of histogram names
      //first bit 1/0 barrel/forward, second 1/0 multi/single, cotbeta, cotalpha
      sprintf(title,
              "cotbeta %.1f-%.1f cotalpha %.2f-%.2f npixel=1 X",
              cotbetaLowEdge_ + ii * cotbetaBinWidth_,
              cotbetaLowEdge_ + (ii + 1) * cotbetaBinWidth_,
              cotalphaLowEdge_ + jj * cotalphaBinWidth_,
              cotalphaLowEdge_ + (jj + 1) * cotalphaBinWidth_);
      //
      resSinglePixelXHist_[ii][jj] = new TH1F(histo, title, 1000, -0.05, 0.05);

      sprintf(histo, "hy%d0%02d%d", detType_, ii + 1, jj + 1);
      sprintf(title,
              "cotbeta %.1f-%.1f cotalpha %.2f-%.2f npixel=1 Y",
              cotbetaLowEdge_ + ii * cotbetaBinWidth_,
              cotbetaLowEdge_ + (ii + 1) * cotbetaBinWidth_,
              cotalphaLowEdge_ + jj * cotalphaBinWidth_,
              cotalphaLowEdge_ + (jj + 1) * cotalphaBinWidth_);
      //
      resSinglePixelYHist_[ii][jj] = new TH1F(histo, title, 1000, -0.05, 0.05);

      sprintf(histo, "hqbin%d%02d%d", detType_, ii + 1, jj + 1);
      sprintf(title,
              "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin",
              cotbetaLowEdge_ + ii * cotbetaBinWidth_,
              cotbetaLowEdge_ + (ii + 1) * cotbetaBinWidth_,
              cotalphaLowEdge_ + jj * cotalphaBinWidth_,
              cotalphaLowEdge_ + (jj + 1) * cotalphaBinWidth_);
      //
      qbinHist_[ii][jj] = new TH1F(histo, title, 4, -0.49, 3.51);
    }
  }
}

//------------------------------------------------------------------------------
//  Another constructor: load the histograms from one file.
//     filename = full path to filename
//     rootdir  = ROOT directory inside the file
//
//  The other parameters are the same (needed later) and must correspond
//  to the histograms we are loading from the file.
//------------------------------------------------------------------------------
PixelResolutionHistograms::PixelResolutionHistograms(std::string filename,
                                                     std::string rootdir,
                                                     int detType,
                                                     bool ignore_multi,
                                                     bool ignore_single,
                                                     bool ignore_qBin)
    : weOwnHistograms_(false),  // resolution histograms are owned by the ROOT file
      detType_(-1),
      cotbetaBinWidth_(0),
      cotbetaLowEdge_(0),
      cotbetaBins_(0),
      cotalphaBinWidth_(0),
      cotalphaLowEdge_(0),
      cotalphaBins_(0),
      qbinWidth_(1),
      qbins_(4),
      binningHisto_(nullptr),
      resMultiPixelXHist_(),
      resSinglePixelXHist_(),  // all to nullptr
      resMultiPixelYHist_(),
      resSinglePixelYHist_(),  // all to nullptr
      qbinHist_(),             // all to nullptr
      file_(nullptr),
      status_(0),
      resMultiPixelXGen_(),
      resSinglePixelXGen_(),
      resMultiPixelYGen_(),
      resSinglePixelYGen_(),
      qbinGen_() {
  Char_t histo[200];        // the name of the histogram
  Char_t title[200];        // histo title, for debugging and sanity checking (compare inside file)
  TH1F* tmphist = nullptr;  // cache for histo pointer

  //--- Open the file for reading.
  file_ = std::make_unique<TFile>(filename.c_str(), "READ");
  if (!file_) {
    status_ = 1;
    LOGERROR << "PixelResolutionHistograms:: Error, file " << filename << " not found.";
    return;  // PixelTemplateSmearerBase will throw an exception upon our return.
  }

  //--- The dummy 2D histogram with the binning of cot\beta and cot\alpha:
  binningHisto_ = (TH2F*)file_->Get(Form("%s%s", rootdir.c_str(), "ResHistoBinning"));
  if (!binningHisto_) {
    status_ = 11;
    LOGERROR << "PixelResolutionHistograms:: Error, binning histogrram ResHistoBinning not found.";
    return;  // PixelTemplateSmearerBase will throw an exception upon our return.
  }

  if (detType == -1) {
    //--- Fish out detType from the underflow bin:
    detType_ = binningHisto_->GetBinContent(0, 0);
  } else {
    detType_ = detType;  // constructor's argument overrides what's in ResHistoBinning histogram.
  }

  //--- Now we fill the binning variables:
  cotbetaAxis_ = binningHisto_->GetXaxis();
  cotbetaBinWidth_ = binningHisto_->GetXaxis()->GetBinWidth(1);  // assume all same width
  cotbetaLowEdge_ = binningHisto_->GetXaxis()->GetXmin();        // low edge of the first bin
  cotbetaBins_ = binningHisto_->GetXaxis()->GetNbins();
  cotalphaAxis_ = binningHisto_->GetYaxis();
  cotalphaBinWidth_ = binningHisto_->GetYaxis()->GetBinWidth(1);  // assume all same width;
  cotalphaLowEdge_ = binningHisto_->GetYaxis()->GetXmin();        // low edge of the first bin;
  cotalphaBins_ = binningHisto_->GetYaxis()->GetNbins();

  if (!ignore_multi) {
    //
    //--- Histograms for clusters with multiple pixels hit in a given direction.
    //
    for (int ii = 0; ii < cotbetaBins_; ii++) {
      for (int jj = 0; jj < cotalphaBins_; jj++) {
        for (int kk = 0; kk < qbins_; kk++) {
          //
          sprintf(histo, "hx%d1%02d%d%d", detType_, ii + 1, jj + 1, kk + 1);  //information of bits of histogram names
          //--- First bit 1/0 barrel/forward, second 1/0 multi/single, cotbeta, cotalpha, qbins
          sprintf(title,
                  "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin %d npixel>1 X",
                  cotbetaLowEdge_ + ii * cotbetaBinWidth_,
                  cotbetaLowEdge_ + (ii + 1) * cotbetaBinWidth_,
                  cotalphaLowEdge_ + jj * cotalphaBinWidth_,
                  cotalphaLowEdge_ + (jj + 1) * cotalphaBinWidth_,
                  kk + 1);
          //
          tmphist = (TH1F*)file_->Get(Form("%s%s", rootdir.c_str(), histo));
          if (!tmphist) {
            status_ = 2;
            LOGERROR << "Failed to find histogram=" << std::string(histo);
            return;
          }
          LOGDEBUG << "Found histo " << std::string(histo) << " with title = " << std::string(tmphist->GetTitle())
                   << std::endl;
          if (tmphist->GetEntries() < 5) {
            LOGINFO << "Histogram " << std::string(histo) << " has only " << tmphist->GetEntries()
                    << " entries. Trouble ahead." << std::endl;
          }
          resMultiPixelXHist_[ii][jj][kk] = tmphist;
          resMultiPixelXGen_[ii][jj][kk] = new SimpleHistogramGenerator(tmphist);

          sprintf(histo, "hy%d1%02d%d%d", detType_, ii + 1, jj + 1, kk + 1);
          sprintf(title,
                  "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin %d npixel>1 Y",
                  cotbetaLowEdge_ + ii * cotbetaBinWidth_,
                  cotbetaLowEdge_ + (ii + 1) * cotbetaBinWidth_,
                  cotalphaLowEdge_ + jj * cotalphaBinWidth_,
                  cotalphaLowEdge_ + (jj + 1) * cotalphaBinWidth_,
                  kk + 1);
          //
          tmphist = (TH1F*)file_->Get(Form("%s%s", rootdir.c_str(), histo));
          if (!tmphist) {
            status_ = 3;
            LOGERROR << "Failed to find histogram=" << std::string(histo);
            return;
          }
          LOGDEBUG << "Found histo " << std::string(histo) << " with title = " << std::string(tmphist->GetTitle())
                   << std::endl;
          if (tmphist->GetEntries() < 5) {
            LOGINFO << "Histogram " << std::string(histo) << " has only " << tmphist->GetEntries()
                    << " entries. Trouble ahead." << std::endl;
          }
          resMultiPixelYHist_[ii][jj][kk] = tmphist;
          resMultiPixelYGen_[ii][jj][kk] = new SimpleHistogramGenerator(tmphist);
        }
      }
    }
    //
  }  // if (not ignore multi)

  //
  //--- Histograms for clusters where only a single pixel was hit in a given direction.
  //
  for (int ii = 0; ii < cotbetaBins_; ii++) {
    for (int jj = 0; jj < cotalphaBins_; jj++) {
      //--- Single pixel, along X.
      //
      sprintf(histo, "hx%d0%02d%d", detType_, ii + 1, jj + 1);  //information of bits of histogram names
      //--- First bit 1/0 barrel/forward, second 1/0 multi/single, cotbeta, cotalpha
      sprintf(title,
              "cotbeta %.1f-%.1f cotalpha %.2f-%.2f npixel=1 X",
              cotbetaLowEdge_ + ii * cotbetaBinWidth_,
              cotbetaLowEdge_ + (ii + 1) * cotbetaBinWidth_,
              cotalphaLowEdge_ + jj * cotalphaBinWidth_,
              cotalphaLowEdge_ + (jj + 1) * cotalphaBinWidth_);
      //
      tmphist = (TH1F*)file_->Get(Form("%s%s", rootdir.c_str(), histo));
      if (!tmphist) {
        if (!ignore_single) {
          LOGERROR << "Failed to find histogram=" << std::string(histo);
          status_ = 4;
          return;
        }
      } else {
        LOGDEBUG << "Found histo " << std::string(histo) << " with title = " << std::string(tmphist->GetTitle())
                 << std::endl;
        LOGDEBUG << "Found histo with title = " << std::string(tmphist->GetTitle()) << std::endl;
        if (tmphist->GetEntries() < 5) {
          LOGINFO << "Histogram " << std::string(histo) << " has only " << tmphist->GetEntries()
                  << " entries. Trouble ahead." << std::endl;
        }
        resSinglePixelXHist_[ii][jj] = tmphist;
        resSinglePixelXGen_[ii][jj] = new SimpleHistogramGenerator(tmphist);
      }

      //--- Single pixel, along Y.
      //
      sprintf(histo, "hy%d0%02d%d", detType_, ii + 1, jj + 1);
      sprintf(title,
              "cotbeta %.1f-%.1f cotalpha %.2f-%.2f npixel=1 Y",
              cotbetaLowEdge_ + ii * cotbetaBinWidth_,
              cotbetaLowEdge_ + (ii + 1) * cotbetaBinWidth_,
              cotalphaLowEdge_ + jj * cotalphaBinWidth_,
              cotalphaLowEdge_ + (jj + 1) * cotalphaBinWidth_);
      //
      tmphist = (TH1F*)file_->Get(Form("%s%s", rootdir.c_str(), histo));
      if (!tmphist) {
        if (!ignore_single) {
          LOGERROR << "Failed to find histogram=" << std::string(histo);
          status_ = 5;
          return;
        }
      } else {
        LOGDEBUG << "Found histo " << std::string(histo) << " with title = " << std::string(tmphist->GetTitle())
                 << std::endl;
        if (tmphist->GetEntries() < 5) {
          LOGINFO << "Histogram " << std::string(histo) << " has only " << tmphist->GetEntries()
                  << " entries. Trouble ahead." << std::endl;
        }
        resSinglePixelYHist_[ii][jj] = tmphist;
        resSinglePixelYGen_[ii][jj] = new SimpleHistogramGenerator(tmphist);
      }

      //--- qBin distribution, for this (cotbeta, cotalpha) bin.
      //
      sprintf(histo, "hqbin%d%02d%d", detType_, ii + 1, jj + 1);
      sprintf(title,
              "cotbeta %.1f-%.1f cotalpha %.2f-%.2f qbin",
              cotbetaLowEdge_ + ii * cotbetaBinWidth_,
              cotbetaLowEdge_ + (ii + 1) * cotbetaBinWidth_,
              cotalphaLowEdge_ + jj * cotalphaBinWidth_,
              cotalphaLowEdge_ + (jj + 1) * cotalphaBinWidth_);
      //
      tmphist = (TH1F*)file_->Get(Form("%s%s", rootdir.c_str(), histo));
      if (!tmphist) {
        if (!ignore_qBin) {
          LOGERROR << "Failed to find histogram=" << std::string(histo);
          status_ = 6;
          return;
        }
      } else {
        LOGDEBUG << "Found histo " << std::string(histo) << " with title = " << std::string(tmphist->GetTitle())
                 << std::endl;
        if (tmphist->GetEntries() < 5) {
          LOGINFO << "Histogram " << std::string(histo) << " has only " << tmphist->GetEntries()
                  << " entries. Trouble ahead." << std::endl;
        }
        qbinHist_[ii][jj] = tmphist;
        qbinGen_[ii][jj] = new SimpleHistogramGenerator(tmphist);
      }
    }
  }
}

//------------------------------------------------------------------------------
//  Destructor.  Use file_ pointer to tell whether we loaded the histograms
//  from a file (and do not own them), or we built them ourselves and thus need
//  to delete them.
//------------------------------------------------------------------------------
PixelResolutionHistograms::~PixelResolutionHistograms() {
  //--- Delete histograms, but only if we own them. If
  //--- they came from a file, let them be.
  //
  if (!weOwnHistograms_) {
    //--- Read the histograms from the TFile, the file will take care of them.
    file_->Close();
    /// delete file_ ;   // no need to delete if unique_ptr<>
    /// file_ = 0;
  } else {
    //--- We made the histograms, so first write them inthe output ROOT file and close it.
    LOGINFO << "PixelResHistoStore: Writing the histograms to the output file. "  // << filename
            << std::endl;
    file_->Write();
    file_->Close();

    // ROOT file has the ownership, and once the file is closed,
    // all of these guys are deleted.  So, we don't need to do anything.
  }  // else

  //--- Delete FastSim generators. (It's safe to delete a nullptr.)
  for (int ii = 0; ii < cotbetaBins_; ii++) {
    for (int jj = 0; jj < cotalphaBins_; jj++) {
      for (int kk = 0; kk < qbins_; kk++) {
        delete resMultiPixelXGen_[ii][jj][kk];
        delete resMultiPixelYGen_[ii][jj][kk];
      }
    }
  }
  for (int ii = 0; ii < cotbetaBins_; ii++) {
    for (int jj = 0; jj < cotalphaBins_; jj++) {
      delete resSinglePixelXGen_[ii][jj];
      delete resSinglePixelYGen_[ii][jj];
      delete qbinGen_[ii][jj];
    }
  }
}

//------------------------------------------------------------------------------
//  Fills the appropriate FastSim histograms.
//  Returns 0 if the relevant histogram(s) were found and filled, 1 if not.
//------------------------------------------------------------------------------
int PixelResolutionHistograms::Fill(
    double dx, double dy, double cotalpha, double cotbeta, int qbin, int nxpix, int nypix) {
  int icotalpha, icotbeta, iqbin;
  icotalpha = (int)floor((cotalpha - cotalphaLowEdge_) / cotalphaBinWidth_);
  icotbeta = (int)floor((cotbeta - cotbetaLowEdge_) / cotbetaBinWidth_);
  iqbin = qbin > 2 ? 3 : qbin;
  if (icotalpha >= 0 && icotalpha < cotalphaBins_ && icotbeta >= 0 && icotbeta < cotbetaBins_) {
    qbinHist_[icotbeta][icotalpha]->Fill((double)iqbin);
    if (nxpix == 1)
      resSinglePixelXHist_[icotbeta][icotalpha]->Fill(dx / cmtomicron);
    else
      resMultiPixelXHist_[icotbeta][icotalpha][iqbin]->Fill(dx / cmtomicron);
    if (nypix == 1)
      resSinglePixelYHist_[icotbeta][icotalpha]->Fill(dy / cmtomicron);
    else
      resMultiPixelYHist_[icotbeta][icotalpha][iqbin]->Fill(dy / cmtomicron);
  }

  return 0;
}

//------------------------------------------------------------------------------
//  Return the histogram generator for resolution in X.  A generator contains
//  both the histogram and knows how to throw a random number off it.  It is
//  called from FastSim (from PixelTemplateSmearerBase).
//  If cotalpha or cotbeta are outside of the range, return the end of the range.
//------------------------------------------------------------------------------
const SimpleHistogramGenerator* PixelResolutionHistograms::getGeneratorX(double cotalpha,
                                                                         double cotbeta,
                                                                         int qbin,
                                                                         bool single) {
  int icotalpha, icotbeta, iqbin;
  icotalpha = (int)floor((cotalpha - cotalphaLowEdge_) / cotalphaBinWidth_);
  icotbeta = (int)floor((cotbeta - cotbetaLowEdge_) / cotbetaBinWidth_);
  iqbin = qbin > 2 ? 3 : qbin;  // if (qbin>2) then = 3, else return qbin
  //
  //if( icotalpha >= 0 && icotalpha < cotalphaBins_ && icotbeta >= 0 && icotbeta < cotbetaBins_ ) {

  if (icotalpha < 0)
    icotalpha = 0;
  if (icotalpha >= cotalphaBins_)
    icotalpha = cotalphaBins_ - 1;

  if (icotbeta < 0)
    icotbeta = 0;
  if (icotbeta >= cotbetaBins_)
    icotbeta = cotbetaBins_ - 1;

  // At this point we are sure to return *some bin* from the 3D histogram

  if (single)
    return resSinglePixelXGen_[icotbeta][icotalpha];
  else
    return resMultiPixelXGen_[icotbeta][icotalpha][iqbin];

  // }
  //else
  //return nullptr;
}

//------------------------------------------------------------------------------
//  Return the histogram generator for resolution in Y.  A generator contains
//  both the histogram and knows how to throw a random number off it.  It is
//  called from FastSim (from PixelTemplateSmearerBase).
//  If cotalpha or cotbeta are outside of the range, return the end of the range.
//------------------------------------------------------------------------------
const SimpleHistogramGenerator* PixelResolutionHistograms::getGeneratorY(double cotalpha,
                                                                         double cotbeta,
                                                                         int qbin,
                                                                         bool single) {
  int icotalpha, icotbeta, iqbin;
  icotalpha = (int)floor((cotalpha - cotalphaLowEdge_) / cotalphaBinWidth_);
  icotbeta = (int)floor((cotbeta - cotbetaLowEdge_) / cotbetaBinWidth_);
  iqbin = qbin > 2 ? 3 : qbin;  // if (qbin>2) then = 3, else return qbin
  //
  //if( icotalpha >= 0 && icotalpha < cotalphaBins_ && icotbeta >= 0 && icotbeta < cotbetaBins_ ) {

  if (icotalpha < 0)
    icotalpha = 0;
  if (icotalpha >= cotalphaBins_)
    icotalpha = cotalphaBins_ - 1;

  if (icotbeta < 0)
    icotbeta = 0;
  if (icotbeta >= cotbetaBins_)
    icotbeta = cotbetaBins_ - 1;

  // At this point we are sure to return *some bin* from the 3D histogram

  if (single)
    return resSinglePixelYGen_[icotbeta][icotalpha];
  else
    return resMultiPixelYGen_[icotbeta][icotalpha][iqbin];

  //}
  //else
  //return nullptr;
}