CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Static Private Member Functions | Private Attributes
BPHWriteSpecificDecay Class Reference

#include <BPHWriteSpecificDecay.h>

Inheritance diagram for BPHWriteSpecificDecay:
BPHAnalyzerWrapper< BPHModuleWrapper::stream_producer > edm::stream::EDProducer< T >

Public Member Functions

 BPHWriteSpecificDecay (const edm::ParameterSet &ps)
 
virtual void fill (edm::Event &ev, const BPHEventSetupWrapper &es)
 
void produce (edm::Event &ev, const edm::EventSetup &es) override
 
 ~BPHWriteSpecificDecay () override=default
 
- Public Member Functions inherited from edm::stream::EDProducer< T >
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Types

typedef edm::Ref< pat::CompositeCandidateCollectioncompcc_ref
 
enum  parType {
  ptMin, etaMax, mPsiMin, mPsiMax,
  mKx0Min, mKx0Max, mPhiMin, mPhiMax,
  mK0sMin, mK0sMax, mLambda0Min, mLambda0Max,
  massMin, massMax, probMin, mFitMin,
  mFitMax, constrMass, constrSigma, requireJPsi,
  constrMJPsi, constrMPsi2, writeCandidate
}
 
enum  recoType {
  Onia, Pmm, Psi1, Psi2,
  Ups, Ups1, Ups2, Ups3,
  Kx0, Pkk, Bu, Bp,
  Bd, Bs, K0s, Lambda0,
  B0, Lambdab, Bc, Psi2S,
  X3872
}
 
typedef edm::Ref< std::vector< reco::Vertex > > vertex_ref
 

Private Member Functions

void setRecoParameters (const edm::ParameterSet &ps)
 
template<class T >
edm::OrphanHandle< pat::CompositeCandidateCollectionwrite (edm::Event &ev, const std::vector< T > &list, const std::string &name)
 

Static Private Member Functions

static void addTrackModes (const std::string &name, const BPHRecoCandidate &cand, std::string &modes, bool &count)
 
static void addTrackModes (const std::string &name, const BPHRecoCandidate &cand, pat::CompositeCandidate &cc)
 

Private Attributes

bool allK0s
 
bool allKx0
 
bool allLambda0
 
bool allPkk
 
std::string b0Name
 
std::string bcName
 
std::string bdName
 
std::string bpName
 
std::string bsName
 
std::string buName
 
std::string ccCandsLabel
 
BPHTokenWrapper< std::vector< pat::CompositeCandidate > > ccCandsToken
 
std::map< const BPHRecoCandidate *, compcc_refccRefMap
 
std::map< const BPHRecoCandidate *, const BPHRecoCandidate * > daughMap
 
std::map< std::string, parTypefMap
 
std::string gpCandsLabel
 
BPHTokenWrapper< std::vector< pat::GenericParticle > > gpCandsToken
 
std::map< const BPHRecoCandidate *, const BPHRecoCandidate * > jPsiOMap
 
std::string k0CandsLabel
 
BPHTokenWrapper< std::vector< reco::VertexCompositeCandidate > > k0CandsToken
 
std::string k0Name
 
std::string kSCandsLabel
 
BPHTokenWrapper< std::vector< reco::VertexCompositePtrCandidate > > kSCandsToken
 
std::string l0CandsLabel
 
BPHTokenWrapper< std::vector< reco::VertexCompositeCandidate > > l0CandsToken
 
std::string l0Name
 
std::vector< BPHRecoConstCandPtrlB0
 
std::vector< BPHRecoConstCandPtrlBc
 
std::vector< BPHRecoConstCandPtrlBd
 
std::string lbName
 
std::vector< BPHRecoConstCandPtrlBp
 
std::vector< BPHRecoConstCandPtrlBs
 
std::vector< BPHRecoConstCandPtrlBu
 
std::vector< BPHPlusMinusConstCandPtrlFull
 
std::vector< BPHPlusMinusConstCandPtrlJPsi
 
std::vector< BPHPlusMinusConstCandPtrlK0
 
std::vector< BPHPlusMinusConstCandPtrlL0
 
std::vector< BPHRecoConstCandPtrlLb
 
std::vector< BPHRecoConstCandPtrlPsi2S
 
std::string lSCandsLabel
 
BPHTokenWrapper< std::vector< reco::VertexCompositePtrCandidate > > lSCandsToken
 
std::vector< BPHRecoConstCandPtrlSd
 
std::vector< BPHRecoConstCandPtrlSs
 
std::vector< BPHRecoConstCandPtrlX3872
 
BPHESTokenWrapper< MagneticField, IdealMagneticFieldRecordmagFieldToken
 
std::string oniaName
 
std::map< recoType, std::map< parType, double > > parMap
 
std::string patMuonLabel
 
BPHTokenWrapper< pat::MuonCollectionpatMuonToken
 
std::string pcCandsLabel
 
BPHTokenWrapper< std::vector< BPHTrackReference::candidate > > pcCandsToken
 
std::string pfCandsLabel
 
BPHTokenWrapper< std::vector< reco::PFCandidate > > pfCandsToken
 
std::map< std::string, parTypepMap
 
std::string psi2SName
 
std::string pVertexLabel
 
BPHTokenWrapper< std::vector< reco::Vertex > > pVertexToken
 
std::map< const BPHRecoCandidate *, vertex_refpvRefMap
 
bool recoB0
 
bool recoBc
 
bool recoBd
 
bool recoBp
 
bool recoBs
 
bool recoBu
 
bool recoK0s
 
bool recoKx0
 
bool recoLambda0
 
bool recoLambdab
 
bool recoOnia
 
bool recoPkk
 
bool recoPsi2S
 
bool recoX3872
 
std::map< std::string, recoTyperMap
 
std::string sdName
 
std::string ssName
 
BPHESTokenWrapper< TransientTrackBuilder, TransientTrackRecordttBToken
 
bool useCC
 
bool useGP
 
bool useK0
 
bool useKS
 
bool useL0
 
bool useLS
 
bool usePC
 
bool usePF
 
bool usePM
 
bool usePV
 
bool writeB0
 
bool writeBc
 
bool writeBd
 
bool writeBp
 
bool writeBs
 
bool writeBu
 
bool writeK0s
 
bool writeKx0
 
bool writeLambda0
 
bool writeLambdab
 
bool writeMomentum
 
bool writeOnia
 
bool writePkk
 
bool writePsi2S
 
bool writeVertex
 
bool writeX3872
 
std::string x3872Name
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer< T >
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 
- Protected Member Functions inherited from BPHAnalyzerWrapper< BPHModuleWrapper::stream_producer >
void consume (BPHTokenWrapper< Obj > &tw, const std::string &label)
 
void consume (BPHTokenWrapper< Obj > &tw, const edm::InputTag &tag)
 
void esConsume (BPHESTokenWrapper< Obj, Rec > &tw)
 
void esConsume (BPHESTokenWrapper< Obj, Rec > &tw, const std::string &label)
 
void esConsume (BPHESTokenWrapper< Obj, Rec > &tw, const edm::ESInputTag &tag)
 

Detailed Description

Definition at line 40 of file BPHWriteSpecificDecay.h.

Member Typedef Documentation

◆ compcc_ref

Definition at line 213 of file BPHWriteSpecificDecay.h.

◆ vertex_ref

Definition at line 211 of file BPHWriteSpecificDecay.h.

Member Enumeration Documentation

◆ parType

Enumerator
ptMin 
etaMax 
mPsiMin 
mPsiMax 
mKx0Min 
mKx0Max 
mPhiMin 
mPhiMax 
mK0sMin 
mK0sMax 
mLambda0Min 
mLambda0Max 
massMin 
massMax 
probMin 
mFitMin 
mFitMax 
constrMass 
constrSigma 
requireJPsi 
constrMJPsi 
constrMPsi2 
writeCandidate 

Definition at line 125 of file BPHWriteSpecificDecay.h.

125  {
126  ptMin,
127  etaMax,
128  mPsiMin,
129  mPsiMax,
130  mKx0Min,
131  mKx0Max,
132  mPhiMin,
133  mPhiMax,
134  mK0sMin,
135  mK0sMax,
136  mLambda0Min,
137  mLambda0Max,
138  massMin,
139  massMax,
140  probMin,
141  mFitMin,
142  mFitMax,
143  constrMass,
144  constrSigma,
145  requireJPsi,
146  constrMJPsi,
147  constrMPsi2,
149  };

◆ recoType

Enumerator
Onia 
Pmm 
Psi1 
Psi2 
Ups 
Ups1 
Ups2 
Ups3 
Kx0 
Pkk 
Bu 
Bp 
Bd 
Bs 
K0s 
Lambda0 
B0 
Lambdab 
Bc 
Psi2S 
X3872 

Definition at line 102 of file BPHWriteSpecificDecay.h.

102  {
103  Onia,
104  Pmm,
105  Psi1,
106  Psi2,
107  Ups,
108  Ups1,
109  Ups2,
110  Ups3,
111  Kx0,
112  Pkk,
113  Bu,
114  Bp,
115  Bd,
116  Bs,
117  K0s,
118  Lambda0,
119  B0,
120  Lambdab,
121  Bc,
122  Psi2S,
123  X3872
124  };

Constructor & Destructor Documentation

◆ BPHWriteSpecificDecay()

BPHWriteSpecificDecay::BPHWriteSpecificDecay ( const edm::ParameterSet ps)
explicit

Definition at line 60 of file BPHWriteSpecificDecay.cc.

References heavyFlavorDQMFirstStep_cff::b0Name, heavyFlavorDQMFirstStep_cff::bcName, heavyFlavorDQMFirstStep_cff::bdName, heavyFlavorDQMFirstStep_cff::bpName, heavyFlavorDQMFirstStep_cff::bsName, heavyFlavorDQMFirstStep_cff::buName, vertexSelectForHeavyFlavorDQM_cfi::constrMass, vertexSelectForHeavyFlavorDQM_cfi::constrMJPsi, vertexSelectForHeavyFlavorDQM_cfi::constrMPsi2, vertexSelectForHeavyFlavorDQM_cfi::constrSigma, ALCARECOTkAlBeamHalo_cff::etaMax, edm::ParameterSet::getParameter(), heavyFlavorDQMFirstStep_cff::k0CandsLabel, heavyFlavorDQMFirstStep_cff::k0Name, heavyFlavorDQMFirstStep_cff::l0CandsLabel, heavyFlavorDQMFirstStep_cff::l0Name, heavyFlavorDQMFirstStep_cff::lbName, vertexSelectForHeavyFlavorDQM_cfi::massMax, vertexSelectForHeavyFlavorDQM_cfi::massMin, vertexSelectForHeavyFlavorDQM_cfi::mK0sMax, vertexSelectForHeavyFlavorDQM_cfi::mK0sMin, vertexSelectForHeavyFlavorDQM_cfi::mKx0Max, vertexSelectForHeavyFlavorDQM_cfi::mKx0Min, vertexSelectForHeavyFlavorDQM_cfi::mLambda0Max, vertexSelectForHeavyFlavorDQM_cfi::mLambda0Min, vertexSelectForHeavyFlavorDQM_cfi::mPhiMax, vertexSelectForHeavyFlavorDQM_cfi::mPhiMin, heavyFlavorDQMFirstStep_cff::oniaName, heavyFlavorDQMFirstStep_cff::patMuonLabel, heavyFlavorDQMFirstStep_cff::pfCandsLabel, vertexSelectForHeavyFlavorDQM_cfi::probMin, heavyFlavorDQMFirstStep_cff::psi2SName, ptMin, heavyFlavorDQMFirstStep_cff::pVertexLabel, heavyFlavorDQMFirstStep_cff::recoSelect, vertexSelectForHeavyFlavorDQM_cfi::requireJPsi, heavyFlavorDQMFirstStep_cff::sdName, SET_PAR, heavyFlavorDQMFirstStep_cff::ssName, trackingRecoMaterialAnalyzer_cfi::usePV, heavyFlavorDQMFirstStep_cff::writeMomentum, and heavyFlavorDQMFirstStep_cff::writeVertex.

60  {
61  usePV = (!SET_PAR(string, pVertexLabel, ps).empty());
62  usePM = (!SET_PAR(string, patMuonLabel, ps).empty());
63  useCC = (!SET_PAR(string, ccCandsLabel, ps).empty());
64  usePF = (!SET_PAR(string, pfCandsLabel, ps).empty());
65  usePC = (!SET_PAR(string, pcCandsLabel, ps).empty());
66  useGP = (!SET_PAR(string, gpCandsLabel, ps).empty());
67  useK0 = (!SET_PAR(string, k0CandsLabel, ps).empty());
68  useL0 = (!SET_PAR(string, l0CandsLabel, ps).empty());
69  useKS = (!SET_PAR(string, kSCandsLabel, ps).empty());
70  useLS = (!SET_PAR(string, lSCandsLabel, ps).empty());
71  SET_PAR(string, oniaName, ps);
72  SET_PAR(string, sdName, ps);
73  SET_PAR(string, ssName, ps);
74  SET_PAR(string, buName, ps);
75  SET_PAR(string, bpName, ps);
76  SET_PAR(string, bdName, ps);
77  SET_PAR(string, bsName, ps);
78  SET_PAR(string, k0Name, ps);
79  SET_PAR(string, l0Name, ps);
80  SET_PAR(string, b0Name, ps);
81  SET_PAR(string, lbName, ps);
82  SET_PAR(string, bcName, ps);
83  SET_PAR(string, psi2SName, ps);
84  SET_PAR(string, x3872Name, ps);
85 
86  SET_PAR(bool, writeMomentum, ps);
87  SET_PAR(bool, writeVertex, ps);
88 
89  rMap["Onia"] = Onia;
90  rMap["PhiMuMu"] = Pmm;
91  rMap["Psi1"] = Psi1;
92  rMap["Psi2"] = Psi2;
93  rMap["Ups"] = Ups;
94  rMap["Ups1"] = Ups1;
95  rMap["Ups2"] = Ups2;
96  rMap["Ups3"] = Ups3;
97  rMap["Kx0"] = Kx0;
98  rMap["PhiKK"] = Pkk;
99  rMap["Bu"] = Bu;
100  rMap["Bp"] = Bp;
101  rMap["Bd"] = Bd;
102  rMap["Bs"] = Bs;
103  rMap["K0s"] = K0s;
104  rMap["Lambda0"] = Lambda0;
105  rMap["B0"] = B0;
106  rMap["Lambdab"] = Lambdab;
107  rMap["Bc"] = Bc;
108  rMap["Psi2S"] = Psi2S;
109  rMap["X3872"] = X3872;
110 
111  pMap["ptMin"] = ptMin;
112  pMap["etaMax"] = etaMax;
113  pMap["mJPsiMin"] = mPsiMin;
114  pMap["mJPsiMax"] = mPsiMax;
115  pMap["mKx0Min"] = mKx0Min;
116  pMap["mKx0Max"] = mKx0Max;
117  pMap["mPhiMin"] = mPhiMin;
118  pMap["mPhiMax"] = mPhiMax;
119  pMap["mK0sMin"] = mK0sMin;
120  pMap["mK0sMax"] = mK0sMax;
121  pMap["mLambda0Min"] = mLambda0Min;
122  pMap["mLambda0Max"] = mLambda0Max;
123  pMap["massMin"] = massMin;
124  pMap["massMax"] = massMax;
125  pMap["probMin"] = probMin;
126  pMap["massFitMin"] = mFitMin;
127  pMap["massFitMax"] = mFitMax;
128  pMap["constrMass"] = constrMass;
129  pMap["constrSigma"] = constrSigma;
130 
131  fMap["requireJPsi"] = requireJPsi;
132  fMap["constrMJPsi"] = constrMJPsi;
133  fMap["constrMPsi2"] = constrMPsi2;
134 
135  fMap["writeCandidate"] = writeCandidate;
136 
140 
141  writeOnia = true;
142  const vector<edm::ParameterSet> recoSelect = ps.getParameter<vector<edm::ParameterSet> >("recoSelect");
143  int iSel;
144  int nSel = recoSelect.size();
145  for (iSel = 0; iSel < nSel; ++iSel)
147  if (!recoOnia)
148  writeOnia = false;
149 
150  if (recoBu)
151  recoOnia = true;
152  if (recoBd)
153  recoOnia = recoKx0 = true;
154  if (recoBs)
155  recoOnia = recoPkk = true;
156  if (recoB0)
157  recoOnia = recoK0s = true;
158  if (recoLambdab)
159  recoOnia = recoLambda0 = true;
160  if (recoBc)
161  recoOnia = true;
162  if (recoPsi2S)
163  recoOnia = true;
164  if (recoX3872)
165  recoOnia = true;
166  if (writeBu)
167  writeOnia = true;
168  if (writeBd)
169  writeOnia = writeKx0 = true;
170  if (writeBs)
171  writeOnia = writePkk = true;
172  if (writeB0)
173  writeOnia = writeK0s = true;
174  if (writeLambdab)
175  writeOnia = writeLambda0 = true;
176  if (writeBc)
177  writeOnia = true;
178  if (writePsi2S)
179  writeOnia = true;
180  if (writeX3872)
181  writeOnia = true;
182  if (recoBp && !recoPsi2S && !recoX3872)
183  recoPsi2S = true;
184  if (writeBp && !writePsi2S && !writeX3872)
185  writePsi2S = true;
186  allKx0 = (parMap[Kx0][requireJPsi] < 0);
187  allPkk = (parMap[Pkk][requireJPsi] < 0);
188  allK0s = (parMap[K0s][requireJPsi] < 0);
190 
191  esConsume<MagneticField, IdealMagneticFieldRecord>(magFieldToken);
192  esConsume<TransientTrackBuilder, TransientTrackRecord>(ttBToken, "TransientTrackBuilder");
193  if (usePV)
194  consume<vector<reco::Vertex> >(pVertexToken, pVertexLabel);
195  if (usePM)
196  consume<pat::MuonCollection>(patMuonToken, patMuonLabel);
197  if (useCC)
198  consume<vector<pat::CompositeCandidate> >(ccCandsToken, ccCandsLabel);
199  if (usePF)
200  consume<vector<reco::PFCandidate> >(pfCandsToken, pfCandsLabel);
201  if (usePC)
202  consume<vector<BPHTrackReference::candidate> >(pcCandsToken, pcCandsLabel);
203  if (useGP)
204  consume<vector<pat::GenericParticle> >(gpCandsToken, gpCandsLabel);
205  if (useK0)
206  consume<vector<reco::VertexCompositeCandidate> >(k0CandsToken, k0CandsLabel);
207  if (useL0)
208  consume<vector<reco::VertexCompositeCandidate> >(l0CandsToken, l0CandsLabel);
209  if (useKS)
210  consume<vector<reco::VertexCompositePtrCandidate> >(kSCandsToken, kSCandsLabel);
211  if (useLS)
212  consume<vector<reco::VertexCompositePtrCandidate> >(lSCandsToken, lSCandsLabel);
213 
214  if (writeOnia)
215  produces<pat::CompositeCandidateCollection>(oniaName);
216  if (writeKx0)
217  produces<pat::CompositeCandidateCollection>(sdName);
218  if (writePkk)
219  produces<pat::CompositeCandidateCollection>(ssName);
220  if (writeBu)
221  produces<pat::CompositeCandidateCollection>(buName);
222  if (writeBp)
223  produces<pat::CompositeCandidateCollection>(bpName);
224  if (writeBd)
225  produces<pat::CompositeCandidateCollection>(bdName);
226  if (writeBs)
227  produces<pat::CompositeCandidateCollection>(bsName);
228  if (writeK0s)
229  produces<pat::CompositeCandidateCollection>(k0Name);
230  if (writeLambda0)
231  produces<pat::CompositeCandidateCollection>(l0Name);
232  if (writeB0)
233  produces<pat::CompositeCandidateCollection>(b0Name);
234  if (writeLambdab)
235  produces<pat::CompositeCandidateCollection>(lbName);
236  if (writeBc)
237  produces<pat::CompositeCandidateCollection>(bcName);
238  if (writePsi2S)
239  produces<pat::CompositeCandidateCollection>(psi2SName);
240  if (writeX3872)
241  produces<pat::CompositeCandidateCollection>(x3872Name);
242 }
BPHTokenWrapper< std::vector< reco::VertexCompositeCandidate > > l0CandsToken
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
BPHTokenWrapper< std::vector< reco::Vertex > > pVertexToken
BPHTokenWrapper< std::vector< reco::VertexCompositePtrCandidate > > kSCandsToken
std::map< std::string, recoType > rMap
BPHTokenWrapper< std::vector< BPHTrackReference::candidate > > pcCandsToken
BPHESTokenWrapper< TransientTrackBuilder, TransientTrackRecord > ttBToken
BPHESTokenWrapper< MagneticField, IdealMagneticFieldRecord > magFieldToken
#define SET_PAR(TYPE, NAME, PSET)
std::map< std::string, parType > pMap
BPHTokenWrapper< std::vector< reco::VertexCompositeCandidate > > k0CandsToken
BPHTokenWrapper< pat::MuonCollection > patMuonToken
BPHTokenWrapper< std::vector< reco::PFCandidate > > pfCandsToken
BPHTokenWrapper< std::vector< reco::VertexCompositePtrCandidate > > lSCandsToken
BPHTokenWrapper< std::vector< pat::CompositeCandidate > > ccCandsToken
BPHTokenWrapper< std::vector< pat::GenericParticle > > gpCandsToken
std::map< recoType, std::map< parType, double > > parMap
void setRecoParameters(const edm::ParameterSet &ps)
std::map< std::string, parType > fMap

◆ ~BPHWriteSpecificDecay()

BPHWriteSpecificDecay::~BPHWriteSpecificDecay ( )
overridedefault

Member Function Documentation

◆ addTrackModes() [1/2]

void BPHWriteSpecificDecay::addTrackModes ( const std::string &  name,
const BPHRecoCandidate cand,
std::string &  modes,
bool &  count 
)
staticprivate

Definition at line 1617 of file BPHWriteSpecificDecay.cc.

References submitPVResolutionJobs::count, mps_splice::entry, and Skims_PA_cff::name.

Referenced by write().

1620  {
1622  if (count)
1623  modes += "#";
1624  modes += (name + entry.first + ":" + cand.getTMode(entry.second));
1625  count = true;
1626  }
1628  addTrackModes(entry.first + "/", *entry.second, modes, count);
1629  }
1630  return;
1631 }
static void addTrackModes(const std::string &name, const BPHRecoCandidate &cand, std::string &modes, bool &count)

◆ addTrackModes() [2/2]

void BPHWriteSpecificDecay::addTrackModes ( const std::string &  name,
const BPHRecoCandidate cand,
pat::CompositeCandidate cc 
)
staticprivate

Definition at line 1633 of file BPHWriteSpecificDecay.cc.

References gpuPixelDoublets::cc, mps_splice::entry, and Skims_PA_cff::name.

1635  {
1637  cc.addUserData(name + entry.first, string(1, cand.getTMode(entry.second)), true);
1639  addTrackModes(name + entry.first + "/", *entry.second, cc);
1640  return;
1641 }
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
static void addTrackModes(const std::string &name, const BPHRecoCandidate &cand, std::string &modes, bool &count)

◆ fill()

void BPHWriteSpecificDecay::fill ( edm::Event ev,
const BPHEventSetupWrapper es 
)
virtual

Definition at line 337 of file BPHWriteSpecificDecay.cc.

References accept(), BPHMassFitSelect::accept(), PVValHelper::add(), b0, reco::Candidate::begin(), cms::cuda::bs, BPHDecayGenericBuilder< ProdType >::build(), TwoTrackMinimumDistance::calculate(), gpuPixelDoublets::cc, edm::HandleBase::clear(), BPHPlusMinusCandidate::composite(), vertexSelectForHeavyFlavorDQM_cfi::constrMass, vertexSelectForHeavyFlavorDQM_cfi::constrMJPsi, vertexSelectForHeavyFlavorDQM_cfi::constrMPsi2, vertexSelectForHeavyFlavorDQM_cfi::constrSigma, BPHRecoBuilder::createCollection(), BPHDecayToFlyingCascadeBuilderBase::daughMap(), BPHBuToPsi2SKBuilder::daughMap(), BPHDecayMomentum::daughters(), l1tTrackerHTMiss_cfi::deltaZ, Calorimetry_cff::dp, MillePedeFileConverter_cfg::e, ALCARECOTkAlBeamHalo_cff::etaMax, makeMEIFBenchmarkPlots::ev, BPHEventSetupWrapper::get(), BPHOniaToMuMuBuilder::getConstrMass(), BPHOniaToMuMuBuilder::getConstrSigma(), BPHDecayMomentum::getDaug(), BPHOniaToMuMuBuilder::getList(), BPHDecayGenericBuilderBase::getMassFitMax(), BPHDecayGenericBuilderBase::getMassFitMin(), mps_fire::i, dqmiolumiharvest::j, BPHParticleMasses::jPsiMass, BPHParticleMasses::jPsiMWidth, MainPageGenerator::l, visualization-live-secondInstance_cfg::m, HLT_2023v12_cff::magneticField, vertexSelectForHeavyFlavorDQM_cfi::massMax, vertexSelectForHeavyFlavorDQM_cfi::massMin, pfMETCorrectionType0_cfi::minDz, vertexSelectForHeavyFlavorDQM_cfi::mK0sMax, vertexSelectForHeavyFlavorDQM_cfi::mK0sMin, vertexSelectForHeavyFlavorDQM_cfi::mKx0Max, vertexSelectForHeavyFlavorDQM_cfi::mKx0Min, vertexSelectForHeavyFlavorDQM_cfi::mLambda0Max, vertexSelectForHeavyFlavorDQM_cfi::mLambda0Min, vertexSelectForHeavyFlavorDQM_cfi::mPhiMax, vertexSelectForHeavyFlavorDQM_cfi::mPhiMin, dqmiodumpmetadata::n, BPHDecayMomentum::originalReco(), BPHOniaToMuMuBuilder::Phi, TwoTrackMinimumDistance::points(), reco::Vertex::position(), funct::pow(), vertexSelectForHeavyFlavorDQM_cfi::probMin, edm::Handle< T >::product(), BPHOniaToMuMuBuilder::Psi1, BPHOniaToMuMuBuilder::Psi2, BPHParticleMasses::psi2Mass, BPHParticleMasses::psi2MWidth, ptMin, MetAnalyzer::pv(), BPHRecoBuilder::sameTrack(), edm::second(), BPHDecayConstrainedBuilderBase::setConstr(), BPHOniaToMuMuBuilder::setConstr(), BPHDecayToChargedXXbarBuilder::setEtaMax(), BPHOniaToMuMuBuilder::setEtaMax(), BPHKx0ToKPiBuilder::setEtaMax(), BPHDecayToV0DiffMassBuilder::setEtaMax(), BPHBdToJPsiKxBuilder::setJPsiMassMax(), BPHLbToJPsiL0Builder::setJPsiMassMax(), BPHBuToJPsiKBuilder::setJPsiMassMax(), BPHDecayToJPsiPiPiBuilder::setJPsiMassMax(), BPHBcToJPsiPiBuilder::setJPsiMassMax(), BPHBdToJPsiKxBuilder::setJPsiMassMin(), BPHLbToJPsiL0Builder::setJPsiMassMin(), BPHBuToJPsiKBuilder::setJPsiMassMin(), BPHBcToJPsiPiBuilder::setJPsiMassMin(), BPHDecayToJPsiPiPiBuilder::setJPsiMassMin(), BPHBuToJPsiKBuilder::setKEtaMax(), BPHBuToPsi2SKBuilder::setKEtaMax(), BPHBuToJPsiKBuilder::setKPtMin(), BPHBuToPsi2SKBuilder::setKPtMin(), BPHBdToJPsiKxBuilder::setKxMassMax(), BPHBdToJPsiKxBuilder::setKxMassMin(), BPHLbToJPsiL0Builder::setLambda0MassMax(), BPHLbToJPsiL0Builder::setLambda0MassMin(), BPHDecayGenericBuilderBase::setMassFitMax(), BPHDecayGenericBuilderBase::setMassFitMin(), BPHDecayConstrainedBuilderBase::setMassFitSelect(), BPHDecayGenericBuilderBase::setMassMax(), BPHOniaToMuMuBuilder::setMassMax(), BPHDecayGenericBuilderBase::setMassMin(), BPHOniaToMuMuBuilder::setMassMin(), BPHBcToJPsiPiBuilder::setPiEtaMax(), BPHDecayToJPsiPiPiBuilder::setPiEtaMax(), BPHDecayToJPsiPiPiBuilder::setPiPtMin(), BPHBcToJPsiPiBuilder::setPiPtMin(), BPHDecayGenericBuilderBase::setProbMin(), BPHOniaToMuMuBuilder::setProbMin(), BPHBuToPsi2SKBuilder::setPsi2SMassMax(), BPHBuToPsi2SKBuilder::setPsi2SMassMin(), BPHDecayToChargedXXbarBuilder::setPtMin(), BPHKx0ToKPiBuilder::setPtMin(), BPHOniaToMuMuBuilder::setPtMin(), BPHDecayToV0DiffMassBuilder::setPtMin(), jetUpdater_cfi::sort, BPHOniaToMuMuBuilder::Ups, BPHOniaToMuMuBuilder::Ups1, BPHOniaToMuMuBuilder::Ups2, BPHOniaToMuMuBuilder::Ups3, pat::PATObject< ObjectType >::userData(), and BPHDecayVertex::vertex().

337  {
338  lFull.clear();
339  lJPsi.clear();
340  lSd.clear();
341  lSs.clear();
342  lBu.clear();
343  lBp.clear();
344  lBd.clear();
345  lBs.clear();
346  lK0.clear();
347  lL0.clear();
348  lB0.clear();
349  lLb.clear();
350  lBc.clear();
351  lPsi2S.clear();
352  lX3872.clear();
353  jPsiOMap.clear();
354  daughMap.clear();
355  pvRefMap.clear();
356  ccRefMap.clear();
357 
358  // get magnetic field
359  // data are got through "BPHESTokenWrapper" interface to allow
360  // uniform access in different CMSSW versions
363 
364  // get object collections
365  // collections are got through "BPHTokenWrapper" interface to allow
366  // uniform access in different CMSSW versions
367 
369  pVertexToken.get(ev, pVertices);
370  int npv = pVertices->size();
371 
372  int nrc = 0;
373 
374  // get reco::PFCandidate collection (in full AOD )
376  if (usePF) {
377  pfCandsToken.get(ev, pfCands);
378  nrc = pfCands->size();
379  }
380 
381  // get pat::PackedCandidate collection (in MiniAOD)
382  // pat::PackedCandidate is not defined in CMSSW_5XY, so a
383  // typedef (BPHTrackReference::candidate) is used, actually referring
384  // to pat::PackedCandidate only for CMSSW versions where it's defined
386  if (usePC) {
387  pcCandsToken.get(ev, pcCands);
388  nrc = pcCands->size();
389  }
390 
391  // get pat::GenericParticle collection (in skimmed data)
393  if (useGP) {
394  gpCandsToken.get(ev, gpCands);
395  nrc = gpCands->size();
396  }
397 
398  // get pat::Muon collection (in full AOD and MiniAOD)
400  if (usePM) {
401  patMuonToken.get(ev, patMuon);
402  }
403 
404  // get K0 reco::VertexCompositeCandidate collection (in full AOD)
406  if (useK0) {
407  k0CandsToken.get(ev, k0Cand);
408  }
409 
410  // get Lambda0 reco::VertexCompositeCandidate collection (in full AOD)
412  if (useL0) {
413  l0CandsToken.get(ev, l0Cand);
414  }
415 
416  // get K0 reco::VertexCompositePtrCandidate collection (in MiniAOD)
418  if (useKS) {
419  kSCandsToken.get(ev, kSCand);
420  }
421 
422  // get Lambda0 reco::VertexCompositePtrCandidate collection (in MiniAOD)
424  if (useLS) {
425  lSCandsToken.get(ev, lSCand);
426  }
427 
428  // get muons from pat::CompositeCandidate objects describing onia;
429  // muons from all composite objects are copied to an unique std::vector
430  vector<const reco::Candidate*> muDaugs;
431  set<const pat::Muon*> muonSet;
432  typedef multimap<const reco::Candidate*, const pat::CompositeCandidate*> mu_cc_map;
433  mu_cc_map muCCMap;
434  if (useCC) {
436  ccCandsToken.get(ev, ccCands);
437  int n = ccCands->size();
438  muDaugs.clear();
439  muDaugs.reserve(n);
440  muonSet.clear();
441  set<const pat::Muon*>::const_iterator iter;
442  set<const pat::Muon*>::const_iterator iend;
443  int i;
444  for (i = 0; i < n; ++i) {
445  const pat::CompositeCandidate& cc = ccCands->at(i);
446  int j;
447  int m = cc.numberOfDaughters();
448  for (j = 0; j < m; ++j) {
449  const reco::Candidate* dp = cc.daughter(j);
450  const pat::Muon* mp = dynamic_cast<const pat::Muon*>(dp);
451  iter = muonSet.begin();
452  iend = muonSet.end();
453  bool add = (mp != nullptr) && (muonSet.find(mp) == iend);
454  while (add && (iter != iend)) {
455  if (BPHRecoBuilder::sameTrack(mp, *iter++, 1.0e-5))
456  add = false;
457  }
458  if (add)
459  muonSet.insert(mp);
460  // associate muon to the CompositeCandidate containing it
461  muCCMap.insert(pair<const reco::Candidate*, const pat::CompositeCandidate*>(dp, &cc));
462  }
463  }
464  iter = muonSet.begin();
465  iend = muonSet.end();
466  while (iter != iend)
467  muDaugs.push_back(*iter++);
468  }
469 
470  map<recoType, map<parType, double> >::const_iterator rIter = parMap.begin();
471  map<recoType, map<parType, double> >::const_iterator rIend = parMap.end();
472 
473  // reconstruct quarkonia
474 
475  BPHOniaToMuMuBuilder* onia = nullptr;
476  if (recoOnia) {
477  if (usePM)
478  onia = new BPHOniaToMuMuBuilder(
479  es, BPHRecoBuilder::createCollection(patMuon, "ingmcf"), BPHRecoBuilder::createCollection(patMuon, "ingmcf"));
480  else if (useCC)
481  onia = new BPHOniaToMuMuBuilder(
482  es, BPHRecoBuilder::createCollection(muDaugs, "ingmcf"), BPHRecoBuilder::createCollection(muDaugs, "ingmcf"));
483  }
484 
485  if (onia != nullptr) {
486  while (rIter != rIend) {
487  const map<recoType, map<parType, double> >::value_type& rEntry = *rIter++;
488  recoType rType = rEntry.first;
489  const map<parType, double>& pMap = rEntry.second;
491  switch (rType) {
492  case Pmm:
494  break;
495  case Psi1:
497  break;
498  case Psi2:
500  break;
501  case Ups:
503  break;
504  case Ups1:
506  break;
507  case Ups2:
509  break;
510  case Ups3:
512  break;
513  default:
514  continue;
515  }
516  map<parType, double>::const_iterator pIter = pMap.begin();
517  map<parType, double>::const_iterator pIend = pMap.end();
518  while (pIter != pIend) {
519  const map<parType, double>::value_type& pEntry = *pIter++;
520  parType id = pEntry.first;
521  double pv = pEntry.second;
522  switch (id) {
523  case ptMin:
524  onia->setPtMin(type, pv);
525  break;
526  case etaMax:
527  onia->setEtaMax(type, pv);
528  break;
529  case massMin:
530  onia->setMassMin(type, pv);
531  break;
532  case massMax:
533  onia->setMassMax(type, pv);
534  break;
535  case probMin:
536  onia->setProbMin(type, pv);
537  break;
538  case constrMass:
539  onia->setConstr(type, pv, onia->getConstrSigma(type));
540  break;
541  case constrSigma:
542  onia->setConstr(type, onia->getConstrMass(type), pv);
543  break;
544  default:
545  break;
546  }
547  }
548  }
549  lFull = onia->build();
550  }
551 
552  // associate onia to primary vertex
553 
554  int iFull;
555  int nFull = lFull.size();
556  map<const BPHRecoCandidate*, const reco::Vertex*> oniaVtxMap;
557 
558  typedef mu_cc_map::const_iterator mu_cc_iter;
559  for (iFull = 0; iFull < nFull; ++iFull) {
560  const reco::Vertex* pVtx = nullptr;
561  int pvId = 0;
562  const BPHPlusMinusCandidate* ptr = lFull[iFull].get();
563  const std::vector<const reco::Candidate*>& daugs = ptr->daughters();
564 
565  // try to recover primary vertex association in skim data:
566  // get the CompositeCandidate containing both muons
567  pair<mu_cc_iter, mu_cc_iter> cc0 = muCCMap.equal_range(ptr->originalReco(daugs[0]));
568  pair<mu_cc_iter, mu_cc_iter> cc1 = muCCMap.equal_range(ptr->originalReco(daugs[1]));
569  mu_cc_iter iter0 = cc0.first;
570  mu_cc_iter iend0 = cc0.second;
571  mu_cc_iter iter1 = cc1.first;
572  mu_cc_iter iend1 = cc1.second;
573  while ((iter0 != iend0) && (pVtx == nullptr)) {
574  const pat::CompositeCandidate* ccp = iter0++->second;
575  while (iter1 != iend1) {
576  if (ccp != iter1++->second)
577  continue;
578  pVtx = ccp->userData<reco::Vertex>("PVwithmuons");
579  const reco::Vertex* sVtx = nullptr;
580  const reco::Vertex::Point& pPos = pVtx->position();
581  float dMin = 999999.;
582  int ipv;
583  for (ipv = 0; ipv < npv; ++ipv) {
584  const reco::Vertex* tVtx = &pVertices->at(ipv);
585  const reco::Vertex::Point& tPos = tVtx->position();
586  float dist = pow(pPos.x() - tPos.x(), 2) + pow(pPos.y() - tPos.y(), 2) + pow(pPos.z() - tPos.z(), 2);
587  if (dist < dMin) {
588  dMin = dist;
589  sVtx = tVtx;
590  pvId = ipv;
591  }
592  }
593  pVtx = sVtx;
594  break;
595  }
596  }
597 
598  // if not found, as for other type of input data,
599  // try to get the nearest primary vertex in z direction
600  if (pVtx == nullptr) {
601  const reco::Vertex::Point& sVtp = ptr->vertex().position();
602  GlobalPoint cPos(sVtp.x(), sVtp.y(), sVtp.z());
603  const pat::CompositeCandidate& sCC = ptr->composite();
604  GlobalVector cDir(sCC.px(), sCC.py(), sCC.pz());
605  GlobalPoint bPos(0.0, 0.0, 0.0);
606  GlobalVector bDir(0.0, 0.0, 1.0);
608  bool state = ttmd.calculate(GlobalTrajectoryParameters(cPos, cDir, TrackCharge(0), &(*magneticField)),
610  float minDz = 999999.;
611  float extrapZ = (state ? ttmd.points().first.z() : -9e20);
612  int ipv;
613  for (ipv = 0; ipv < npv; ++ipv) {
614  const reco::Vertex& tVtx = pVertices->at(ipv);
615  float deltaZ = fabs(extrapZ - tVtx.position().z());
616  if (deltaZ < minDz) {
617  minDz = deltaZ;
618  pVtx = &tVtx;
619  pvId = ipv;
620  }
621  }
622  }
623 
624  oniaVtxMap[ptr] = pVtx;
625  pvRefMap[ptr] = vertex_ref(pVertices, pvId);
626  }
627  pVertexToken.get(ev, pVertices);
628 
629  // get JPsi subsample and associate JPsi candidate to original
630  // generic onia candidate
631  if (nFull)
633 
634  bool jPsiFound = !lJPsi.empty();
635  delete onia;
636 
637  if (!nrc)
638  return;
639 
640  int ij;
641  int io;
642  int nj = lJPsi.size();
643  int no = lFull.size();
644  for (ij = 0; ij < nj; ++ij) {
645  const BPHRecoCandidate* jp = lJPsi[ij].get();
646  for (io = 0; io < no; ++io) {
647  const BPHRecoCandidate* oc = lFull[io].get();
648  if ((jp->originalReco(jp->getDaug("MuPos")) == oc->originalReco(oc->getDaug("MuPos"))) &&
649  (jp->originalReco(jp->getDaug("MuNeg")) == oc->originalReco(oc->getDaug("MuNeg")))) {
650  jPsiOMap[jp] = oc;
651  break;
652  }
653  }
654  }
655 
656  // build and dump Bu
657 
658  BPHBuToJPsiKBuilder* bu = nullptr;
659  if (recoBu && jPsiFound) {
660  if (usePF)
661  bu = new BPHBuToJPsiKBuilder(es, lJPsi, BPHRecoBuilder::createCollection(pfCands, "f"));
662  else if (usePC)
663  bu = new BPHBuToJPsiKBuilder(es, lJPsi, BPHRecoBuilder::createCollection(pcCands, "p"));
664  else if (useGP)
665  bu = new BPHBuToJPsiKBuilder(es, lJPsi, BPHRecoBuilder::createCollection(gpCands, "h"));
666  }
667 
668  if (bu != nullptr) {
669  rIter = parMap.find(Bu);
670  if (rIter != rIend) {
671  const map<parType, double>& pMap = rIter->second;
672  map<parType, double>::const_iterator pIter = pMap.begin();
673  map<parType, double>::const_iterator pIend = pMap.end();
674  while (pIter != pIend) {
675  const map<parType, double>::value_type& pEntry = *pIter++;
676  parType id = pEntry.first;
677  double pv = pEntry.second;
678  switch (id) {
679  case ptMin:
680  bu->setKPtMin(pv);
681  break;
682  case etaMax:
683  bu->setKEtaMax(pv);
684  break;
685  case mPsiMin:
686  bu->setJPsiMassMin(pv);
687  break;
688  case mPsiMax:
689  bu->setJPsiMassMax(pv);
690  break;
691  case massMin:
692  bu->setMassMin(pv);
693  break;
694  case massMax:
695  bu->setMassMax(pv);
696  break;
697  case probMin:
698  bu->setProbMin(pv);
699  break;
700  case mFitMin:
701  bu->setMassFitMin(pv);
702  break;
703  case mFitMax:
704  bu->setMassFitMax(pv);
705  break;
706  case constrMJPsi:
707  bu->setConstr(pv > 0);
708  break;
709  case writeCandidate:
710  writeBu = (pv > 0);
711  break;
712  default:
713  break;
714  }
715  }
716  }
717  lBu = bu->build();
718  delete bu;
719  }
720 
721  // build and dump Kx0
722 
723  vector<BPHPlusMinusConstCandPtr> lKx0;
724  BPHKx0ToKPiBuilder* kx0 = nullptr;
725  if (recoKx0 && (jPsiFound || allKx0)) {
726  if (usePF)
727  kx0 = new BPHKx0ToKPiBuilder(
728  es, BPHRecoBuilder::createCollection(pfCands, "f"), BPHRecoBuilder::createCollection(pfCands, "f"));
729  else if (usePC)
730  kx0 = new BPHKx0ToKPiBuilder(
731  es, BPHRecoBuilder::createCollection(pcCands, "p"), BPHRecoBuilder::createCollection(pcCands, "p"));
732  else if (useGP)
733  kx0 = new BPHKx0ToKPiBuilder(
734  es, BPHRecoBuilder::createCollection(gpCands, "h"), BPHRecoBuilder::createCollection(gpCands, "h"));
735  }
736 
737  set<BPHRecoConstCandPtr> sKx0;
738 
739  if (kx0 != nullptr) {
740  rIter = parMap.find(Kx0);
741  if (rIter != rIend) {
742  const map<parType, double>& pMap = rIter->second;
743  map<parType, double>::const_iterator pIter = pMap.begin();
744  map<parType, double>::const_iterator pIend = pMap.end();
745  while (pIter != pIend) {
746  const map<parType, double>::value_type& pEntry = *pIter++;
747  parType id = pEntry.first;
748  double pv = pEntry.second;
749  switch (id) {
750  case ptMin:
751  kx0->setPtMin(pv);
752  break;
753  case etaMax:
754  kx0->setEtaMax(pv);
755  break;
756  case massMin:
757  kx0->setMassMin(pv);
758  break;
759  case massMax:
760  kx0->setMassMax(pv);
761  break;
762  case probMin:
763  kx0->setProbMin(pv);
764  break;
765  case writeCandidate:
766  writeKx0 = (pv > 0);
767  break;
768  default:
769  break;
770  }
771  }
772  }
773  lKx0 = kx0->build();
774  if (allKx0)
775  sKx0.insert(lKx0.begin(), lKx0.end());
776  delete kx0;
777  }
778 
779  bool kx0Found = !lKx0.empty();
780 
781  // build and dump Bd -> JPsi Kx0
782 
783  if (recoBd && jPsiFound && kx0Found) {
784  BPHBdToJPsiKxBuilder* bd = new BPHBdToJPsiKxBuilder(es, lJPsi, lKx0);
785  rIter = parMap.find(Bd);
786  if (rIter != rIend) {
787  const map<parType, double>& pMap = rIter->second;
788  map<parType, double>::const_iterator pIter = pMap.begin();
789  map<parType, double>::const_iterator pIend = pMap.end();
790  while (pIter != pIend) {
791  const map<parType, double>::value_type& pEntry = *pIter++;
792  parType id = pEntry.first;
793  double pv = pEntry.second;
794  switch (id) {
795  case mPsiMin:
796  bd->setJPsiMassMin(pv);
797  break;
798  case mPsiMax:
799  bd->setJPsiMassMax(pv);
800  break;
801  case mKx0Min:
802  bd->setKxMassMin(pv);
803  break;
804  case mKx0Max:
805  bd->setKxMassMax(pv);
806  break;
807  case massMin:
808  bd->setMassMin(pv);
809  break;
810  case massMax:
811  bd->setMassMax(pv);
812  break;
813  case probMin:
814  bd->setProbMin(pv);
815  break;
816  case mFitMin:
817  bd->setMassFitMin(pv);
818  break;
819  case mFitMax:
820  bd->setMassFitMax(pv);
821  break;
822  case constrMJPsi:
823  bd->setConstr(pv > 0);
824  break;
825  case writeCandidate:
826  writeBd = (pv > 0);
827  break;
828  default:
829  break;
830  }
831  }
832  }
833 
834  lBd = bd->build();
835  delete bd;
836 
837  int iBd;
838  int nBd = lBd.size();
839  for (iBd = 0; iBd < nBd; ++iBd)
840  sKx0.insert(lBd[iBd]->getComp("Kx0"));
841  }
842  set<BPHRecoConstCandPtr>::const_iterator kx0_iter = sKx0.begin();
843  set<BPHRecoConstCandPtr>::const_iterator kx0_iend = sKx0.end();
844  lSd.reserve(sKx0.size());
845  while (kx0_iter != kx0_iend)
846  lSd.push_back(*kx0_iter++);
847 
848  // build and dump Phi
849 
850  vector<BPHPlusMinusConstCandPtr> lPhi;
851  BPHPhiToKKBuilder* phi = nullptr;
852  if (recoPkk && (jPsiFound || allPkk)) {
853  if (usePF)
854  phi = new BPHPhiToKKBuilder(
855  es, BPHRecoBuilder::createCollection(pfCands, "f"), BPHRecoBuilder::createCollection(pfCands, "f"));
856  else if (usePC)
857  phi = new BPHPhiToKKBuilder(
858  es, BPHRecoBuilder::createCollection(pcCands, "p"), BPHRecoBuilder::createCollection(pcCands, "p"));
859  else if (useGP)
860  phi = new BPHPhiToKKBuilder(
861  es, BPHRecoBuilder::createCollection(gpCands, "h"), BPHRecoBuilder::createCollection(gpCands, "h"));
862  }
863 
864  set<BPHRecoConstCandPtr> sPhi;
865 
866  if (phi != nullptr) {
867  rIter = parMap.find(Pkk);
868  if (rIter != rIend) {
869  const map<parType, double>& pMap = rIter->second;
870  map<parType, double>::const_iterator pIter = pMap.begin();
871  map<parType, double>::const_iterator pIend = pMap.end();
872  while (pIter != pIend) {
873  const map<parType, double>::value_type& pEntry = *pIter++;
874  parType id = pEntry.first;
875  double pv = pEntry.second;
876  switch (id) {
877  case ptMin:
878  phi->setPtMin(pv);
879  break;
880  case etaMax:
881  phi->setEtaMax(pv);
882  break;
883  case massMin:
884  phi->setMassMin(pv);
885  break;
886  case massMax:
887  phi->setMassMax(pv);
888  break;
889  case probMin:
890  phi->setProbMin(pv);
891  break;
892  case writeCandidate:
893  writePkk = (pv > 0);
894  break;
895  default:
896  break;
897  }
898  }
899  }
900  lPhi = phi->build();
901  if (allPkk)
902  sPhi.insert(lPhi.begin(), lPhi.end());
903  delete phi;
904  }
905 
906  bool phiFound = !lPhi.empty();
907 
908  // build and dump Bs
909 
910  if (recoBs && jPsiFound && phiFound) {
912  rIter = parMap.find(Bs);
913  if (rIter != rIend) {
914  const map<parType, double>& pMap = rIter->second;
915  map<parType, double>::const_iterator pIter = pMap.begin();
916  map<parType, double>::const_iterator pIend = pMap.end();
917  while (pIter != pIend) {
918  const map<parType, double>::value_type& pEntry = *pIter++;
919  parType id = pEntry.first;
920  double pv = pEntry.second;
921  switch (id) {
922  case mPsiMin:
923  bs->setJPsiMassMin(pv);
924  break;
925  case mPsiMax:
926  bs->setJPsiMassMax(pv);
927  break;
928  case mPhiMin:
929  bs->setPhiMassMin(pv);
930  break;
931  case mPhiMax:
932  bs->setPhiMassMax(pv);
933  break;
934  case massMin:
935  bs->setMassMin(pv);
936  break;
937  case massMax:
938  bs->setMassMax(pv);
939  break;
940  case probMin:
941  bs->setProbMin(pv);
942  break;
943  case mFitMin:
944  bs->setMassFitMin(pv);
945  break;
946  case mFitMax:
947  bs->setMassFitMax(pv);
948  break;
949  case constrMJPsi:
950  bs->setConstr(pv > 0);
951  break;
952  case writeCandidate:
953  writeBs = (pv > 0);
954  break;
955  default:
956  break;
957  }
958  }
959  }
960 
961  lBs = bs->build();
962  delete bs;
963 
964  int iBs;
965  int nBs = lBs.size();
966  for (iBs = 0; iBs < nBs; ++iBs)
967  sPhi.insert(lBs[iBs]->getComp("Phi"));
968  }
969  set<BPHRecoConstCandPtr>::const_iterator phi_iter = sPhi.begin();
970  set<BPHRecoConstCandPtr>::const_iterator phi_iend = sPhi.end();
971  lSs.reserve(sPhi.size());
972  while (phi_iter != phi_iend)
973  lSs.push_back(*phi_iter++);
974 
975  // build K0
976 
977  BPHK0sToPiPiBuilder* k0s = nullptr;
978  if (recoK0s && (jPsiFound || allK0s)) {
979  if (useK0)
980  k0s = new BPHK0sToPiPiBuilder(es, k0Cand.product(), "cfp");
981  else if (useKS)
982  k0s = new BPHK0sToPiPiBuilder(es, kSCand.product(), "cfp");
983  }
984  if (k0s != nullptr) {
985  rIter = parMap.find(K0s);
986  if (rIter != rIend) {
987  const map<parType, double>& pMap = rIter->second;
988  map<parType, double>::const_iterator pIter = pMap.begin();
989  map<parType, double>::const_iterator pIend = pMap.end();
990  while (pIter != pIend) {
991  const map<parType, double>::value_type& pEntry = *pIter++;
992  parType id = pEntry.first;
993  double pv = pEntry.second;
994  switch (id) {
995  case ptMin:
996  k0s->setPtMin(pv);
997  break;
998  case etaMax:
999  k0s->setEtaMax(pv);
1000  break;
1001  case massMin:
1002  k0s->setMassMin(pv);
1003  break;
1004  case massMax:
1005  k0s->setMassMax(pv);
1006  break;
1007  case probMin:
1008  k0s->setProbMin(pv);
1009  break;
1010  case writeCandidate:
1011  writeK0s = (pv > 0);
1012  break;
1013  default:
1014  break;
1015  }
1016  }
1017  }
1018  lK0 = k0s->build();
1019  delete k0s;
1020  }
1021 
1022  bool k0Found = !lK0.empty();
1023 
1024  // build Lambda0
1025 
1026  BPHLambda0ToPPiBuilder* l0s = nullptr;
1027  if (recoLambda0 && (jPsiFound || allLambda0)) {
1028  if (useL0)
1029  l0s = new BPHLambda0ToPPiBuilder(es, l0Cand.product(), "cfp");
1030  else if (useLS)
1031  l0s = new BPHLambda0ToPPiBuilder(es, lSCand.product(), "cfp");
1032  }
1033  if (l0s != nullptr) {
1034  rIter = parMap.find(Lambda0);
1035  if (rIter != rIend) {
1036  const map<parType, double>& pMap = rIter->second;
1037  map<parType, double>::const_iterator pIter = pMap.begin();
1038  map<parType, double>::const_iterator pIend = pMap.end();
1039  while (pIter != pIend) {
1040  const map<parType, double>::value_type& pEntry = *pIter++;
1041  parType id = pEntry.first;
1042  double pv = pEntry.second;
1043  switch (id) {
1044  case ptMin:
1045  l0s->setPtMin(pv);
1046  break;
1047  case etaMax:
1048  l0s->setEtaMax(pv);
1049  break;
1050  case massMin:
1051  l0s->setMassMin(pv);
1052  break;
1053  case massMax:
1054  l0s->setMassMax(pv);
1055  break;
1056  case probMin:
1057  l0s->setProbMin(pv);
1058  break;
1059  case writeCandidate:
1060  writeLambda0 = (pv > 0);
1061  break;
1062  default:
1063  break;
1064  }
1065  }
1066  }
1067  lL0 = l0s->build();
1068  delete l0s;
1069  }
1070 
1071  bool l0Found = !lL0.empty();
1072 
1073  // build and dump Bd -> JPsi K0s
1074 
1075  if (recoB0 && jPsiFound && k0Found) {
1077  rIter = parMap.find(B0);
1078  if (rIter != rIend) {
1079  const map<parType, double>& pMap = rIter->second;
1080  map<parType, double>::const_iterator pIter = pMap.begin();
1081  map<parType, double>::const_iterator pIend = pMap.end();
1082  while (pIter != pIend) {
1083  const map<parType, double>::value_type& pEntry = *pIter++;
1084  parType id = pEntry.first;
1085  double pv = pEntry.second;
1086  switch (id) {
1087  case mPsiMin:
1088  b0->setJPsiMassMin(pv);
1089  break;
1090  case mPsiMax:
1091  b0->setJPsiMassMax(pv);
1092  break;
1093  case mK0sMin:
1094  b0->setK0MassMin(pv);
1095  break;
1096  case mK0sMax:
1097  b0->setK0MassMax(pv);
1098  break;
1099  case massMin:
1100  b0->setMassMin(pv);
1101  break;
1102  case massMax:
1103  b0->setMassMax(pv);
1104  break;
1105  case probMin:
1106  b0->setProbMin(pv);
1107  break;
1108  case mFitMin:
1109  b0->setMassFitMin(pv);
1110  break;
1111  case mFitMax:
1112  b0->setMassFitMax(pv);
1113  break;
1114  case constrMJPsi:
1115  b0->setConstr(pv > 0);
1116  break;
1117  case writeCandidate:
1118  writeB0 = (pv > 0);
1119  break;
1120  default:
1121  break;
1122  }
1123  }
1124  }
1125 
1126  lB0 = b0->build();
1127  const map<const BPHRecoCandidate*, const BPHRecoCandidate*>& b0Map = b0->daughMap();
1128  daughMap.insert(b0Map.begin(), b0Map.end());
1129  delete b0;
1130  }
1131 
1132  // build and dump Lambdab -> JPsi Lambda0
1133 
1134  if (recoLambdab && jPsiFound && l0Found) {
1136  rIter = parMap.find(Lambdab);
1137  if (rIter != rIend) {
1138  const map<parType, double>& pMap = rIter->second;
1139  map<parType, double>::const_iterator pIter = pMap.begin();
1140  map<parType, double>::const_iterator pIend = pMap.end();
1141  while (pIter != pIend) {
1142  const map<parType, double>::value_type& pEntry = *pIter++;
1143  parType id = pEntry.first;
1144  double pv = pEntry.second;
1145  switch (id) {
1146  case mPsiMin:
1147  lb->setJPsiMassMin(pv);
1148  break;
1149  case mPsiMax:
1150  lb->setJPsiMassMax(pv);
1151  break;
1152  case mLambda0Min:
1153  lb->setLambda0MassMin(pv);
1154  break;
1155  case mLambda0Max:
1156  lb->setLambda0MassMax(pv);
1157  break;
1158  case massMin:
1159  lb->setMassMin(pv);
1160  break;
1161  case massMax:
1162  lb->setMassMax(pv);
1163  break;
1164  case probMin:
1165  lb->setProbMin(pv);
1166  break;
1167  case mFitMin:
1168  lb->setMassFitMin(pv);
1169  break;
1170  case mFitMax:
1171  lb->setMassFitMax(pv);
1172  break;
1173  case constrMJPsi:
1174  lb->setConstr(pv > 0);
1175  break;
1176  case writeCandidate:
1177  writeLambdab = (pv > 0);
1178  break;
1179  default:
1180  break;
1181  }
1182  }
1183  }
1184 
1185  lLb = lb->build();
1186  const map<const BPHRecoCandidate*, const BPHRecoCandidate*>& ldMap = lb->daughMap();
1187  daughMap.insert(ldMap.begin(), ldMap.end());
1188  delete lb;
1189  }
1190 
1191  // build and dump Bc
1192 
1193  BPHBcToJPsiPiBuilder* bc = nullptr;
1194  if (recoBc && jPsiFound) {
1195  if (usePF)
1196  bc = new BPHBcToJPsiPiBuilder(es, lJPsi, BPHRecoBuilder::createCollection(pfCands, "f"));
1197  else if (usePC)
1198  bc = new BPHBcToJPsiPiBuilder(es, lJPsi, BPHRecoBuilder::createCollection(pcCands, "p"));
1199  else if (useGP)
1200  bc = new BPHBcToJPsiPiBuilder(es, lJPsi, BPHRecoBuilder::createCollection(gpCands, "h"));
1201  }
1202 
1203  if (bc != nullptr) {
1204  rIter = parMap.find(Bc);
1205  if (rIter != rIend) {
1206  const map<parType, double>& pMap = rIter->second;
1207  map<parType, double>::const_iterator pIter = pMap.begin();
1208  map<parType, double>::const_iterator pIend = pMap.end();
1209  while (pIter != pIend) {
1210  const map<parType, double>::value_type& pEntry = *pIter++;
1211  parType id = pEntry.first;
1212  double pv = pEntry.second;
1213  switch (id) {
1214  case ptMin:
1215  bc->setPiPtMin(pv);
1216  break;
1217  case etaMax:
1218  bc->setPiEtaMax(pv);
1219  break;
1220  case mPsiMin:
1221  bc->setJPsiMassMin(pv);
1222  break;
1223  case mPsiMax:
1224  bc->setJPsiMassMax(pv);
1225  break;
1226  case massMin:
1227  bc->setMassMin(pv);
1228  break;
1229  case massMax:
1230  bc->setMassMax(pv);
1231  break;
1232  case probMin:
1233  bc->setProbMin(pv);
1234  break;
1235  case mFitMin:
1236  bc->setMassFitMin(pv);
1237  break;
1238  case mFitMax:
1239  bc->setMassFitMax(pv);
1240  break;
1241  case constrMJPsi:
1242  bc->setConstr(pv > 0);
1243  break;
1244  case writeCandidate:
1245  writeBc = (pv > 0);
1246  break;
1247  default:
1248  break;
1249  }
1250  }
1251  }
1252  lBc = bc->build();
1253  delete bc;
1254  }
1255 
1256  // build and dump Psi2S
1257 
1258  BPHPsi2SToJPsiPiPiBuilder* psi2S = nullptr;
1259  if (recoPsi2S && jPsiFound) {
1260  if (usePF)
1261  psi2S = new BPHPsi2SToJPsiPiPiBuilder(
1262  es, lJPsi, BPHRecoBuilder::createCollection(pfCands, "f"), BPHRecoBuilder::createCollection(pfCands, "f"));
1263  else if (usePC)
1264  psi2S = new BPHPsi2SToJPsiPiPiBuilder(
1265  es, lJPsi, BPHRecoBuilder::createCollection(pcCands, "p"), BPHRecoBuilder::createCollection(pcCands, "p"));
1266  else if (useGP)
1267  psi2S = new BPHPsi2SToJPsiPiPiBuilder(
1268  es, lJPsi, BPHRecoBuilder::createCollection(gpCands, "h"), BPHRecoBuilder::createCollection(gpCands, "h"));
1269  }
1270 
1271  if (psi2S != nullptr) {
1272  rIter = parMap.find(Psi2S);
1273  if (rIter != rIend) {
1274  const map<parType, double>& pMap = rIter->second;
1275  map<parType, double>::const_iterator pIter = pMap.begin();
1276  map<parType, double>::const_iterator pIend = pMap.end();
1277  while (pIter != pIend) {
1278  const map<parType, double>::value_type& pEntry = *pIter++;
1279  parType id = pEntry.first;
1280  double pv = pEntry.second;
1281  switch (id) {
1282  case ptMin:
1283  psi2S->setPiPtMin(pv);
1284  break;
1285  case etaMax:
1286  psi2S->setPiEtaMax(pv);
1287  break;
1288  case mPsiMin:
1289  psi2S->setJPsiMassMin(pv);
1290  break;
1291  case mPsiMax:
1292  psi2S->setJPsiMassMax(pv);
1293  break;
1294  case massMin:
1295  psi2S->setMassMin(pv);
1296  break;
1297  case massMax:
1298  psi2S->setMassMax(pv);
1299  break;
1300  case probMin:
1301  psi2S->setProbMin(pv);
1302  break;
1303  case mFitMin:
1304  psi2S->setMassFitMin(pv);
1305  break;
1306  case mFitMax:
1307  psi2S->setMassFitMax(pv);
1308  break;
1309  case constrMJPsi:
1310  psi2S->setConstr(pv > 0);
1311  break;
1312  case writeCandidate:
1313  writePsi2S = (pv > 0);
1314  break;
1315  default:
1316  break;
1317  }
1318  }
1319  }
1320  lPsi2S = psi2S->build();
1321  delete psi2S;
1322  }
1323 
1324  // build and dump X3872
1325 
1326  BPHX3872ToJPsiPiPiBuilder* x3872 = nullptr;
1327  if (recoX3872 && jPsiFound) {
1328  if (usePF)
1329  x3872 = new BPHX3872ToJPsiPiPiBuilder(
1330  es, lJPsi, BPHRecoBuilder::createCollection(pfCands, "f"), BPHRecoBuilder::createCollection(pfCands, "f"));
1331  else if (usePC)
1332  x3872 = new BPHX3872ToJPsiPiPiBuilder(
1333  es, lJPsi, BPHRecoBuilder::createCollection(pcCands, "p"), BPHRecoBuilder::createCollection(pcCands, "p"));
1334  else if (useGP)
1335  x3872 = new BPHX3872ToJPsiPiPiBuilder(
1336  es, lJPsi, BPHRecoBuilder::createCollection(gpCands, "h"), BPHRecoBuilder::createCollection(gpCands, "h"));
1337  }
1338 
1339  if (x3872 != nullptr) {
1340  rIter = parMap.find(X3872);
1341  if (rIter != rIend) {
1342  const map<parType, double>& pMap = rIter->second;
1343  map<parType, double>::const_iterator pIter = pMap.begin();
1344  map<parType, double>::const_iterator pIend = pMap.end();
1345  while (pIter != pIend) {
1346  const map<parType, double>::value_type& pEntry = *pIter++;
1347  parType id = pEntry.first;
1348  double pv = pEntry.second;
1349  switch (id) {
1350  case ptMin:
1351  x3872->setPiPtMin(pv);
1352  break;
1353  case etaMax:
1354  x3872->setPiEtaMax(pv);
1355  break;
1356  case mPsiMin:
1357  x3872->setJPsiMassMin(pv);
1358  break;
1359  case mPsiMax:
1360  x3872->setJPsiMassMax(pv);
1361  break;
1362  case massMin:
1363  x3872->setMassMin(pv);
1364  break;
1365  case massMax:
1366  x3872->setMassMax(pv);
1367  break;
1368  case probMin:
1369  x3872->setProbMin(pv);
1370  break;
1371  case mFitMin:
1372  x3872->setMassFitMin(pv);
1373  break;
1374  case mFitMax:
1375  x3872->setMassFitMax(pv);
1376  break;
1377  case constrMJPsi:
1378  x3872->setConstr(pv > 0);
1379  break;
1380  case writeCandidate:
1381  writeX3872 = (pv > 0);
1382  break;
1383  default:
1384  break;
1385  }
1386  }
1387  }
1388  lX3872 = x3872->build();
1389  delete x3872;
1390  }
1391 
1392  // merge Psi2S and X3872
1393  class ResTrkTrkCompare {
1394  public:
1395  bool operator()(const BPHRecoConstCandPtr& l, const BPHRecoConstCandPtr& r) const {
1396  vector<const reco::Track*> tl = l->tracks();
1397  vector<const reco::Track*> tr = r->tracks();
1398  if (tl.size() < tr.size())
1399  return true;
1400  sort(tl.begin(), tl.end());
1401  sort(tr.begin(), tr.end());
1402  int n = tr.size();
1403  int i;
1404  for (i = 0; i < n; ++i) {
1405  if (tl[i] < tr[i])
1406  return true;
1407  if (tl[i] > tr[i])
1408  return false;
1409  }
1410  return false;
1411  }
1412  } rttc;
1413  set<BPHRecoConstCandPtr, ResTrkTrkCompare> sjpPiPi(rttc);
1414  sjpPiPi.insert(lPsi2S.begin(), lPsi2S.end());
1415  sjpPiPi.insert(lX3872.begin(), lX3872.end());
1416  vector<BPHRecoConstCandPtr> ljpPiPi;
1417  ljpPiPi.insert(ljpPiPi.end(), sjpPiPi.begin(), sjpPiPi.end());
1418  bool jpPiPiFound = !ljpPiPi.empty();
1419 
1420  // build and dump Bp
1421 
1422  BPHBuToPsi2SKBuilder* bp = nullptr;
1423  if (recoBp && jpPiPiFound) {
1424  if (usePF)
1425  bp = new BPHBuToPsi2SKBuilder(es, ljpPiPi, BPHRecoBuilder::createCollection(pfCands, "f"));
1426  else if (usePC)
1427  bp = new BPHBuToPsi2SKBuilder(es, ljpPiPi, BPHRecoBuilder::createCollection(pcCands, "p"));
1428  else if (useGP)
1429  bp = new BPHBuToPsi2SKBuilder(es, ljpPiPi, BPHRecoBuilder::createCollection(gpCands, "h"));
1430  }
1431 
1432  if (bp != nullptr) {
1433  class BPHBuToPsi2SSelect : public BPHMassFitSelect {
1434  public:
1435  BPHBuToPsi2SSelect()
1436  : BPHMassFitSelect("Psi2S", BPHParticleMasses::psi2Mass, BPHParticleMasses::psi2MWidth, 5.0, 6.0) {}
1437  ~BPHBuToPsi2SSelect() override = default;
1438  bool accept(const BPHKinematicFit& cand) const override {
1439  const_cast<BPHRecoCandidate*>(cand.getComp("Psi2S").get())
1440  ->setIndependentFit("JPsi", true, BPHParticleMasses::jPsiMass, BPHParticleMasses::jPsiMWidth);
1442  }
1443  };
1444  bool mcJPsi = false;
1445  bool mcPsi2 = true;
1446  rIter = parMap.find(Bp);
1447  if (rIter != rIend) {
1448  const map<parType, double>& pMap = rIter->second;
1449  map<parType, double>::const_iterator pIter = pMap.begin();
1450  map<parType, double>::const_iterator pIend = pMap.end();
1451  while (pIter != pIend) {
1452  const map<parType, double>::value_type& pEntry = *pIter++;
1453  parType id = pEntry.first;
1454  double pv = pEntry.second;
1455  switch (id) {
1456  case ptMin:
1457  bp->setKPtMin(pv);
1458  break;
1459  case etaMax:
1460  bp->setKEtaMax(pv);
1461  break;
1462  case mPsiMin:
1463  bp->setPsi2SMassMin(pv);
1464  break;
1465  case mPsiMax:
1466  bp->setPsi2SMassMax(pv);
1467  break;
1468  case massMin:
1469  bp->setMassMin(pv);
1470  break;
1471  case massMax:
1472  bp->setMassMax(pv);
1473  break;
1474  case probMin:
1475  bp->setProbMin(pv);
1476  break;
1477  case mFitMin:
1478  bp->setMassFitMin(pv);
1479  break;
1480  case mFitMax:
1481  bp->setMassFitMax(pv);
1482  break;
1483  case constrMJPsi:
1484  mcJPsi = (pv > 0);
1485  break;
1486  case constrMPsi2:
1487  mcPsi2 = (pv > 0);
1488  break;
1489  case writeCandidate:
1490  writeBp = (pv > 0);
1491  break;
1492  default:
1493  break;
1494  }
1495  }
1496  }
1497  if (mcJPsi)
1498  bp->setMassFitSelect(mcPsi2 ? new BPHBuToPsi2SSelect
1499  : new BPHMassFitSelect("Psi2S/JPsi",
1502  bp->getMassFitMin(),
1503  bp->getMassFitMax()));
1504  else
1505  bp->setConstr(mcPsi2);
1506  lBp = bp->build();
1507  const map<const BPHRecoCandidate*, const BPHRecoCandidate*>& bpMap = bp->daughMap();
1508  daughMap.insert(bpMap.begin(), bpMap.end());
1509  delete bp;
1510  }
1511 
1512  return;
1513 }
BPHTokenWrapper< std::vector< reco::VertexCompositeCandidate > > l0CandsToken
virtual const reco::Candidate * originalReco(const reco::Candidate *daug) const
get the original particle from the clone
Analysis-level particle class.
static bool sameTrack(const reco::Candidate *lCand, const reco::Candidate *rCand, double minPDifference)
void setJPsiMassMax(double m)
BPHTokenWrapper< std::vector< reco::Vertex > > pVertexToken
std::vector< BPHPlusMinusConstCandPtr > lJPsi
void setJPsiMassMin(double m)
set cuts
BPHTokenWrapper< std::vector< reco::VertexCompositePtrCandidate > > kSCandsToken
std::vector< BPHPlusMinusConstCandPtr > lK0
edm::Ref< std::vector< reco::Vertex > > vertex_ref
void setMassMax(oniaType type, double m)
std::map< const BPHRecoCandidate *, compcc_ref > ccRefMap
BPHGenericPtr< const BPHRecoCandidate >::type BPHRecoConstCandPtr
static const double jPsiMWidth
std::vector< BPHPlusMinusConstCandPtr > lFull
const Point & position() const
position
Definition: Vertex.h:127
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
void setJPsiMassMin(double m)
T const * product() const
Definition: Handle.h:70
void setKPtMin(double pt)
set cuts
std::vector< BPHRecoConstCandPtr > lBs
std::vector< BPHRecoConstCandPtr > lBp
bool get(const edm::Event &ev, edm::Handle< Obj > &obj)
BPHTokenWrapper< std::vector< BPHTrackReference::candidate > > pcCandsToken
void setKEtaMax(double eta)
virtual std::vector< prod_ptr > build()
build candidates
static const double jPsiMass
void setPtMin(double pt)
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
BPHESTokenWrapper< MagneticField, IdealMagneticFieldRecord > magFieldToken
const pat::CompositeCandidate & composite() const override
get a composite by the simple sum of simple particles
void setLambda0MassMin(double m)
void setJPsiMassMax(double m)
std::map< const BPHRecoCandidate *, const BPHRecoCandidate * > daughMap
void setProbMin(oniaType type, double p)
bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb) override
std::map< std::string, parType > pMap
void setMassFitSelect(BPHMassFitSelect *mfs)
static BPHGenericCollection * createCollection(const edm::Handle< T > &collection, const std::string &list="cfhpmig")
std::vector< BPHRecoConstCandPtr > lSs
U second(std::pair< T, U > const &p)
double getConstrSigma(oniaType type) const
int TrackCharge
Definition: TrackCharge.h:4
virtual const reco::Candidate * getDaug(const std::string &name) const
BPHTokenWrapper< std::vector< reco::VertexCompositeCandidate > > k0CandsToken
std::vector< BPHPlusMinusConstCandPtr > getList(oniaType type, BPHRecoSelect *dSel=nullptr, BPHMomentumSelect *mSel=nullptr, BPHVertexSelect *vSel=nullptr, BPHFitSelect *kSel=nullptr)
BPHTokenWrapper< pat::MuonCollection > patMuonToken
void setEtaMax(double eta)
BPHTokenWrapper< std::vector< reco::PFCandidate > > pfCandsToken
def pv(vc)
Definition: MetAnalyzer.py:7
void setJPsiMassMin(double m)
set cuts
BPHTokenWrapper< std::vector< reco::VertexCompositePtrCandidate > > lSCandsToken
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
BPHTokenWrapper< std::vector< pat::CompositeCandidate > > ccCandsToken
std::pair< GlobalPoint, GlobalPoint > points() const override
std::vector< BPHRecoConstCandPtr > lSd
double getConstrMass(oniaType type) const
const T * userData(const std::string &key) const
Returns user-defined data. Returns NULL if the data is not present, or not of type T...
Definition: PATObject.h:322
virtual const reco::Vertex & vertex(VertexFitter< 5 > *fitter=nullptr, const reco::BeamSpot *bs=nullptr, const GlobalPoint *priorPos=nullptr, const GlobalError *priorError=nullptr) const
get reconstructed vertex
void setPsi2SMassMax(double m)
void setLambda0MassMax(double m)
BPHTokenWrapper< std::vector< pat::GenericParticle > > gpCandsToken
void setKEtaMax(double eta)
void setPiEtaMax(double eta)
std::map< const BPHRecoCandidate *, vertex_ref > pvRefMap
std::vector< BPHPlusMinusConstCandPtr > lL0
std::vector< BPHRecoConstCandPtr > lBu
std::vector< BPHRecoConstCandPtr > lLb
void setPiPtMin(double pt)
set cuts
void add(std::map< std::string, TH1 *> &h, TH1 *hist)
const edm::EventSetup * get() const
void setJPsiMassMin(double m)
void setConstr(oniaType type, double mass, double sigma)
void setEtaMax(oniaType type, double eta)
std::map< recoType, std::map< parType, double > > parMap
std::vector< BPHRecoConstCandPtr > lPsi2S
virtual const std::vector< const reco::Candidate * > & daughters() const
const std::map< const BPHRecoCandidate *, const BPHRecoCandidate * > & daughMap() const
get original daughters map
std::map< const BPHRecoCandidate *, const BPHRecoCandidate * > jPsiOMap
static constexpr float b0
std::vector< BPHRecoConstCandPtr > lB0
void setKPtMin(double pt)
set cuts
void setMassMin(oniaType type, double m)
std::vector< BPHRecoConstCandPtr > lBc
std::vector< BPHRecoConstCandPtr > lX3872
const std::map< const BPHRecoCandidate *, const BPHRecoCandidate * > & daughMap() const
get original daughters map
bool accept(const BPHKinematicFit &cand) const override
select particle
std::vector< BPHRecoConstCandPtr > lBd
Analysis-level muon class.
Definition: Muon.h:51
bool get(const edm::EventSetup &es, edm::ESHandle< Obj > &obj)
void setJPsiMassMax(double m)
void setPsi2SMassMin(double m)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
void setPtMin(oniaType type, double pt)
set cuts
void setJPsiMassMax(double m)
const_iterator begin() const
first daughter const_iterator
Definition: Candidate.h:143
void setPiPtMin(double pt)
set cuts

◆ fillDescriptions()

void BPHWriteSpecificDecay::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 244 of file BPHWriteSpecificDecay.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), and submitPVResolutionJobs::desc.

244  {
246  desc.add<string>("pVertexLabel", "");
247  desc.add<string>("patMuonLabel", "");
248  desc.add<string>("ccCandsLabel", "");
249  desc.add<string>("pfCandsLabel", "");
250  desc.add<string>("pcCandsLabel", "");
251  desc.add<string>("gpCandsLabel", "");
252  desc.add<string>("k0CandsLabel", "");
253  desc.add<string>("l0CandsLabel", "");
254  desc.add<string>("kSCandsLabel", "");
255  desc.add<string>("lSCandsLabel", "");
256  desc.add<string>("oniaName", "oniaCand");
257  desc.add<string>("sdName", "kx0Cand");
258  desc.add<string>("ssName", "phiCand");
259  desc.add<string>("buName", "buFitted");
260  desc.add<string>("bpName", "bpFitted");
261  desc.add<string>("bdName", "bdFitted");
262  desc.add<string>("bsName", "bsFitted");
263  desc.add<string>("k0Name", "k0Fitted");
264  desc.add<string>("l0Name", "l0Fitted");
265  desc.add<string>("b0Name", "b0Fitted");
266  desc.add<string>("lbName", "lbFitted");
267  desc.add<string>("bcName", "bcFitted");
268  desc.add<string>("psi2SName", "psi2SFitted");
269  desc.add<string>("x3872Name", "x3872Fitted");
270  desc.add<bool>("writeVertex", true);
271  desc.add<bool>("writeMomentum", true);
273  dpar.add<string>("name");
274  dpar.add<double>("ptMin", -2.0e35);
275  dpar.add<double>("etaMax", -2.0e35);
276  dpar.add<double>("mJPsiMin", -2.0e35);
277  dpar.add<double>("mJPsiMax", -2.0e35);
278  dpar.add<double>("mKx0Min", -2.0e35);
279  dpar.add<double>("mKx0Max", -2.0e35);
280  dpar.add<double>("mPhiMin", -2.0e35);
281  dpar.add<double>("mPhiMax", -2.0e35);
282  dpar.add<double>("mK0sMin", -2.0e35);
283  dpar.add<double>("mK0sMax", -2.0e35);
284  dpar.add<double>("mLambda0Min", -2.0e35);
285  dpar.add<double>("mLambda0Max", -2.0e35);
286  dpar.add<double>("massMin", -2.0e35);
287  dpar.add<double>("massMax", -2.0e35);
288  dpar.add<double>("probMin", -2.0e35);
289  dpar.add<double>("massFitMin", -2.0e35);
290  dpar.add<double>("massFitMax", -2.0e35);
291  dpar.add<double>("constrMass", -2.0e35);
292  dpar.add<double>("constrSigma", -2.0e35);
293  dpar.add<bool>("requireJPsi", true);
294  dpar.add<bool>("constrMJPsi", true);
295  dpar.add<bool>("constrMPsi2", true);
296  dpar.add<bool>("writeCandidate", true);
297  vector<edm::ParameterSet> rpar;
298  desc.addVPSet("recoSelect", dpar, rpar);
299  descriptions.add("bphWriteSpecificDecay", desc);
300  return;
301 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ produce()

void BPHWriteSpecificDecay::produce ( edm::Event ev,
const edm::EventSetup es 
)
override

Definition at line 303 of file BPHWriteSpecificDecay.cc.

References heavyFlavorDQMFirstStep_cff::b0Name, heavyFlavorDQMFirstStep_cff::bcName, heavyFlavorDQMFirstStep_cff::bdName, heavyFlavorDQMFirstStep_cff::bpName, heavyFlavorDQMFirstStep_cff::bsName, heavyFlavorDQMFirstStep_cff::buName, makeMEIFBenchmarkPlots::ev, ntuplemaker::fill, heavyFlavorDQMFirstStep_cff::k0Name, heavyFlavorDQMFirstStep_cff::l0Name, heavyFlavorDQMFirstStep_cff::lbName, heavyFlavorDQMFirstStep_cff::oniaName, heavyFlavorDQMFirstStep_cff::psi2SName, heavyFlavorDQMFirstStep_cff::sdName, heavyFlavorDQMFirstStep_cff::ssName, BPHRecoCandidate::transientTrackBuilder, and writeEcalDQMStatus::write.

303  {
305  fill(ev, ew);
306  if (writeOnia)
307  write(ev, lFull, oniaName);
308  if (writeKx0)
309  write(ev, lSd, sdName);
310  if (writePkk)
311  write(ev, lSs, ssName);
312  if (writeBu)
313  write(ev, lBu, buName);
314  if (writeBp)
315  write(ev, lBp, bpName);
316  if (writeBd)
317  write(ev, lBd, bdName);
318  if (writeBs)
319  write(ev, lBs, bsName);
320  if (writeK0s)
321  write(ev, lK0, k0Name);
322  if (writeLambda0)
323  write(ev, lL0, l0Name);
324  if (writeB0)
325  write(ev, lB0, b0Name);
326  if (writeLambdab)
327  write(ev, lLb, lbName);
328  if (writeBc)
329  write(ev, lBc, bcName);
330  if (writePsi2S)
332  if (writeX3872)
334  return;
335 }
std::vector< BPHPlusMinusConstCandPtr > lK0
std::vector< BPHPlusMinusConstCandPtr > lFull
std::vector< BPHRecoConstCandPtr > lBs
std::vector< BPHRecoConstCandPtr > lBp
BPHESTokenWrapper< TransientTrackBuilder, TransientTrackRecord > ttBToken
std::vector< BPHRecoConstCandPtr > lSs
std::vector< BPHRecoConstCandPtr > lSd
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::vector< BPHPlusMinusConstCandPtr > lL0
std::vector< BPHRecoConstCandPtr > lBu
std::vector< BPHRecoConstCandPtr > lLb
std::vector< BPHRecoConstCandPtr > lPsi2S
std::vector< BPHRecoConstCandPtr > lB0
std::vector< BPHRecoConstCandPtr > lBc
std::vector< BPHRecoConstCandPtr > lX3872
std::vector< BPHRecoConstCandPtr > lBd

◆ setRecoParameters()

void BPHWriteSpecificDecay::setRecoParameters ( const edm::ParameterSet ps)
private

Definition at line 1515 of file BPHWriteSpecificDecay.cc.

References mps_splice::entry, personalPlayback::fn, edm::ParameterSet::getParameter(), l1ctLayer2EG_cff::id, Skims_PA_cff::name, and MetAnalyzer::pv().

1515  {
1516  const string& name = ps.getParameter<string>("name");
1517  bool writeCandidate = ps.getParameter<bool>("writeCandidate");
1518  switch (rMap[name]) {
1519  case Onia:
1520  recoOnia = true;
1522  break;
1523  case Pmm:
1524  case Psi1:
1525  case Psi2:
1526  case Ups:
1527  case Ups1:
1528  case Ups2:
1529  case Ups3:
1530  recoOnia = true;
1531  break;
1532  case Kx0:
1533  recoKx0 = true;
1534  allKx0 = false;
1536  break;
1537  case Pkk:
1538  recoPkk = true;
1539  allPkk = false;
1541  break;
1542  case Bu:
1543  recoBu = true;
1545  break;
1546  case Bp:
1547  recoBp = true;
1549  break;
1550  case Bd:
1551  recoBd = true;
1553  break;
1554  case Bs:
1555  recoBs = true;
1557  break;
1558  case K0s:
1559  recoK0s = true;
1560  allK0s = false;
1562  break;
1563  case Lambda0:
1564  recoLambda0 = true;
1565  allLambda0 = false;
1567  break;
1568  case B0:
1569  recoB0 = true;
1571  break;
1572  case Lambdab:
1573  recoLambdab = true;
1575  break;
1576  case Bc:
1577  recoBc = true;
1579  break;
1580  case Psi2S:
1581  recoPsi2S = true;
1583  break;
1584  case X3872:
1585  recoX3872 = true;
1587  break;
1588  }
1589 
1590  map<string, parType>::const_iterator pIter = pMap.begin();
1591  map<string, parType>::const_iterator pIend = pMap.end();
1592  while (pIter != pIend) {
1593  const map<string, parType>::value_type& entry = *pIter++;
1594  const string& pn = entry.first;
1595  parType id = entry.second;
1596  double pv = ps.getParameter<double>(pn);
1597  if (pv > -1.0e35)
1598  edm::LogVerbatim("Configuration") << "BPHWriteSpecificDecay::setRecoParameters: set " << pn << " for " << name
1599  << " : " << (parMap[rMap[name]][id] = pv);
1600  }
1601 
1602  map<string, parType>::const_iterator fIter = fMap.begin();
1603  map<string, parType>::const_iterator fIend = fMap.end();
1604  while (fIter != fIend) {
1605  const map<string, parType>::value_type& entry = *fIter++;
1606  const string& fn = entry.first;
1607  parType id = entry.second;
1608  double pv = (ps.getParameter<bool>(fn) ? 1 : -1);
1609  if (pv > -1.0e35)
1610  edm::LogVerbatim("Configuration") << "BPHWriteSpecificDecay::setRecoParameters: set " << fn << " for " << name
1611  << " : " << (parMap[rMap[name]][id] = pv);
1612  }
1613 
1614  return;
1615 }
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::map< std::string, recoType > rMap
std::map< std::string, parType > pMap
def pv(vc)
Definition: MetAnalyzer.py:7
std::map< recoType, std::map< parType, double > > parMap
std::map< std::string, parType > fMap

◆ write()

template<class T >
edm::OrphanHandle<pat::CompositeCandidateCollection> BPHWriteSpecificDecay::write ( edm::Event ev,
const std::vector< T > &  list,
const std::string &  name 
)
inlineprivate

Definition at line 223 of file BPHWriteSpecificDecay.h.

References addTrackModes(), gpuPixelDoublets::cc, ccRefMap, submitPVResolutionJobs::count, daughMap, makeMEIFBenchmarkPlots::ev, BPHDecayMomentum::getComp(), mps_fire::i, BPHPlusMinusCandidate::isCowboy(), edm::Ref< C, T, F >::isNonnull(), dqmiolumiharvest::j, jPsiOMap, KinematicState::kinematicParameters(), visualization-live-secondInstance_cfg::m, KinematicState::mass(), KinematicParameters::momentum(), dqmiodumpmetadata::n, Skims_PA_cff::name, pvRefMap, AlCaHLTBitMon_QueryRunRegistry::string, writeMomentum, and writeVertex.

225  {
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  }
Analysis-level particle class.
std::map< const BPHRecoCandidate *, compcc_ref > ccRefMap
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
edm::Ref< pat::CompositeCandidateCollection > compcc_ref
bool isCowboy() const
get cowboy/sailor classification
std::map< const BPHRecoCandidate *, const BPHRecoCandidate * > daughMap
GlobalVector momentum() const
KinematicParameters const & kinematicParameters() const
virtual BPHRecoConstCandPtr getComp(const std::string &name) const
std::vector< CompositeCandidate > CompositeCandidateCollection
std::map< const BPHRecoCandidate *, vertex_ref > pvRefMap
ParticleMass mass() const
std::map< const BPHRecoCandidate *, const BPHRecoCandidate * > jPsiOMap
long double T
static void addTrackModes(const std::string &name, const BPHRecoCandidate &cand, std::string &modes, bool &count)

Member Data Documentation

◆ allK0s

bool BPHWriteSpecificDecay::allK0s
private

Definition at line 172 of file BPHWriteSpecificDecay.h.

◆ allKx0

bool BPHWriteSpecificDecay::allKx0
private

Definition at line 170 of file BPHWriteSpecificDecay.h.

◆ allLambda0

bool BPHWriteSpecificDecay::allLambda0
private

Definition at line 173 of file BPHWriteSpecificDecay.h.

◆ allPkk

bool BPHWriteSpecificDecay::allPkk
private

Definition at line 171 of file BPHWriteSpecificDecay.h.

◆ b0Name

std::string BPHWriteSpecificDecay::b0Name
private

Definition at line 96 of file BPHWriteSpecificDecay.h.

◆ bcName

std::string BPHWriteSpecificDecay::bcName
private

Definition at line 98 of file BPHWriteSpecificDecay.h.

◆ bdName

std::string BPHWriteSpecificDecay::bdName
private

Definition at line 92 of file BPHWriteSpecificDecay.h.

◆ bpName

std::string BPHWriteSpecificDecay::bpName
private

Definition at line 91 of file BPHWriteSpecificDecay.h.

◆ bsName

std::string BPHWriteSpecificDecay::bsName
private

Definition at line 93 of file BPHWriteSpecificDecay.h.

◆ buName

std::string BPHWriteSpecificDecay::buName
private

Definition at line 90 of file BPHWriteSpecificDecay.h.

◆ ccCandsLabel

std::string BPHWriteSpecificDecay::ccCandsLabel
private

Definition at line 53 of file BPHWriteSpecificDecay.h.

◆ ccCandsToken

BPHTokenWrapper<std::vector<pat::CompositeCandidate> > BPHWriteSpecificDecay::ccCandsToken
private

Definition at line 67 of file BPHWriteSpecificDecay.h.

◆ ccRefMap

std::map<const BPHRecoCandidate*, compcc_ref> BPHWriteSpecificDecay::ccRefMap
private

Definition at line 214 of file BPHWriteSpecificDecay.h.

Referenced by write().

◆ daughMap

std::map<const BPHRecoCandidate*, const BPHRecoCandidate*> BPHWriteSpecificDecay::daughMap
private

Definition at line 210 of file BPHWriteSpecificDecay.h.

Referenced by write().

◆ fMap

std::map<std::string, parType> BPHWriteSpecificDecay::fMap
private

Definition at line 152 of file BPHWriteSpecificDecay.h.

◆ gpCandsLabel

std::string BPHWriteSpecificDecay::gpCandsLabel
private

Definition at line 56 of file BPHWriteSpecificDecay.h.

◆ gpCandsToken

BPHTokenWrapper<std::vector<pat::GenericParticle> > BPHWriteSpecificDecay::gpCandsToken
private

Definition at line 70 of file BPHWriteSpecificDecay.h.

◆ jPsiOMap

std::map<const BPHRecoCandidate*, const BPHRecoCandidate*> BPHWriteSpecificDecay::jPsiOMap
private

Definition at line 209 of file BPHWriteSpecificDecay.h.

Referenced by write().

◆ k0CandsLabel

std::string BPHWriteSpecificDecay::k0CandsLabel
private

Definition at line 57 of file BPHWriteSpecificDecay.h.

◆ k0CandsToken

BPHTokenWrapper<std::vector<reco::VertexCompositeCandidate> > BPHWriteSpecificDecay::k0CandsToken
private

Definition at line 71 of file BPHWriteSpecificDecay.h.

◆ k0Name

std::string BPHWriteSpecificDecay::k0Name
private

Definition at line 94 of file BPHWriteSpecificDecay.h.

◆ kSCandsLabel

std::string BPHWriteSpecificDecay::kSCandsLabel
private

Definition at line 59 of file BPHWriteSpecificDecay.h.

◆ kSCandsToken

BPHTokenWrapper<std::vector<reco::VertexCompositePtrCandidate> > BPHWriteSpecificDecay::kSCandsToken
private

Definition at line 73 of file BPHWriteSpecificDecay.h.

◆ l0CandsLabel

std::string BPHWriteSpecificDecay::l0CandsLabel
private

Definition at line 58 of file BPHWriteSpecificDecay.h.

◆ l0CandsToken

BPHTokenWrapper<std::vector<reco::VertexCompositeCandidate> > BPHWriteSpecificDecay::l0CandsToken
private

Definition at line 72 of file BPHWriteSpecificDecay.h.

◆ l0Name

std::string BPHWriteSpecificDecay::l0Name
private

Definition at line 95 of file BPHWriteSpecificDecay.h.

◆ lB0

std::vector<BPHRecoConstCandPtr> BPHWriteSpecificDecay::lB0
private

Definition at line 203 of file BPHWriteSpecificDecay.h.

◆ lBc

std::vector<BPHRecoConstCandPtr> BPHWriteSpecificDecay::lBc
private

Definition at line 205 of file BPHWriteSpecificDecay.h.

◆ lBd

std::vector<BPHRecoConstCandPtr> BPHWriteSpecificDecay::lBd
private

Definition at line 199 of file BPHWriteSpecificDecay.h.

◆ lbName

std::string BPHWriteSpecificDecay::lbName
private

Definition at line 97 of file BPHWriteSpecificDecay.h.

◆ lBp

std::vector<BPHRecoConstCandPtr> BPHWriteSpecificDecay::lBp
private

Definition at line 198 of file BPHWriteSpecificDecay.h.

◆ lBs

std::vector<BPHRecoConstCandPtr> BPHWriteSpecificDecay::lBs
private

Definition at line 200 of file BPHWriteSpecificDecay.h.

◆ lBu

std::vector<BPHRecoConstCandPtr> BPHWriteSpecificDecay::lBu
private

Definition at line 197 of file BPHWriteSpecificDecay.h.

◆ lFull

std::vector<BPHPlusMinusConstCandPtr> BPHWriteSpecificDecay::lFull
private

Definition at line 193 of file BPHWriteSpecificDecay.h.

◆ lJPsi

std::vector<BPHPlusMinusConstCandPtr> BPHWriteSpecificDecay::lJPsi
private

Definition at line 194 of file BPHWriteSpecificDecay.h.

◆ lK0

std::vector<BPHPlusMinusConstCandPtr> BPHWriteSpecificDecay::lK0
private

Definition at line 201 of file BPHWriteSpecificDecay.h.

◆ lL0

std::vector<BPHPlusMinusConstCandPtr> BPHWriteSpecificDecay::lL0
private

Definition at line 202 of file BPHWriteSpecificDecay.h.

◆ lLb

std::vector<BPHRecoConstCandPtr> BPHWriteSpecificDecay::lLb
private

Definition at line 204 of file BPHWriteSpecificDecay.h.

◆ lPsi2S

std::vector<BPHRecoConstCandPtr> BPHWriteSpecificDecay::lPsi2S
private

Definition at line 206 of file BPHWriteSpecificDecay.h.

◆ lSCandsLabel

std::string BPHWriteSpecificDecay::lSCandsLabel
private

Definition at line 60 of file BPHWriteSpecificDecay.h.

◆ lSCandsToken

BPHTokenWrapper<std::vector<reco::VertexCompositePtrCandidate> > BPHWriteSpecificDecay::lSCandsToken
private

Definition at line 74 of file BPHWriteSpecificDecay.h.

◆ lSd

std::vector<BPHRecoConstCandPtr> BPHWriteSpecificDecay::lSd
private

Definition at line 195 of file BPHWriteSpecificDecay.h.

◆ lSs

std::vector<BPHRecoConstCandPtr> BPHWriteSpecificDecay::lSs
private

Definition at line 196 of file BPHWriteSpecificDecay.h.

◆ lX3872

std::vector<BPHRecoConstCandPtr> BPHWriteSpecificDecay::lX3872
private

Definition at line 207 of file BPHWriteSpecificDecay.h.

◆ magFieldToken

BPHESTokenWrapper<MagneticField, IdealMagneticFieldRecord> BPHWriteSpecificDecay::magFieldToken
private

Definition at line 63 of file BPHWriteSpecificDecay.h.

◆ oniaName

std::string BPHWriteSpecificDecay::oniaName
private

Definition at line 87 of file BPHWriteSpecificDecay.h.

◆ parMap

std::map<recoType, std::map<parType, double> > BPHWriteSpecificDecay::parMap
private

Definition at line 153 of file BPHWriteSpecificDecay.h.

◆ patMuonLabel

std::string BPHWriteSpecificDecay::patMuonLabel
private

Definition at line 52 of file BPHWriteSpecificDecay.h.

◆ patMuonToken

BPHTokenWrapper<pat::MuonCollection> BPHWriteSpecificDecay::patMuonToken
private

Definition at line 66 of file BPHWriteSpecificDecay.h.

◆ pcCandsLabel

std::string BPHWriteSpecificDecay::pcCandsLabel
private

Definition at line 55 of file BPHWriteSpecificDecay.h.

◆ pcCandsToken

BPHTokenWrapper<std::vector<BPHTrackReference::candidate> > BPHWriteSpecificDecay::pcCandsToken
private

Definition at line 69 of file BPHWriteSpecificDecay.h.

◆ pfCandsLabel

std::string BPHWriteSpecificDecay::pfCandsLabel
private

Definition at line 54 of file BPHWriteSpecificDecay.h.

◆ pfCandsToken

BPHTokenWrapper<std::vector<reco::PFCandidate> > BPHWriteSpecificDecay::pfCandsToken
private

Definition at line 68 of file BPHWriteSpecificDecay.h.

◆ pMap

std::map<std::string, parType> BPHWriteSpecificDecay::pMap
private

Definition at line 151 of file BPHWriteSpecificDecay.h.

◆ psi2SName

std::string BPHWriteSpecificDecay::psi2SName
private

Definition at line 99 of file BPHWriteSpecificDecay.h.

◆ pVertexLabel

std::string BPHWriteSpecificDecay::pVertexLabel
private

Definition at line 51 of file BPHWriteSpecificDecay.h.

◆ pVertexToken

BPHTokenWrapper<std::vector<reco::Vertex> > BPHWriteSpecificDecay::pVertexToken
private

Definition at line 65 of file BPHWriteSpecificDecay.h.

◆ pvRefMap

std::map<const BPHRecoCandidate*, vertex_ref> BPHWriteSpecificDecay::pvRefMap
private

Definition at line 212 of file BPHWriteSpecificDecay.h.

Referenced by write().

◆ recoB0

bool BPHWriteSpecificDecay::recoB0
private

Definition at line 164 of file BPHWriteSpecificDecay.h.

◆ recoBc

bool BPHWriteSpecificDecay::recoBc
private

Definition at line 166 of file BPHWriteSpecificDecay.h.

◆ recoBd

bool BPHWriteSpecificDecay::recoBd
private

Definition at line 160 of file BPHWriteSpecificDecay.h.

◆ recoBp

bool BPHWriteSpecificDecay::recoBp
private

Definition at line 159 of file BPHWriteSpecificDecay.h.

◆ recoBs

bool BPHWriteSpecificDecay::recoBs
private

Definition at line 161 of file BPHWriteSpecificDecay.h.

◆ recoBu

bool BPHWriteSpecificDecay::recoBu
private

Definition at line 158 of file BPHWriteSpecificDecay.h.

◆ recoK0s

bool BPHWriteSpecificDecay::recoK0s
private

Definition at line 162 of file BPHWriteSpecificDecay.h.

◆ recoKx0

bool BPHWriteSpecificDecay::recoKx0
private

Definition at line 156 of file BPHWriteSpecificDecay.h.

◆ recoLambda0

bool BPHWriteSpecificDecay::recoLambda0
private

Definition at line 163 of file BPHWriteSpecificDecay.h.

◆ recoLambdab

bool BPHWriteSpecificDecay::recoLambdab
private

Definition at line 165 of file BPHWriteSpecificDecay.h.

◆ recoOnia

bool BPHWriteSpecificDecay::recoOnia
private

Definition at line 155 of file BPHWriteSpecificDecay.h.

◆ recoPkk

bool BPHWriteSpecificDecay::recoPkk
private

Definition at line 157 of file BPHWriteSpecificDecay.h.

◆ recoPsi2S

bool BPHWriteSpecificDecay::recoPsi2S
private

Definition at line 167 of file BPHWriteSpecificDecay.h.

◆ recoX3872

bool BPHWriteSpecificDecay::recoX3872
private

Definition at line 168 of file BPHWriteSpecificDecay.h.

◆ rMap

std::map<std::string, recoType> BPHWriteSpecificDecay::rMap
private

Definition at line 150 of file BPHWriteSpecificDecay.h.

◆ sdName

std::string BPHWriteSpecificDecay::sdName
private

Definition at line 88 of file BPHWriteSpecificDecay.h.

◆ ssName

std::string BPHWriteSpecificDecay::ssName
private

Definition at line 89 of file BPHWriteSpecificDecay.h.

◆ ttBToken

BPHESTokenWrapper<TransientTrackBuilder, TransientTrackRecord> BPHWriteSpecificDecay::ttBToken
private

Definition at line 64 of file BPHWriteSpecificDecay.h.

◆ useCC

bool BPHWriteSpecificDecay::useCC
private

Definition at line 78 of file BPHWriteSpecificDecay.h.

◆ useGP

bool BPHWriteSpecificDecay::useGP
private

Definition at line 81 of file BPHWriteSpecificDecay.h.

◆ useK0

bool BPHWriteSpecificDecay::useK0
private

Definition at line 82 of file BPHWriteSpecificDecay.h.

◆ useKS

bool BPHWriteSpecificDecay::useKS
private

Definition at line 84 of file BPHWriteSpecificDecay.h.

◆ useL0

bool BPHWriteSpecificDecay::useL0
private

Definition at line 83 of file BPHWriteSpecificDecay.h.

◆ useLS

bool BPHWriteSpecificDecay::useLS
private

Definition at line 85 of file BPHWriteSpecificDecay.h.

◆ usePC

bool BPHWriteSpecificDecay::usePC
private

Definition at line 80 of file BPHWriteSpecificDecay.h.

◆ usePF

bool BPHWriteSpecificDecay::usePF
private

Definition at line 79 of file BPHWriteSpecificDecay.h.

◆ usePM

bool BPHWriteSpecificDecay::usePM
private

Definition at line 77 of file BPHWriteSpecificDecay.h.

◆ usePV

bool BPHWriteSpecificDecay::usePV
private

Definition at line 76 of file BPHWriteSpecificDecay.h.

◆ writeB0

bool BPHWriteSpecificDecay::writeB0
private

Definition at line 184 of file BPHWriteSpecificDecay.h.

◆ writeBc

bool BPHWriteSpecificDecay::writeBc
private

Definition at line 186 of file BPHWriteSpecificDecay.h.

◆ writeBd

bool BPHWriteSpecificDecay::writeBd
private

Definition at line 180 of file BPHWriteSpecificDecay.h.

◆ writeBp

bool BPHWriteSpecificDecay::writeBp
private

Definition at line 179 of file BPHWriteSpecificDecay.h.

◆ writeBs

bool BPHWriteSpecificDecay::writeBs
private

Definition at line 181 of file BPHWriteSpecificDecay.h.

◆ writeBu

bool BPHWriteSpecificDecay::writeBu
private

Definition at line 178 of file BPHWriteSpecificDecay.h.

◆ writeK0s

bool BPHWriteSpecificDecay::writeK0s
private

Definition at line 182 of file BPHWriteSpecificDecay.h.

◆ writeKx0

bool BPHWriteSpecificDecay::writeKx0
private

Definition at line 176 of file BPHWriteSpecificDecay.h.

◆ writeLambda0

bool BPHWriteSpecificDecay::writeLambda0
private

Definition at line 183 of file BPHWriteSpecificDecay.h.

◆ writeLambdab

bool BPHWriteSpecificDecay::writeLambdab
private

Definition at line 185 of file BPHWriteSpecificDecay.h.

◆ writeMomentum

bool BPHWriteSpecificDecay::writeMomentum
private

Definition at line 191 of file BPHWriteSpecificDecay.h.

Referenced by write().

◆ writeOnia

bool BPHWriteSpecificDecay::writeOnia
private

Definition at line 175 of file BPHWriteSpecificDecay.h.

◆ writePkk

bool BPHWriteSpecificDecay::writePkk
private

Definition at line 177 of file BPHWriteSpecificDecay.h.

◆ writePsi2S

bool BPHWriteSpecificDecay::writePsi2S
private

Definition at line 187 of file BPHWriteSpecificDecay.h.

◆ writeVertex

bool BPHWriteSpecificDecay::writeVertex
private

Definition at line 190 of file BPHWriteSpecificDecay.h.

Referenced by write().

◆ writeX3872

bool BPHWriteSpecificDecay::writeX3872
private

Definition at line 188 of file BPHWriteSpecificDecay.h.

◆ x3872Name

std::string BPHWriteSpecificDecay::x3872Name
private

Definition at line 100 of file BPHWriteSpecificDecay.h.