CMS 3D CMS Logo

IsoDepositVetoFactory.cc
Go to the documentation of this file.
2 
5 #include <regex>
6 
7 // ---------- FIRST DEFINE NEW VETOS ------------
8 namespace reco {
9  namespace isodeposit {
10 
11  class SwitchingEcalVeto : public AbsVeto {
12  public:
13  // creates SwitchingEcalVeto from another AbsVeto (which becomes owned by this veto)
14  SwitchingEcalVeto(AbsVeto *veto, bool isBarrel) : veto_(veto), barrel_(isBarrel) {}
15  bool veto(double eta, double phi, float value) const override {
16  return (fabs(eta) < 1.479) == (barrel_) ? veto_->veto(eta, phi, value) : false;
17  }
18  void centerOn(double eta, double phi) override { veto_->centerOn(eta, phi); }
19 
20  private:
21  std::unique_ptr<AbsVeto> veto_;
22  bool barrel_;
23  };
24 
25  class NumCrystalVeto : public AbsVeto {
26  public:
27  NumCrystalVeto(Direction dir, double iR) : vetoDir_(dir), iR_(iR) {}
28  bool veto(double eta, double phi, float value) const override {
29  if (fabs(vetoDir_.eta()) < 1.479) {
30  return (vetoDir_.deltaR(Direction(eta, phi)) < 0.0174 * iR_);
31  } else {
32  return (vetoDir_.deltaR(Direction(eta, phi)) < 0.00864 * fabs(sinh(eta)) * iR_);
33  }
34  }
35  void centerOn(double eta, double phi) override { vetoDir_ = Direction(eta, phi); }
36 
37  private:
39  float iR_;
40  };
41 
42  class NumCrystalEtaPhiVeto : public AbsVeto {
43  public:
44  NumCrystalEtaPhiVeto(const math::XYZVectorD &dir, double iEta, double iPhi)
45  : vetoDir_(dir.eta(), dir.phi()), iEta_(iEta), iPhi_(iPhi) {}
46  NumCrystalEtaPhiVeto(Direction dir, double iEta, double iPhi)
47  : vetoDir_(dir.eta(), dir.phi()), iEta_(iEta), iPhi_(iPhi) {}
48  bool veto(double eta, double phi, float value) const override {
49  double dPhi = phi - vetoDir_.phi();
50  double dEta = eta - vetoDir_.eta();
51  while (dPhi < -M_PI)
52  dPhi += 2 * M_PI;
53  while (dPhi >= M_PI)
54  dPhi -= 2 * M_PI;
55  if (fabs(vetoDir_.eta()) < 1.479) {
56  return ((fabs(dEta) < 0.0174 * iEta_) && (fabs(dPhi) < 0.0174 * iPhi_));
57  } else {
58  return ((fabs(dEta) < 0.00864 * fabs(sinh(eta)) * iEta_) && (fabs(dPhi) < 0.00864 * fabs(sinh(eta)) * iPhi_));
59  }
60  }
61  void centerOn(double eta, double phi) override { vetoDir_ = Direction(eta, phi); }
62 
63  private:
65  double iEta_, iPhi_;
66  };
67 
68  } // namespace isodeposit
69 } // namespace reco
70 
71 // ---------- THEN THE ACTUAL FACTORY CODE ------------
74  std::unique_ptr<reco::isodeposit::AbsVeto> ret(make(string, evdep, iC));
75  if (evdep != nullptr) {
76  throw cms::Exception("Configuration") << "The resulting AbsVeto depends on the edm::Event.\n"
77  << "Please use the two-arguments IsoDepositVetoFactory::make.\n";
78  }
79  return ret.release();
80 }
81 
85  using namespace reco::isodeposit;
86  static const std::regex ecalSwitch("^Ecal(Barrel|Endcaps):(.*)"), threshold("Threshold\\((\\d+\\.\\d+)\\)"),
87  thresholdtransverse("ThresholdFromTransverse\\((\\d+\\.\\d+)\\)"),
88  absthreshold("AbsThreshold\\((\\d+\\.\\d+)\\)"),
89  absthresholdtransverse("AbsThresholdFromTransverse\\((\\d+\\.\\d+)\\)"), cone("ConeVeto\\((\\d+\\.\\d+)\\)"),
90  angleCone("AngleCone\\((\\d+\\.\\d+)\\)"), angleVeto("AngleVeto\\((\\d+\\.\\d+)\\)"),
91  rectangularEtaPhiVeto(
92  "RectangularEtaPhiVeto\\(([+-]?\\d+\\.\\d+),([+-]?\\d+\\.\\d+),([+-]?\\d+\\.\\d+),([+-]?\\d+\\.\\d+)\\)"),
93  numCrystal("NumCrystalVeto\\((\\d+\\.\\d+)\\)"),
94  numCrystalEtaPhi("NumCrystalEtaPhiVeto\\((\\d+\\.\\d+),(\\d+\\.\\d+)\\)"),
95  otherCandidatesDR("OtherCandidatesByDR\\((\\w+:?\\w*:?\\w*),\\s*(\\d+\\.?|\\d*\\.\\d*)\\)"),
96  otherJetConstituentsDR(
97  "OtherJetConstituentsDeltaRVeto\\((\\w+:?\\w*:?\\w*),\\s*(\\d+\\.?|\\d*\\.\\d*),\\s*(\\w+:?\\w*:?\\w*),\\s*("
98  "\\d+\\.?|\\d*\\.\\d*)\\)"),
99  otherCand("^(.*?):(.*)"), number("^(\\d+\\.?|\\d*\\.\\d*)$");
100  std::cmatch match;
101 
102  evdep = nullptr; // by default it does not depend on this
103  if (regex_match(string, match, ecalSwitch)) {
104  return new SwitchingEcalVeto(make(match[2].first, iC), (match[1] == "Barrel"));
105  } else if (regex_match(string, match, threshold)) {
106  return new ThresholdVeto(atof(match[1].first));
107  } else if (regex_match(string, match, thresholdtransverse)) {
108  return new ThresholdVetoFromTransverse(atof(((std::string)match[1]).c_str()));
109  } else if (regex_match(string, match, absthreshold)) {
110  return new AbsThresholdVeto(atof(match[1].first));
111  } else if (regex_match(string, match, absthresholdtransverse)) {
112  return new AbsThresholdVetoFromTransverse(atof(((std::string)match[1]).c_str()));
113  } else if (regex_match(string, match, cone)) {
114  return new ConeVeto(Direction(), atof(match[1].first));
115  } else if (regex_match(string, match, number)) {
116  return new ConeVeto(Direction(), atof(match[1].first));
117  } else if (regex_match(string, match, angleCone)) {
118  return new AngleCone(Direction(), atof(match[1].first));
119  } else if (regex_match(string, match, angleVeto)) {
120  return new AngleConeVeto(Direction(), atof(match[1].first));
121  } else if (regex_match(string, match, rectangularEtaPhiVeto)) {
122  return new RectangularEtaPhiVeto(
123  Direction(), atof(match[1].first), atof(match[2].first), atof(match[3].first), atof(match[4].first));
124  } else if (regex_match(string, match, numCrystal)) {
125  return new NumCrystalVeto(Direction(), atof(match[1].first));
126  } else if (regex_match(string, match, numCrystalEtaPhi)) {
127  return new NumCrystalEtaPhiVeto(Direction(), atof(match[1].first), atof(match[2].first));
128  } else if (regex_match(string, match, otherCandidatesDR)) {
129  OtherCandidatesDeltaRVeto *ret = new OtherCandidatesDeltaRVeto(edm::InputTag(match[1]), atof(match[2].first), iC);
130  evdep = ret;
131  return ret;
132  } else if (regex_match(string, match, otherJetConstituentsDR)) {
134  Direction(), edm::InputTag(match[1]), atof(match[2].first), edm::InputTag(match[3]), atof(match[4].first), iC);
135  evdep = ret;
136  return ret;
137  } else if (regex_match(string, match, otherCand)) {
138  OtherCandVeto *ret = new OtherCandVeto(edm::InputTag(match[1]), make(match[2].first, iC), iC);
139  evdep = ret;
140  return ret;
141  } else {
142  throw cms::Exception("Not Implemented") << "Veto " << string << " not implemented yet...";
143  }
144 }
static reco::isodeposit::AbsVeto * make(const char *string, edm::ConsumesCollector &iC)
void centerOn(double eta, double phi) override
void centerOn(double eta, double phi) override
ret
prodAgent to be discontinued
S make(const edm::ParameterSet &cfg)
NumCrystalEtaPhiVeto(const math::XYZVectorD &dir, double iEta, double iPhi)
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
bool veto(double eta, double phi, float value) const override
Return "true" if a deposit at specific (eta,phi) with that value must be vetoed in the sum...
Definition: value.py:1
#define M_PI
NumCrystalEtaPhiVeto(Direction dir, double iEta, double iPhi)
void centerOn(double eta, double phi) override
bool veto(double eta, double phi, float value) const override
Return "true" if a deposit at specific (eta,phi) with that value must be vetoed in the sum...
fixed size matrix
bool veto(double eta, double phi, float value) const override
Return "true" if a deposit at specific (eta,phi) with that value must be vetoed in the sum...
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
NumCrystalVeto(Direction dir, double iR)
SwitchingEcalVeto(AbsVeto *veto, bool isBarrel)