CMS 3D CMS Logo

BPHWriteSpecificDecay.h
Go to the documentation of this file.
1 #ifndef HeavyFlavorAnalysis_SpecificDecay_BPHWriteSpecificDecay_h
2 #define HeavyFlavorAnalysis_SpecificDecay_BPHWriteSpecificDecay_h
3 
5 
10 
14 
17 
23 
26 
32 
33 #include <string>
34 #include <vector>
35 #include <map>
36 
37 class TH1F;
38 class BPHRecoCandidate;
39 
40 class BPHWriteSpecificDecay : public BPHAnalyzerWrapper<BPHModuleWrapper::stream_producer> {
41 public:
42  explicit BPHWriteSpecificDecay(const edm::ParameterSet& ps);
43  ~BPHWriteSpecificDecay() override = default;
44 
45  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
46 
47  void produce(edm::Event& ev, const edm::EventSetup& es) override;
48  virtual void fill(edm::Event& ev, const BPHEventSetupWrapper& es);
49 
50 private:
61 
62  // token wrappers to allow running both on "old" and "new" CMSSW versions
75 
76  bool usePV;
77  bool usePM;
78  bool useCC;
79  bool usePF;
80  bool usePC;
81  bool useGP;
82  bool useK0;
83  bool useL0;
84  bool useKS;
85  bool useLS;
86 
101 
102  enum recoType {
113  Bu,
114  Bp,
115  Bd,
116  Bs,
119  B0,
121  Bc,
124  };
125  enum parType {
149  };
150  std::map<std::string, recoType> rMap;
151  std::map<std::string, parType> pMap;
152  std::map<std::string, parType> fMap;
153  std::map<recoType, std::map<parType, double> > parMap;
154 
155  bool recoOnia;
156  bool recoKx0;
157  bool recoPkk;
158  bool recoBu;
159  bool recoBp;
160  bool recoBd;
161  bool recoBs;
162  bool recoK0s;
164  bool recoB0;
166  bool recoBc;
167  bool recoPsi2S;
168  bool recoX3872;
169 
170  bool allKx0;
171  bool allPkk;
172  bool allK0s;
174 
175  bool writeOnia;
176  bool writeKx0;
177  bool writePkk;
178  bool writeBu;
179  bool writeBp;
180  bool writeBd;
181  bool writeBs;
182  bool writeK0s;
184  bool writeB0;
186  bool writeBc;
189 
192 
193  std::vector<BPHPlusMinusConstCandPtr> lFull;
194  std::vector<BPHPlusMinusConstCandPtr> lJPsi;
195  std::vector<BPHRecoConstCandPtr> lSd;
196  std::vector<BPHRecoConstCandPtr> lSs;
197  std::vector<BPHRecoConstCandPtr> lBu;
198  std::vector<BPHRecoConstCandPtr> lBp;
199  std::vector<BPHRecoConstCandPtr> lBd;
200  std::vector<BPHRecoConstCandPtr> lBs;
201  std::vector<BPHPlusMinusConstCandPtr> lK0;
202  std::vector<BPHPlusMinusConstCandPtr> lL0;
203  std::vector<BPHRecoConstCandPtr> lB0;
204  std::vector<BPHRecoConstCandPtr> lLb;
205  std::vector<BPHRecoConstCandPtr> lBc;
206  std::vector<BPHRecoConstCandPtr> lPsi2S;
207  std::vector<BPHRecoConstCandPtr> lX3872;
208 
209  std::map<const BPHRecoCandidate*, const BPHRecoCandidate*> jPsiOMap;
210  std::map<const BPHRecoCandidate*, const BPHRecoCandidate*> daughMap;
212  std::map<const BPHRecoCandidate*, vertex_ref> pvRefMap;
214  std::map<const BPHRecoCandidate*, compcc_ref> ccRefMap;
215 
216  void setRecoParameters(const edm::ParameterSet& ps);
217 
218  static void addTrackModes(const std::string& name, const BPHRecoCandidate& cand, std::string& modes, bool& count);
219 
221 
222  template <class T>
224  const std::vector<T>& list,
225  const std::string& name) {
227  int i;
228  int n = list.size();
229  std::map<const BPHRecoCandidate*, const BPHRecoCandidate*>::const_iterator dauIter;
230  std::map<const BPHRecoCandidate*, const BPHRecoCandidate*>::const_iterator dauIend = daughMap.end();
231  std::map<const BPHRecoCandidate*, const BPHRecoCandidate*>::const_iterator jpoIter;
232  std::map<const BPHRecoCandidate*, const BPHRecoCandidate*>::const_iterator jpoIend = jPsiOMap.end();
233  std::map<const BPHRecoCandidate*, vertex_ref>::const_iterator pvrIter;
234  std::map<const BPHRecoCandidate*, vertex_ref>::const_iterator pvrIend = pvRefMap.end();
235  std::map<const BPHRecoCandidate*, compcc_ref>::const_iterator ccrIter;
236  std::map<const BPHRecoCandidate*, compcc_ref>::const_iterator ccrIend = ccRefMap.end();
237  for (i = 0; i < n; ++i) {
238  const T& ptr = list[i];
239  ccList->push_back(ptr->composite());
240  pat::CompositeCandidate& cc = ccList->back();
241  std::string modes;
242  bool count = false;
243  addTrackModes("", *ptr, modes, count);
244  cc.addUserData("trackModes", modes, true);
245  addTrackModes("trackMode_", *ptr, cc);
246  if ((pvrIter = pvRefMap.find(ptr.get())) != pvrIend)
247  cc.addUserData("primaryVertex", pvrIter->second);
248  const std::vector<std::string>& cNames = ptr->compNames();
249  int j = 0;
250  int m = cNames.size();
251  while (j < m) {
252  const std::string& compName = cNames[j++];
253  const BPHRecoCandidate* cptr = ptr->getComp(compName).get();
254  if ((ccrIter = ccRefMap.find(cptr)) == ccrIend) {
255  if ((dauIter = daughMap.find(cptr)) != dauIend)
256  cptr = dauIter->second;
257  if ((jpoIter = jPsiOMap.find(cptr)) != jpoIend)
258  cptr = jpoIter->second;
259  }
260  if ((ccrIter = ccRefMap.find(cptr)) != ccrIend) {
261  compcc_ref cref = ccrIter->second;
262  if (cref.isNonnull())
263  cc.addUserData("refTo" + compName, cref);
264  }
265  }
266  const BPHPlusMinusCandidate* pmp = dynamic_cast<const BPHPlusMinusCandidate*>(ptr.get());
267  if (pmp != nullptr) {
268  cc.addUserInt("cowboy", (pmp->isCowboy() ? +1 : -1));
269  // cc.addUserFloat( "dca", pmp->cAppInRPhi().distance() );
270  }
271  if (writeVertex)
272  cc.addUserData("vertex", ptr->vertex());
273  if (ptr->isEmpty())
274  continue;
275  if (writeVertex)
276  cc.addUserData("fitVertex", reco::Vertex(*ptr->topDecayVertex()));
277  if (ptr->isValidFit()) {
278  const RefCountedKinematicParticle kinPart = ptr->topParticle();
279  const KinematicState kinStat = kinPart->currentState();
280  cc.addUserFloat("fitMass", kinStat.mass());
281  if (writeMomentum)
282  cc.addUserData("fitMomentum", kinStat.kinematicParameters().momentum());
283  }
284  }
285  typedef std::unique_ptr<pat::CompositeCandidateCollection> ccc_pointer;
286  edm::OrphanHandle<pat::CompositeCandidateCollection> ccHandle = ev.put(ccc_pointer(ccList), name);
287  for (i = 0; i < n; ++i) {
288  const BPHRecoCandidate* ptr = list[i].get();
290  ccRefMap[ptr] = ccRef;
291  }
292  return ccHandle;
293  }
294 };
295 
296 #endif
BPHTokenWrapper< std::vector< reco::VertexCompositeCandidate > > l0CandsToken
Analysis-level particle class.
BPHTokenWrapper< std::vector< reco::Vertex > > pVertexToken
std::vector< BPHPlusMinusConstCandPtr > lJPsi
BPHTokenWrapper< std::vector< reco::VertexCompositePtrCandidate > > kSCandsToken
std::vector< BPHPlusMinusConstCandPtr > lK0
edm::Ref< std::vector< reco::Vertex > > vertex_ref
std::map< const BPHRecoCandidate *, compcc_ref > ccRefMap
std::vector< BPHPlusMinusConstCandPtr > lFull
std::map< std::string, recoType > rMap
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
edm::Ref< pat::CompositeCandidateCollection > compcc_ref
std::vector< BPHRecoConstCandPtr > lBs
std::vector< BPHRecoConstCandPtr > lBp
BPHTokenWrapper< std::vector< BPHTrackReference::candidate > > pcCandsToken
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
BPHESTokenWrapper< TransientTrackBuilder, TransientTrackRecord > ttBToken
BPHESTokenWrapper< MagneticField, IdealMagneticFieldRecord > magFieldToken
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool isCowboy() const
get cowboy/sailor classification
std::map< const BPHRecoCandidate *, const BPHRecoCandidate * > daughMap
std::map< std::string, parType > pMap
std::vector< BPHRecoConstCandPtr > lSs
BPHTokenWrapper< std::vector< reco::VertexCompositeCandidate > > k0CandsToken
BPHWriteSpecificDecay(const edm::ParameterSet &ps)
BPHTokenWrapper< pat::MuonCollection > patMuonToken
GlobalVector momentum() const
BPHTokenWrapper< std::vector< reco::PFCandidate > > pfCandsToken
KinematicParameters const & kinematicParameters() const
BPHTokenWrapper< std::vector< reco::VertexCompositePtrCandidate > > lSCandsToken
BPHTokenWrapper< std::vector< pat::CompositeCandidate > > ccCandsToken
std::vector< BPHRecoConstCandPtr > lSd
~BPHWriteSpecificDecay() override=default
virtual BPHRecoConstCandPtr getComp(const std::string &name) const
BPHTokenWrapper< std::vector< pat::GenericParticle > > gpCandsToken
std::vector< CompositeCandidate > CompositeCandidateCollection
edm::OrphanHandle< pat::CompositeCandidateCollection > write(edm::Event &ev, const std::vector< T > &list, const std::string &name)
virtual void fill(edm::Event &ev, const BPHEventSetupWrapper &es)
std::map< const BPHRecoCandidate *, vertex_ref > pvRefMap
std::vector< BPHPlusMinusConstCandPtr > lL0
std::vector< BPHRecoConstCandPtr > lBu
std::vector< BPHRecoConstCandPtr > lLb
std::map< recoType, std::map< parType, double > > parMap
std::vector< BPHRecoConstCandPtr > lPsi2S
ParticleMass mass() const
std::map< const BPHRecoCandidate *, const BPHRecoCandidate * > jPsiOMap
std::vector< BPHRecoConstCandPtr > lB0
std::vector< BPHRecoConstCandPtr > lBc
std::vector< BPHRecoConstCandPtr > lX3872
long double T
void produce(edm::Event &ev, const edm::EventSetup &es) override
static void addTrackModes(const std::string &name, const BPHRecoCandidate &cand, std::string &modes, bool &count)
std::vector< BPHRecoConstCandPtr > lBd
void setRecoParameters(const edm::ParameterSet &ps)
std::map< std::string, parType > fMap