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
//
// F.Ratnikov (UMd), Aug. 9, 2005
// Adapted for CASTOR by L. Mundim
//

#include "FWCore/Utilities/interface/typelookup.h"

#include "CalibFormats/CastorObjects/interface/CastorDbService.h"
#include "CalibFormats/CastorObjects/interface/CastorCoderDb.h"
#include "CalibFormats/CastorObjects/interface/QieShape.h"

#include "DataFormats/HcalDetId/interface/HcalGenericDetId.h"

#include <cmath>

CastorDbService::CastorDbService()
    : mPedestals(nullptr),
      mPedestalWidths(nullptr),
      mGains(nullptr),
      mGainWidths(nullptr),
      mQIEData(nullptr),
      mChannelQuality(nullptr),
      mElectronicsMap(nullptr) {}

bool CastorDbService::makeCastorCalibration(const HcalGenericDetId& fId,
                                            CastorCalibrations* fObject,
                                            bool pedestalInADC) const {
  if (fObject) {
    const CastorPedestal* pedestal = getPedestal(fId);
    const CastorGain* gain = getGain(fId);

    if (pedestalInADC) {
      const CastorQIEShape* shape = getCastorShape();
      const CastorQIECoder* coder = getCastorCoder(fId);
      if (pedestal && gain && shape && coder) {
        float pedTrue[4];
        for (int i = 0; i < 4; i++) {
          float x = pedestal->getValues()[i];
          int x1 = (int)std::floor(x);
          int x2 = (int)std::floor(x + 1);
          // y = (y2-y1)/(x2-x1) * (x - x1) + y1  [note: x2-x1=1]
          float y2 = coder->charge(*shape, x2, i);
          float y1 = coder->charge(*shape, x1, i);
          pedTrue[i] = (y2 - y1) * (x - x1) + y1;
        }
        *fObject = CastorCalibrations(gain->getValues(), pedTrue);
        return true;
      }
    } else {
      if (pedestal && gain) {
        *fObject = CastorCalibrations(gain->getValues(), pedestal->getValues());
        return true;
      }
    }
  }
  return false;
}

void CastorDbService::buildCalibrations() {
  // we use the set of ids for pedestals as the master list
  if ((!mPedestals) || (!mGains) || (!mQIEData))
    return;
  std::vector<DetId> ids = mPedestals->getAllChannels();
  bool pedsInADC = mPedestals->isADC();
  // clear the calibrations set
  mCalibSet.clear();
  // loop!
  CastorCalibrations tool;

  //  std::cout << " length of id-vector: " << ids.size() << std::endl;
  for (std::vector<DetId>::const_iterator id = ids.begin(); id != ids.end(); ++id) {
    // make
    bool ok = makeCastorCalibration(*id, &tool, pedsInADC);
    // store
    if (ok)
      mCalibSet.setCalibrations(*id, tool);
    //    std::cout << "Castor calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
  }
  mCalibSet.sort();
}

void CastorDbService::buildCalibWidths() {
  // we use the set of ids for pedestal widths as the master list
  if ((!mPedestalWidths) || (!mGainWidths) || (!mQIEData))
    return;

  std::vector<DetId> ids = mPedestalWidths->getAllChannels();
  bool pedsInADC = mPedestalWidths->isADC();
  // clear the calibrations set
  mCalibWidthSet.clear();
  // loop!
  CastorCalibrationWidths tool;

  //  std::cout << " length of id-vector: " << ids.size() << std::endl;
  for (std::vector<DetId>::const_iterator id = ids.begin(); id != ids.end(); ++id) {
    // make
    bool ok = makeCastorCalibrationWidth(*id, &tool, pedsInADC);
    // store
    if (ok)
      mCalibWidthSet.setCalibrationWidths(*id, tool);
    //    std::cout << "Castor calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
  }
  mCalibWidthSet.sort();
}

bool CastorDbService::makeCastorCalibrationWidth(const HcalGenericDetId& fId,
                                                 CastorCalibrationWidths* fObject,
                                                 bool pedestalInADC) const {
  if (fObject) {
    const CastorPedestalWidth* pedestalwidth = getPedestalWidth(fId);
    const CastorGainWidth* gainwidth = getGainWidth(fId);
    if (pedestalInADC) {
      const CastorQIEShape* shape = getCastorShape();
      const CastorQIECoder* coder = getCastorCoder(fId);
      if (pedestalwidth && gainwidth && shape && coder) {
        float pedTrueWidth[4];
        for (int i = 0; i < 4; i++) {
          float x = pedestalwidth->getWidth(i);
          // assume QIE is linear in low range and use x1=0 and x2=1
          // y = (y2-y1) * (x) [do not add any constant, only scale!]
          float y2 = coder->charge(*shape, 1, i);
          float y1 = coder->charge(*shape, 0, i);
          pedTrueWidth[i] = (y2 - y1) * x;
        }
        *fObject = CastorCalibrationWidths(gainwidth->getValues(), pedTrueWidth);
        return true;
      }
    } else {
      if (pedestalwidth && gainwidth) {
        float pedestalWidth[4];
        for (int i = 0; i < 4; i++)
          pedestalWidth[i] = pedestalwidth->getWidth(i);
        *fObject = CastorCalibrationWidths(gainwidth->getValues(), pedestalWidth);
        return true;
      }
    }
  }
  return false;
}

const CastorPedestal* CastorDbService::getPedestal(const HcalGenericDetId& fId) const {
  if (mPedestals) {
    return mPedestals->getValues(fId);
  }
  return nullptr;
}

const CastorPedestalWidth* CastorDbService::getPedestalWidth(const HcalGenericDetId& fId) const {
  if (mPedestalWidths) {
    return mPedestalWidths->getValues(fId);
  }
  return nullptr;
}

const CastorGain* CastorDbService::getGain(const HcalGenericDetId& fId) const {
  if (mGains) {
    return mGains->getValues(fId);
  }
  return nullptr;
}

const CastorGainWidth* CastorDbService::getGainWidth(const HcalGenericDetId& fId) const {
  if (mGainWidths) {
    return mGainWidths->getValues(fId);
  }
  return nullptr;
}

const CastorQIECoder* CastorDbService::getCastorCoder(const HcalGenericDetId& fId) const {
  if (mQIEData) {
    return mQIEData->getCoder(fId);
  }
  return nullptr;
}

const CastorQIEShape* CastorDbService::getCastorShape() const {
  if (mQIEData) {
    return &mQIEData->getShape();
  }
  return nullptr;
}

const CastorElectronicsMap* CastorDbService::getCastorMapping() const { return mElectronicsMap; }

const CastorChannelStatus* CastorDbService::getCastorChannelStatus(const HcalGenericDetId& fId) const {
  return mChannelQuality->getValues(fId);
}

TYPELOOKUP_DATA_REG(CastorDbService);