Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:02:16

0001 #include <algorithm>
0002 #include <iostream>
0003 #include <set>
0004 
0005 #include "CondFormats/HcalObjects/interface/HcalSiPMCharacteristics.h"
0006 #include "CondFormats/HcalObjects/interface/HcalObjectAddons.h"
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008 
0009 HcalSiPMCharacteristics::HcalSiPMCharacteristics(const HcalSiPMCharacteristicsAddons::Helper& helper)
0010     : mPItems(helper.mPItems.begin(), helper.mPItems.end()) {
0011   initialize();
0012 }
0013 
0014 HcalSiPMCharacteristics::~HcalSiPMCharacteristics() {}
0015 
0016 // copy-ctor
0017 HcalSiPMCharacteristics::HcalSiPMCharacteristics(const HcalSiPMCharacteristics& src)
0018     : mPItems(src.mPItems), mPItemsByType(src.mPItemsByType) {}
0019 
0020 // copy assignment operator
0021 HcalSiPMCharacteristics& HcalSiPMCharacteristics::operator=(const HcalSiPMCharacteristics& rhs) {
0022   HcalSiPMCharacteristics temp(rhs);
0023   temp.swap(*this);
0024   return *this;
0025 }
0026 
0027 // public swap function
0028 void HcalSiPMCharacteristics::swap(HcalSiPMCharacteristics& other) {
0029   std::swap(mPItems, other.mPItems);
0030   std::swap(mPItemsByType, other.mPItemsByType);
0031 }
0032 
0033 // move constructor
0034 HcalSiPMCharacteristics::HcalSiPMCharacteristics(HcalSiPMCharacteristics&& other) : HcalSiPMCharacteristics() {
0035   other.swap(*this);
0036 }
0037 
0038 const HcalSiPMCharacteristics::PrecisionItem* HcalSiPMCharacteristics::findByType(int type) const {
0039   const HcalSiPMCharacteristics::PrecisionItem* retItem = nullptr;
0040 
0041   for (unsigned int i = 0; i < getTypes(); i++) {
0042     auto iter = &mPItems.at(i);
0043     if (type == iter->type_)
0044       retItem = iter;
0045   }
0046   return retItem;
0047 
0048   //NOT WORKING:
0049   //PrecisionItem target(type, 0, 0, 0, 0, 0, 0, 0);
0050   //return HcalObjectAddons::findByT<PrecisionItem,HcalSiPMCharacteristicsAddons::LessByType>(&target,mPItemsByType);
0051 }
0052 
0053 HcalSiPMCharacteristicsAddons::Helper::Helper() {}
0054 
0055 bool HcalSiPMCharacteristicsAddons::Helper::loadObject(
0056     int type, int pixels, float parLin1, float parLin2, float parLin3, float crossTalk, int auxi1, float auxi2) {
0057   HcalSiPMCharacteristics::PrecisionItem target(type, pixels, parLin1, parLin2, parLin3, crossTalk, auxi1, auxi2);
0058   auto iter = mPItems.find(target);
0059   if (iter != mPItems.end()) {
0060     edm::LogWarning("HCAL") << "HcalSiPMCharacteristics::loadObject type " << type << " already exists with pixels "
0061                             << iter->pixels_ << " NoLinearity parameters " << iter->parLin1_ << ":" << iter->parLin2_
0062                             << ":" << iter->parLin3_ << " CrossTalk parameter " << iter->crossTalk_ << " new values "
0063                             << pixels << ", " << parLin1 << ", " << parLin2 << ", " << parLin3 << ", " << crossTalk
0064                             << ", " << auxi1 << " and " << auxi2 << " are ignored";
0065     return false;
0066   } else {
0067     mPItems.insert(target);
0068     return true;
0069   }
0070 }
0071 
0072 int HcalSiPMCharacteristics::getPixels(int type) const {
0073   const HcalSiPMCharacteristics::PrecisionItem* item = findByType(type);
0074   return (item ? item->pixels_ : 0);
0075 }
0076 
0077 std::vector<float> HcalSiPMCharacteristics::getNonLinearities(int type) const {
0078   const HcalSiPMCharacteristics::PrecisionItem* item = findByType(type);
0079   std::vector<float> pars;
0080   pars.reserve(3);
0081   if (item) {
0082     pars.push_back(item->parLin1_);
0083     pars.push_back(item->parLin2_);
0084     pars.push_back(item->parLin3_);
0085   }
0086   return pars;
0087 }
0088 
0089 float HcalSiPMCharacteristics::getCrossTalk(int type) const {
0090   const PrecisionItem* item = findByType(type);
0091   return (item ? item->crossTalk_ : 0);
0092 }
0093 
0094 int HcalSiPMCharacteristics::getAuxi1(int type) const {
0095   const HcalSiPMCharacteristics::PrecisionItem* item = findByType(type);
0096   return (item ? item->auxi1_ : 0);
0097 }
0098 
0099 float HcalSiPMCharacteristics::getAuxi2(int type) const {
0100   const HcalSiPMCharacteristics::PrecisionItem* item = findByType(type);
0101   return (item ? item->auxi2_ : 0);
0102 }
0103 
0104 void HcalSiPMCharacteristics::sortByType() {
0105   HcalObjectAddons::sortByT<PrecisionItem, HcalSiPMCharacteristicsAddons::LessByType>(mPItems, mPItemsByType);
0106 }
0107 
0108 void HcalSiPMCharacteristics::initialize() { HcalSiPMCharacteristics::sortByType(); }