CMS 3D CMS Logo

BPHRecoBuilder.h
Go to the documentation of this file.
1 #ifndef HeavyFlavorAnalysis_RecoDecay_BPHRecoBuilder_h
2 #define HeavyFlavorAnalysis_RecoDecay_BPHRecoBuilder_h
3 
14 //----------------------
15 // Base Class Headers --
16 //----------------------
17 
18 //------------------------------------
19 // Collaborating Class Declarations --
20 //------------------------------------
25 class BPHRecoSelect;
26 class BPHMomentumSelect;
27 class BPHVertexSelect;
28 class BPHFitSelect;
29 
30 namespace edm {
31  class EventSetup;
32 }
33 
34 namespace reco {
35  class RecoCandidate;
36 }
37 
38 //---------------
39 // C++ Headers --
40 //---------------
41 #include <string>
42 #include <vector>
43 #include <map>
44 #include <set>
45 
46 // ---------------------
47 // -- Class Interface --
48 // ---------------------
49 
51  friend class BPHRecoSelect;
52 
53 public:
57 
60  virtual ~BPHRecoBuilder();
61 
65  // common object to interface with edm collections
67  public:
68  BPHGenericCollection(const std::string& list) : sList(list) {}
69  virtual ~BPHGenericCollection() {}
70  virtual const reco::Candidate& get(int i) const = 0;
71  virtual int size() const = 0;
72  const std::string& searchList() const { return sList; }
73 
74  private:
76  };
77  template <class T>
78  static BPHGenericCollection* createCollection(const edm::Handle<T>& collection, const std::string& list = "cfhpmig");
79  static BPHGenericCollection* createCollection(const std::vector<const reco::Candidate*>& candList,
80  const std::string& list = "cfhpmig");
81 
90  void add(const std::string& name, const BPHGenericCollection* collection, double mass = -1.0, double msig = -1.0);
91  template <class T>
92  void add(const std::string& name, const edm::Handle<T>& collection, double mass = -1.0, double msig = -1.0);
93  void add(const std::string& name, const std::vector<BPHRecoConstCandPtr>& collection);
94  template <class T>
95  void add(const std::string& name, const std::vector<T>& collection);
96 
99  void filter(const std::string& name, const BPHRecoSelect& sel) const;
101  void filter(const std::string& name, const BPHMomentumSelect& sel) const;
103  void filter(const std::string& name, const BPHVertexSelect& sel) const;
105  void filter(const std::string& name, const BPHFitSelect& sel) const;
106 
109  void filter(const BPHMomentumSelect& sel);
111  void filter(const BPHVertexSelect& sel);
113  void filter(const BPHFitSelect& sel);
114  // apply selection to reconstructed candidate
115  bool accept(const BPHRecoCandidate& cand) const;
116 
119  void setMinPDiffererence(double pMin);
120 
123  struct ComponentSet {
124  std::map<std::string, BPHDecayMomentum::Component> daugMap;
125  std::map<std::string, BPHRecoConstCandPtr> compMap;
126  };
127 
129  std::vector<ComponentSet> build() const;
130 
132  const edm::EventSetup* eventSetup() const;
133 
134  // compare two particles with their track reference and return
135  // true or false for same or different particles, including a
136  // check with momentum difference
137  static bool sameTrack(const reco::Candidate* lCand, const reco::Candidate* rCand, double minPDifference);
138 
139 private:
140  // private copy and assigment constructors
141  BPHRecoBuilder(const BPHRecoBuilder& x) = delete;
142  BPHRecoBuilder& operator=(const BPHRecoBuilder& x) = delete;
143 
144  // object to interface with a specific edm collection
145  template <class T>
147  public:
148  BPHSpecificCollection(const T& c, const std::string& list) : BPHGenericCollection(list), cPtr(&c) {}
150  const reco::Candidate& get(int i) const override { return (*cPtr)[i]; }
151  int size() const override { return cPtr->size(); }
152 
153  private:
154  const T* cPtr;
155  };
156 
157  // object to contain a list of simple particles
158  // with their names, selections, masses and sigma
159  struct BPHRecoSource {
162  std::vector<const BPHRecoSelect*> selector;
163  double mass;
164  double msig;
165  };
166 
167  // object to contain a list of previously reconstructed particles
168  // with their names and selections
169  struct BPHCompSource {
171  const std::vector<BPHRecoConstCandPtr>* collection;
172  std::vector<const BPHMomentumSelect*> momSelector;
173  std::vector<const BPHVertexSelect*> vtxSelector;
174  std::vector<const BPHFitSelect*> fitSelector;
175  };
176 
177  // map of names to simple or previously recontructed particles
178  // for currently tested combination
179  mutable std::map<std::string, const reco::Candidate*> daugMap;
180  mutable std::map<std::string, BPHRecoConstCandPtr> compMap;
181 
183  double minPDiff;
184 
185  // list of simple and previously recontructed particles in the decay
186  std::vector<BPHRecoSource*> sourceList;
187  std::vector<BPHCompSource*> srCompList;
188 
189  // set of copies of previously reconstructed particles list
190  // for bookkeeping and cleanup
191  std::set<const std::vector<BPHRecoConstCandPtr>*> compCollectList;
192 
193  // list fo selections to reconstructed particle
194  std::vector<const BPHMomentumSelect*> msList;
195  std::vector<const BPHVertexSelect*> vsList;
196  std::vector<const BPHFitSelect*> fsList;
197 
198  // map linking particles names to position in list position
199  std::map<std::string, int> sourceId;
200  std::map<std::string, int> srCompId;
201 
202  // recursive function to build particles combinations
203  void build(std::vector<ComponentSet>& compList,
204  ComponentSet& compSet,
205  std::vector<BPHRecoSource*>::const_iterator r_iter,
206  std::vector<BPHRecoSource*>::const_iterator r_iend,
207  std::vector<BPHCompSource*>::const_iterator c_iter,
208  std::vector<BPHCompSource*>::const_iterator c_iend) const;
209 
210  // check for already used particles in a combination
211  // previously recontructed particles are assumed to be included
212  // after simple particles
213  bool contained(ComponentSet& compSet, const reco::Candidate* cand) const;
214  bool contained(ComponentSet& compSet, BPHRecoConstCandPtr cand) const;
215  // compare two particles with their track reference and return
216  // true or false for same or different particles, including a
217  // check with momentum difference
218  bool sameTrack(const reco::Candidate* lCand, const reco::Candidate* rCand) const;
219 };
220 
221 template <class T>
223  const std::string& list) {
225 }
226 
227 template <class T>
228 void BPHRecoBuilder::add(const std::string& name, const edm::Handle<T>& collection, double mass, double msig) {
229  // forward call after creating an interface to the collection
230  add(name, new BPHSpecificCollection<T>(*collection, "cfhpmig"), mass, msig);
231  return;
232 }
233 
234 template <class T>
235 void BPHRecoBuilder::add(const std::string& name, const std::vector<T>& collection) {
236  // forward call after converting the list of pointer to a list
237  // of pointer to base objects
238  int i;
239  int n = collection.size();
240  std::vector<BPHRecoConstCandPtr>* compCandList = new std::vector<BPHRecoConstCandPtr>(n);
241  for (i = 0; i < n; ++i)
242  (*compCandList)[i] = collection[i];
243  // save the converted list for cleanup
244  compCollectList.insert(compCandList);
245  add(name, *compCandList);
246  return;
247 }
248 
249 #endif
std::vector< const BPHRecoSelect * > selector
size
Write out results.
std::map< std::string, BPHRecoConstCandPtr > compMap
BPHGenericPtr< const BPHRecoCandidate >::type BPHRecoConstCandPtr
std::vector< const BPHFitSelect * > fsList
const std::vector< BPHRecoConstCandPtr > * collection
std::map< std::string, BPHDecayMomentum::Component > daugMap
const edm::EventSetup * evSetup
BPHSpecificCollection(const T &c, const std::string &list)
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
BPHGenericCollection(const std::string &list)
std::vector< BPHCompSource * > srCompList
const BPHGenericCollection * collection
static BPHGenericCollection * createCollection(const edm::Handle< T > &collection, const std::string &list="cfhpmig")
std::vector< const BPHVertexSelect * > vsList
std::map< std::string, const reco::Candidate * > daugMap
std::vector< const BPHFitSelect * > fitSelector
std::vector< BPHRecoSource * > sourceList
void add(std::map< std::string, TH1 * > &h, TH1 *hist)
void add(const std::string &name, const BPHGenericCollection *collection, double mass=-1.0, double msig=-1.0)
const std::string & searchList() const
std::vector< const BPHMomentumSelect * > momSelector
std::map< std::string, BPHRecoConstCandPtr > compMap
fixed size matrix
HLT enums.
std::vector< const BPHMomentumSelect * > msList
std::vector< const BPHVertexSelect * > vtxSelector
std::set< const std::vector< BPHRecoConstCandPtr > * > compCollectList
long double T
std::map< std::string, int > srCompId
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger list("!*","!HLTx*"if it matches 2 triggers or more) will accept the event if all the matching triggers are FAIL.It will reject the event if any of the triggers are PASS or EXCEPTION(this matches the behavior of"!*"before the partial wildcard feature was incorporated).Triggers which are in the READY state are completely ignored.(READY should never be returned since the trigger paths have been run
std::map< std::string, int > sourceId