Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:51

0001 #ifndef DQMOffline_Trigger_HLTDQMHist_h
0002 #define DQMOffline_Trigger_HLTDQMHist_h
0003 
0004 //********************************************************************************
0005 //
0006 // Description:
0007 //   A MonitorElement together with the necessary information to fill it when pass an
0008 //   object and minimal selection cuts. These selection cuts are intended to be
0009 //   simple selections on kinematic variables. There are two levels of these
0010 //   selection cuts, global (which apply to all MonitorElements in the group and passed
0011 //   by the calling function) and local (which are stored with the MonitorElement
0012 //   and are specific to that MonitorElement). Global selection cuts on the variable being
0013 //   plotted are ignored, for example Et cuts are removed for turn on plots
0014 //
0015 // Implementation:
0016 //   std::function holds the function which generates the vs variable
0017 //   the name of the vs variable is also stored so we can determine if we should not apply
0018 //   a given selection cut
0019 //   selection is done by VarRangeCutColl
0020 //
0021 // Author : Sam Harper , RAL, May 2017
0022 //
0023 //***********************************************************************************
0024 
0025 #include "DQMOffline/Trigger/interface/FunctionDefs.h"
0026 #include "DQMOffline/Trigger/interface/VarRangeCutColl.h"
0027 #include "DQMServices/Core/interface/MonitorElement.h"
0028 #include "FWCore/Framework/interface/Frameworkfwd.h"
0029 
0030 //our base class for our MonitorElements
0031 //takes an object, edm::Event,edm::EventSetup and fills the MonitorElement
0032 //with the predetermined variable (or variables)
0033 template <typename ObjType>
0034 class HLTDQMHist {
0035 public:
0036   HLTDQMHist() = default;
0037   virtual ~HLTDQMHist() = default;
0038   virtual void fill(const ObjType& objType,
0039                     const edm::Event& event,
0040                     const edm::EventSetup& setup,
0041                     const VarRangeCutColl<ObjType>& rangeCuts) = 0;
0042 };
0043 
0044 //this class is a specific implimentation of a HLTDQMHist
0045 //it has the value with which to fill the MonitorElement
0046 //and the MonitorElement itself
0047 //we do not own the MonitorElement
0048 template <typename ObjType, typename ValType>
0049 class HLTDQMHist1D : public HLTDQMHist<ObjType> {
0050 public:
0051   typedef dqm::legacy::MonitorElement MonitorElement;
0052 
0053   HLTDQMHist1D(MonitorElement* me_ptr,
0054                std::string varName,
0055                std::function<ValType(const ObjType&)> func,
0056                VarRangeCutColl<ObjType> rangeCuts)
0057       : var_(std::move(func)), varName_(std::move(varName)), localRangeCuts_(std::move(rangeCuts)), me_ptr_(me_ptr) {}
0058 
0059   void fill(const ObjType& obj,
0060             const edm::Event& event,
0061             const edm::EventSetup& setup,
0062             const VarRangeCutColl<ObjType>& globalRangeCuts) override {
0063     if (me_ptr_ == nullptr)
0064       return;
0065     //local range cuts are specific to a MonitorElement so dont ignore variables
0066     //like global ones (all local cuts should be appropriate to the MonitorElement in question
0067     if (globalRangeCuts(obj, varName_) && localRangeCuts_(obj)) {
0068       me_ptr_->Fill(var_(obj));
0069     }
0070   }
0071 
0072 private:
0073   std::function<ValType(const ObjType&)> var_;
0074   std::string varName_;
0075   VarRangeCutColl<ObjType> localRangeCuts_;
0076   MonitorElement* const me_ptr_;  // we do not own this
0077 };
0078 
0079 template <typename ObjType, typename XValType, typename YValType = XValType>
0080 class HLTDQMHist2D : public HLTDQMHist<ObjType> {
0081 public:
0082   typedef dqm::legacy::MonitorElement MonitorElement;
0083 
0084   HLTDQMHist2D(MonitorElement* me_ptr,
0085                std::string xVarName,
0086                std::string yVarName,
0087                std::function<XValType(const ObjType&)> xFunc,
0088                std::function<YValType(const ObjType&)> yFunc,
0089                VarRangeCutColl<ObjType> rangeCuts)
0090       : xVar_(std::move(xFunc)),
0091         yVar_(std::move(yFunc)),
0092         xVarName_(std::move(xVarName)),
0093         yVarName_(std::move(yVarName)),
0094         localRangeCuts_(std::move(rangeCuts)),
0095         me_ptr_(me_ptr) {}
0096 
0097   void fill(const ObjType& obj,
0098             const edm::Event& event,
0099             const edm::EventSetup& setup,
0100             const VarRangeCutColl<ObjType>& globalRangeCuts) override {
0101     if (me_ptr_ == nullptr)
0102       return;
0103     //local range cuts are specific to a MonitorElement so dont ignore variables
0104     //like global ones (all local cuts should be appropriate to the MonitorElement in question
0105     if (globalRangeCuts(obj, std::vector<std::string>{xVarName_, yVarName_}) && localRangeCuts_(obj)) {
0106       me_ptr_->Fill(xVar_(obj), yVar_(obj));
0107     }
0108   }
0109 
0110 private:
0111   std::function<XValType(const ObjType&)> xVar_;
0112   std::function<YValType(const ObjType&)> yVar_;
0113   std::string xVarName_;
0114   std::string yVarName_;
0115   VarRangeCutColl<ObjType> localRangeCuts_;
0116   MonitorElement* const me_ptr_;  // we do not own this
0117 };
0118 
0119 #endif  // DQMOffline_Trigger_HLTDQMHist_h