SiPixelVCalDB

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
// system includes
#include <fstream>
#include <iostream>
#include <map>
#include <memory>
#include <string>

// user includes
#include "CondFormats/SiPixelObjects/interface/SiPixelVCal.h"
#include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
#include "DataFormats/DetId/interface/DetId.h"
#include "DataFormats/TrackerCommon/interface/PixelBarrelName.h"
#include "DataFormats/TrackerCommon/interface/PixelEndcapName.h"
#include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ServiceRegistry/interface/Service.h"
#include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
#include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
#include "Geometry/Records/interface/TrackerTopologyRcd.h"
#include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"

class SiPixelVCalDB : public edm::one::EDAnalyzer<> {
public:
  explicit SiPixelVCalDB(const edm::ParameterSet& conf);
  explicit SiPixelVCalDB();
  ~SiPixelVCalDB() override;
  void analyze(const edm::Event&, const edm::EventSetup&) override;

private:
  const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> tkGeomToken_;
  const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> tkTopoToken_;
  std::string recordName_;
  typedef std::vector<edm::ParameterSet> Parameters;
  Parameters BPixParameters_;
  Parameters FPixParameters_;
};

using namespace std;
using namespace edm;

SiPixelVCalDB::SiPixelVCalDB(edm::ParameterSet const& iConfig)
    : tkGeomToken_(esConsumes()), tkTopoToken_(esConsumes()) {
  recordName_ = iConfig.getUntrackedParameter<std::string>("record", "SiPixelVCalRcd");
  BPixParameters_ = iConfig.getUntrackedParameter<Parameters>("BPixParameters");
  FPixParameters_ = iConfig.getUntrackedParameter<Parameters>("FPixParameters");
}

SiPixelVCalDB::~SiPixelVCalDB() = default;

// Analyzer: Functions that gets called by framework every event
void SiPixelVCalDB::analyze(const edm::Event& e, const edm::EventSetup& iSetup) {
  SiPixelVCal vcal;
  bool phase1 = true;

  // Retrieve tracker topology from geometry
  const TrackerTopology* const tTopo = &iSetup.getData(tkTopoToken_);

  // Retrieve old style tracker geometry from geometry
  const TrackerGeometry* pDD = &iSetup.getData(tkGeomToken_);
  edm::LogPrint("SiPixelVCalDB") << " There are " << pDD->detUnits().size() << " modules" << std::endl;

  for (const auto& it : pDD->detUnits()) {
    if (dynamic_cast<PixelGeomDetUnit const*>(it) != nullptr) {
      const DetId detid = it->geographicalId();
      const unsigned int rawDetId = detid.rawId();
      int subid = detid.subdetId();

      // FILL BPIX
      if (subid == static_cast<int>(PixelSubdetector::PixelBarrel)) {
        int layer = tTopo->pxbLayer(detid);    // 1, 2, 3, 4
        int ladder = tTopo->pxbLadder(detid);  // 1-12/28/44/64
        edm::LogPrint("SiPixelVCalDB") << " pixel barrel:"
                                       << " detId=" << rawDetId << ", layer=" << layer << ", ladder=" << ladder;
        for (Parameters::iterator it = BPixParameters_.begin(); it != BPixParameters_.end(); ++it) {
          if (it->getParameter<int>("layer") == layer && it->getParameter<int>("ladder") == ladder) {
            float slope = (float)it->getParameter<double>("slope");
            float offset = (float)it->getParameter<double>("offset");
            edm::LogPrint("SiPixelVCalDB") << ";  VCal slope " << slope << ", offset " << offset;
            // edm::LogInfo("SiPixelVCalDB")  << "  detId " << rawDetId << " \t
            // VCal slope " << slope << ", offset " << offset;
            vcal.putSlopeAndOffset(detid, slope, offset);
          }
        }
        edm::LogPrint("SiPixelVCalDB") << std::endl;

        // FILL FPIX
      } else if (subid == static_cast<int>(PixelSubdetector::PixelEndcap)) {
        PixelEndcapName fpix(detid, tTopo, phase1);
        int side = tTopo->pxfSide(detid);   // 1 (-z), 2 for (+z)
        int disk = fpix.diskName();         // 1, 2, 3
        int disk2 = tTopo->pxfDisk(detid);  // 1, 2, 3
        int ring = fpix.ringName();         // 1 (lower), 2 (upper)
        if (disk != disk2) {
          edm::LogError("SiPixelVCalDB::analyze")
              << "Found contradicting FPIX disk number: " << disk << " vs." << disk2 << std::endl;
        }
        edm::LogPrint("SiPixelVCalDB") << " pixel endcap:"
                                       << " detId=" << rawDetId << ", side=" << side << ", disk=" << disk
                                       << ", ring=" << ring;
        for (Parameters::iterator it = FPixParameters_.begin(); it != FPixParameters_.end(); ++it) {
          if (it->getParameter<int>("side") == side && it->getParameter<int>("disk") == disk &&
              it->getParameter<int>("ring") == ring) {
            float slope = (float)it->getParameter<double>("slope");
            float offset = (float)it->getParameter<double>("offset");
            edm::LogPrint("SiPixelVCalDB") << ";  VCal slope " << slope << ", offset " << offset;
            // edm::LogInfo("SiPixelVCalDB")  << "  detId " << rawDetId << " \t
            // VCal slope " << slope << ", offset " << offset;
            vcal.putSlopeAndOffset(rawDetId, slope, offset);
          }
        }
        edm::LogPrint("SiPixelVCalDB") << std::endl;

      } else {
        edm::LogError("SiPixelVCalDB::analyze") << "detid is Pixel but neither bpix nor fpix" << std::endl;
      }
    }
  }

  // Save to DB
  edm::Service<cond::service::PoolDBOutputService> mydbservice;
  if (mydbservice.isAvailable()) {
    try {
      if (mydbservice->isNewTagRequest(recordName_)) {
        mydbservice->createOneIOV<SiPixelVCal>(vcal, mydbservice->beginOfTime(), recordName_);
      } else {
        mydbservice->appendOneIOV<SiPixelVCal>(vcal, mydbservice->currentTime(), recordName_);
      }
    } catch (const cond::Exception& er) {
      edm::LogError("SiPixelVCalDB") << er.what() << std::endl;
    } catch (const std::exception& er) {
      edm::LogError("SiPixelVCalDB") << "caught std::exception " << er.what() << std::endl;
    } catch (...) {
      edm::LogError("SiPixelVCalDB") << "Funny error" << std::endl;
    }
  } else {
    edm::LogError("SiPixelVCalDB") << "Service is unavailable" << std::endl;
  }
}
DEFINE_FWK_MODULE(SiPixelVCalDB);