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 //------------------------------------
20 // Collaborating Class Declarations --
21 //------------------------------------
26 class BPHRecoSelect;
27 class BPHMomentumSelect;
28 class BPHVertexSelect;
29 class BPHFitSelect;
30 
31 namespace edm {
32  class EventSetup;
33 }
34 
35 namespace reco {
36  class RecoCandidate;
37 }
38 
39 //---------------
40 // C++ Headers --
41 //---------------
42 #include <string>
43 #include <vector>
44 #include <map>
45 #include <set>
46 
47 // ---------------------
48 // -- Class Interface --
49 // ---------------------
50 
52 
53  friend class BPHRecoSelect;
54 
55  public:
56 
59  BPHRecoBuilder( const edm::EventSetup& es );
60 
63  virtual ~BPHRecoBuilder();
64 
68  // common object to interface with edm collections
70  public:
71  BPHGenericCollection( const std::string& list ): sList( list ) {}
72  virtual ~BPHGenericCollection() {}
73  virtual const reco::Candidate& get( int i ) const = 0;
74  virtual int size() const = 0;
75  const std::string& searchList() const { return sList; }
76  private:
78  };
79  template <class T>
80  static BPHGenericCollection* createCollection(
82  const std::string& list = "cfhpmig" );
83  static BPHGenericCollection* createCollection(
84  const std::vector<const reco::Candidate*>&
85  candList,
86  const std::string& list = "cfhpmig" );
87 
96  void add( const std::string& name,
98  double mass = -1.0,
99  double msig = -1.0 );
100  template <class T>
101  void add( const std::string& name,
102  const edm::Handle<T>& collection,
103  double mass = -1.0,
104  double msig = -1.0 );
105  void add( const std::string& name,
106  const std::vector<BPHRecoConstCandPtr>& collection );
107  template <class T>
108  void add( const std::string& name,
109  const std::vector<T>& collection );
110 
113  void filter( const std::string& name, const BPHRecoSelect & sel ) const;
115  void filter( const std::string& name, const BPHMomentumSelect& sel ) const;
117  void filter( const std::string& name, const BPHVertexSelect & sel ) const;
119  void filter( const std::string& name, const BPHFitSelect & sel ) const;
120 
123  void filter( const BPHMomentumSelect& sel );
125  void filter( const BPHVertexSelect& sel );
127  void filter( const BPHFitSelect& sel );
128  // apply selection to reconstructed candidate
129  bool accept( const BPHRecoCandidate& cand ) const;
130 
133  void setMinPDiffererence( double pMin );
134 
137  struct ComponentSet {
138  std::map<std::string,BPHDecayMomentum::Component> daugMap;
139  std::map<std::string,BPHRecoConstCandPtr > compMap;
140  };
141 
143  std::vector<ComponentSet> build() const;
144 
146  const edm::EventSetup* eventSetup() const;
147 
148  // compare two particles with their track reference and return
149  // true or false for same or different particles, including a
150  // check with momentum difference
151  static
152  bool sameTrack( const reco::Candidate* lCand,
153  const reco::Candidate* rCand, double minPDifference );
154 
155  private:
156 
157  // private copy and assigment constructors
158  BPHRecoBuilder ( const BPHRecoBuilder& x );
159  BPHRecoBuilder& operator=( const BPHRecoBuilder& x );
160 
161  // object to interface with a specific edm collection
162  template <class T>
164  public:
166  BPHGenericCollection( list ), cPtr( &c ) {}
168  virtual const reco::Candidate& get( int i ) const { return (*cPtr)[i]; }
169  virtual int size() const { return cPtr->size(); }
170  private:
171  const T* cPtr;
172  };
173 
174  // object to contain a list of simple particles
175  // with their names, selections, masses and sigma
176  struct BPHRecoSource {
179  std::vector<const BPHRecoSelect*> selector;
180  double mass;
181  double msig;
182  };
183 
184  // object to contain a list of previously reconstructed particles
185  // with their names and selections
186  struct BPHCompSource {
188  const std::vector<BPHRecoConstCandPtr>* collection;
189  std::vector<const BPHMomentumSelect*> momSelector;
190  std::vector<const BPHVertexSelect* > vtxSelector;
191  std::vector<const BPHFitSelect* > fitSelector;
192  };
193 
194  // map of names to simple or previously recontructed particles
195  // for currently tested combination
196  mutable std::map<std::string,const reco::Candidate*> daugMap;
197  mutable std::map<std::string,BPHRecoConstCandPtr > compMap;
198 
200  double minPDiff;
201 
202  // list of simple and previously recontructed particles in the decay
203  std::vector<BPHRecoSource*> sourceList;
204  std::vector<BPHCompSource*> srCompList;
205 
206  // set of copies of previously reconstructed particles list
207  // for bookkeeping and cleanup
208  std::set<const std::vector<BPHRecoConstCandPtr>*> compCollectList;
209 
210  // list fo selections to reconstructed particle
211  std::vector<const BPHMomentumSelect*> msList;
212  std::vector<const BPHVertexSelect *> vsList;
213  std::vector<const BPHFitSelect *> fsList;
214 
215  // map linking particles names to position in list position
216  std::map<std::string,int> sourceId;
217  std::map<std::string,int> srCompId;
218 
219  // recursive function to build particles combinations
220  void build( std::vector<ComponentSet>& compList,
221  ComponentSet& compSet,
222  std::vector<BPHRecoSource*>::const_iterator r_iter,
223  std::vector<BPHRecoSource*>::const_iterator r_iend,
224  std::vector<BPHCompSource*>::const_iterator c_iter,
225  std::vector<BPHCompSource*>::const_iterator c_iend ) const;
226 
227  // check for already used particles in a combination
228  // previously recontructed particles are assumed to be included
229  // after simple particles
230  bool contained( ComponentSet& compSet,
231  const reco::Candidate* cand ) const;
232  bool contained( ComponentSet& compSet,
233  BPHRecoConstCandPtr cand ) const;
234  // compare two particles with their track reference and return
235  // true or false for same or different particles, including a
236  // check with momentum difference
237  bool sameTrack( const reco::Candidate* lCand,
238  const reco::Candidate* rCand ) const;
239 
240 };
241 
242 
243 template <class T>
245  const edm::Handle<T>& collection,
246  const std::string& list ) {
247  return new BPHSpecificCollection<T>( *collection, list );
248 }
249 
250 
251 template <class T>
253  const edm::Handle<T>& collection,
254  double mass,
255  double msig ) {
256  // forward call after creating an interface to the collection
257  add( name, new BPHSpecificCollection<T>( *collection,
258  "cfhpmig" ), mass, msig );
259  return;
260 }
261 
262 
263 template <class T>
265  const std::vector<T>& collection ) {
266  // forward call after converting the list of pointer to a list
267  // of pointer to base objects
268  int i;
269  int n = collection.size();
270  std::vector<BPHRecoConstCandPtr>* compCandList =
271  new std::vector<BPHRecoConstCandPtr>( n );
272  for ( i = 0; i < n; ++i ) (*compCandList)[i] = collection[i];
273  // save the converted list for cleanup
274  compCollectList.insert( compCandList );
275  add( name, *compCandList );
276  return;
277 }
278 
279 
280 #endif
281 
std::vector< const BPHRecoSelect * > selector
size
Write out results.
std::vector< const BPHVertexSelect * > vsList
BPHGenericPtr< const BPHRecoCandidate >::type BPHRecoConstCandPtr
std::map< std::string, BPHRecoConstCandPtr > compMap
const std::vector< BPHRecoConstCandPtr > * collection
std::map< std::string, int > sourceId
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
std::vector< const BPHFitSelect * > fsList
const BPHGenericCollection * collection
static BPHGenericCollection * createCollection(const edm::Handle< T > &collection, const std::string &list="cfhpmig")
std::vector< const BPHFitSelect * > fitSelector
std::vector< const BPHVertexSelect * > vtxSelector
std::map< std::string, const reco::Candidate * > daugMap
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)
std::map< std::string, BPHDecayMomentum::Component > daugMap
const std::string & searchList() const
std::vector< const BPHMomentumSelect * > momSelector
std::map< std::string, int > srCompId
fixed size matrix
HLT enums.
std::vector< const BPHMomentumSelect * > msList
std::set< const std::vector< BPHRecoConstCandPtr > * > compCollectList
long double T
std::map< std::string, BPHRecoConstCandPtr > compMap
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