APVGainHistograms

APVmon

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
#ifndef CALIBTRACKER_SISTRIPCHANNELGAIN_APVGAINHELPERS_H
#define CALIBTRACKER_SISTRIPCHANNELGAIN_APVGAINHELPERS_H

#include "CalibTracker/SiStripChannelGain/interface/APVGainStruct.h"
#include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "DQMServices/Core/interface/DQMStore.h"

#include <string>
#include <vector>
#include <utility>
#include <cstdint>
#include <unordered_map>

namespace APVGain {

  typedef dqm::legacy::MonitorElement MonitorElement;

  int subdetectorId(uint32_t);
  int subdetectorId(const std::string&);
  int subdetectorSide(uint32_t, const TrackerTopology*);
  int thickness(uint32_t);
  int thickness(const std::string& tag);
  int subdetectorSide(const std::string&);
  int subdetectorPlane(uint32_t, const TrackerTopology*);
  int subdetectorPlane(const std::string&);

  std::vector<std::pair<std::string, std::string>> monHnames(std::vector<std::string>, bool, const char* tag);

  struct APVmon {
  public:
    APVmon(int v0, int v1, int v2, int v3, MonitorElement* v4)
        : m_thickness(v0), m_subdetectorId(v1), m_subdetectorSide(v2), m_subdetectorPlane(v3), m_monitor(v4) {}

    int getSubdetectorId() { return m_subdetectorId; }

    int getSubdetectorSide() { return m_subdetectorSide; }

    int getSubdetectorPlane() { return m_subdetectorPlane; }

    int getThickness() { return m_thickness; }

    MonitorElement* getMonitor() { return m_monitor; }

    void printAll() {
      LogDebug("APVGainHelpers") << "subDetectorID:" << m_subdetectorId << std::endl;
      LogDebug("APVGainHelpers") << "subDetectorSide:" << m_subdetectorSide << std::endl;
      LogDebug("APVGainHelpers") << "subDetectorPlane:" << m_subdetectorPlane << std::endl;
      LogDebug("APVGainHelpers") << "thickness:" << m_thickness << std::endl;
      LogDebug("APVGainHelpers") << "histoName:" << m_monitor->getName() << std::endl;
      return;
    }

  private:
    int m_thickness;
    int m_subdetectorId;
    int m_subdetectorSide;
    int m_subdetectorPlane;
    MonitorElement* m_monitor;
  };

  struct APVGainHistograms {
  public:
    APVGainHistograms()
        : EventStats(),
          Charge_Vs_Index(7),
          Charge_1(),
          Charge_2(),
          Charge_3(),
          Charge_4(),
          Charge_Vs_PathlengthTIB(7),
          Charge_Vs_PathlengthTOB(7),
          Charge_Vs_PathlengthTIDP(7),
          Charge_Vs_PathlengthTIDM(7),
          Charge_Vs_PathlengthTECP1(7),
          Charge_Vs_PathlengthTECP2(7),
          Charge_Vs_PathlengthTECM1(7),
          Charge_Vs_PathlengthTECM2(7),
          NStripAPVs(0),
          NPixelDets(0),
          APVsCollOrdered(),
          APVsColl() {}

    dqm::reco::MonitorElement* EventStats;
    std::vector<dqm::reco::MonitorElement*> Charge_Vs_Index;         /*!< Charge per cm for each detector id */
    std::array<std::vector<dqm::reco::MonitorElement*>, 7> Charge_1; /*!< Charge per cm per layer / wheel */
    std::array<std::vector<dqm::reco::MonitorElement*>, 7> Charge_2; /*!< Charge per cm per layer / wheel without G2 */
    std::array<std::vector<dqm::reco::MonitorElement*>, 7> Charge_3; /*!< Charge per cm per layer / wheel without G1 */
    std::array<std::vector<dqm::reco::MonitorElement*>, 7>
        Charge_4; /*!< Charge per cm per layer / wheel without G1 and G1*/

    std::vector<dqm::reco::MonitorElement*> Charge_Vs_PathlengthTIB;   /*!< Charge vs pathlength in TIB */
    std::vector<dqm::reco::MonitorElement*> Charge_Vs_PathlengthTOB;   /*!< Charge vs pathlength in TOB */
    std::vector<dqm::reco::MonitorElement*> Charge_Vs_PathlengthTIDP;  /*!< Charge vs pathlength in TIDP */
    std::vector<dqm::reco::MonitorElement*> Charge_Vs_PathlengthTIDM;  /*!< Charge vs pathlength in TIDM */
    std::vector<dqm::reco::MonitorElement*> Charge_Vs_PathlengthTECP1; /*!< Charge vs pathlength in TECP thin */
    std::vector<dqm::reco::MonitorElement*> Charge_Vs_PathlengthTECP2; /*!< Charge vs pathlength in TECP thick */
    std::vector<dqm::reco::MonitorElement*> Charge_Vs_PathlengthTECM1; /*!< Charge vs pathlength in TECP thin */
    std::vector<dqm::reco::MonitorElement*> Charge_Vs_PathlengthTECM2; /*!< Charge vs pathlength in TECP thick */
    mutable std::atomic<unsigned int> NStripAPVs;
    mutable std::atomic<unsigned int> NPixelDets;
    std::vector<std::shared_ptr<stAPVGain>> APVsCollOrdered;
    std::unordered_map<unsigned int, std::shared_ptr<stAPVGain>> APVsColl;
  };

  std::vector<MonitorElement*> FetchMonitor(std::vector<APVmon>, uint32_t, const TrackerTopology* topo = nullptr);
  std::vector<unsigned int> FetchIndices(std::map<unsigned int, APVloc>,
                                         uint32_t,
                                         const TrackerTopology* topo = nullptr);

};  // namespace APVGain

#endif