CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Private Attributes
BPHWriteSpecificDecay Class Reference

#include <BPHWriteSpecificDecay.h>

Inheritance diagram for BPHWriteSpecificDecay:
BPHAnalyzerWrapper< BPHModuleWrapper::one_producer > edm::one::EDProducer< T > edm::one::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

void beginJob () override
 
 BPHWriteSpecificDecay (const edm::ParameterSet &ps)
 
void endJob () override
 
virtual void fill (edm::Event &ev, const edm::EventSetup &es)
 
void produce (edm::Event &ev, const edm::EventSetup &es) override
 
 ~BPHWriteSpecificDecay () override
 
- Public Member Functions inherited from edm::one::EDProducer< T >
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
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
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsInputProcessBlocks () const final
 
bool wantsProcessBlocks () const final
 
- Public Member Functions inherited from edm::one::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDProducerBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector
< edm::ProductResolverIndex >
const & 
indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector
< edm::ProductResolverIndex >
const & 
putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
 ~ProducerBase () noexcept(false) override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESProxyIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex >
const & 
esGetTokenRecordIndicesVector (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector
< ProductResolverIndexAndSkipBit >
const & 
itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::array< std::vector< ModuleDescription const * > *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from edm::one::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 

Private Types

typedef edm::Ref
< pat::CompositeCandidateCollection
compcc_ref
 
enum  parType {
  ptMin, etaMax, mPsiMin, mPsiMax,
  mKx0Min, mKx0Max, mPhiMin, mPhiMax,
  mK0sMin, mK0sMax, mLambda0Min, mLambda0Max,
  massMin, massMax, probMin, mFitMin,
  mFitMax, constrMass, constrSigma, constrMJPsi,
  writeCandidate
}
 
enum  recoType {
  Onia, Pmm, Psi1, Psi2,
  Ups, Ups1, Ups2, Ups3,
  Kx0, Pkk, Bu, Bd,
  Bs, K0s, Lambda0, B0,
  Lambdab, Bc, 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::CompositeCandidateCollection
write (edm::Event &ev, const std::vector< T > &list, const std::string &name)
 

Private Attributes

std::string b0Name
 
std::string bcName
 
std::string bdName
 
std::string bsName
 
std::string buName
 
std::string ccCandsLabel
 
BPHTokenWrapper< std::vector
< pat::CompositeCandidate > > 
ccCandsToken
 
std::map< const
BPHRecoCandidate *, compcc_ref
ccRefMap
 
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< BPHRecoConstCandPtrlBs
 
std::vector< BPHRecoConstCandPtrlBu
 
std::vector
< BPHPlusMinusConstCandPtr
lFull
 
std::vector
< BPHPlusMinusConstCandPtr
lJPsi
 
std::vector
< BPHPlusMinusConstCandPtr
lK0
 
std::vector
< BPHPlusMinusConstCandPtr
lL0
 
std::vector< BPHRecoConstCandPtrlLb
 
std::string lSCandsLabel
 
BPHTokenWrapper< std::vector
< reco::VertexCompositePtrCandidate > > 
lSCandsToken
 
std::vector< BPHRecoConstCandPtrlSd
 
std::vector< BPHRecoConstCandPtrlSs
 
std::vector< BPHRecoConstCandPtrlX3872
 
std::string oniaName
 
std::map< recoType, std::map
< parType, double > > 
parMap
 
std::string patMuonLabel
 
BPHTokenWrapper
< pat::MuonCollection
patMuonToken
 
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 pVertexLabel
 
BPHTokenWrapper< std::vector
< reco::Vertex > > 
pVertexToken
 
std::map< const
BPHRecoCandidate *, vertex_ref
pvRefMap
 
bool recoB0
 
bool recoBc
 
bool recoBd
 
bool recoBs
 
bool recoBu
 
bool recoK0s
 
bool recoKx0
 
bool recoLambda0
 
bool recoLambdab
 
bool recoOnia
 
bool recoPkk
 
bool recoX3872
 
std::map< std::string, recoTyperMap
 
std::string sdName
 
std::string ssName
 
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 writeBs
 
bool writeBu
 
bool writeK0s
 
bool writeKx0
 
bool writeLambda0
 
bool writeLambdab
 
bool writeMomentum
 
bool writeOnia
 
bool writePkk
 
bool writeVertex
 
bool writeX3872
 
std::string x3872Name
 

Additional Inherited Members

- Public Types inherited from edm::one::EDProducerBase
typedef EDProducerBase ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Protected Member Functions inherited from BPHAnalyzerWrapper< BPHModuleWrapper::one_producer >
void consume (BPHTokenWrapper< Obj > &tw, const std::string &label)
 
void consume (BPHTokenWrapper< Obj > &tw, const edm::InputTag &tag)
 
- Protected Member Functions inherited from edm::ProducerBase
ProducesCollector producesCollector ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
void resetItemsToGetFrom (BranchType iType)
 

Detailed Description

Definition at line 34 of file BPHWriteSpecificDecay.h.

Member Typedef Documentation

Definition at line 190 of file BPHWriteSpecificDecay.h.

Definition at line 188 of file BPHWriteSpecificDecay.h.

Member Enumeration Documentation

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

Definition at line 115 of file BPHWriteSpecificDecay.h.

115  {
116  ptMin,
117  etaMax,
118  mPsiMin,
119  mPsiMax,
120  mKx0Min,
121  mKx0Max,
122  mPhiMin,
123  mPhiMax,
124  mK0sMin,
125  mK0sMax,
126  mLambda0Min,
127  mLambda0Max,
128  massMin,
129  massMax,
130  probMin,
131  mFitMin,
132  mFitMax,
133  constrMass,
134  constrSigma,
135  constrMJPsi,
137  };
Enumerator
Onia 
Pmm 
Psi1 
Psi2 
Ups 
Ups1 
Ups2 
Ups3 
Kx0 
Pkk 
Bu 
Bd 
Bs 
K0s 
Lambda0 
B0 
Lambdab 
Bc 
X3872 

Definition at line 94 of file BPHWriteSpecificDecay.h.

Constructor & Destructor Documentation

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

Definition at line 56 of file BPHWriteSpecificDecay.cc.

References Puppi_cff::etaMax, edm::ParameterSet::getParameter(), ptMin, recoSelectForWrite_cfi::recoSelect, SET_PAR, and HLT_FULL_cff::usePV.

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

Definition at line 215 of file BPHWriteSpecificDecay.cc.

215 {}

Member Function Documentation

void BPHWriteSpecificDecay::beginJob ( void  )
overridevirtual

Reimplemented from edm::one::EDProducerBase.

Definition at line 272 of file BPHWriteSpecificDecay.cc.

272 { return; }
void BPHWriteSpecificDecay::endJob ( void  )
overridevirtual

Reimplemented from edm::one::EDProducerBase.

Definition at line 1283 of file BPHWriteSpecificDecay.cc.

1283 { return; }
void BPHWriteSpecificDecay::fill ( edm::Event ev,
const edm::EventSetup es 
)
virtual

Definition at line 303 of file BPHWriteSpecificDecay.cc.

References PVValHelper::add(), b0, reco::Candidate::begin(), cms::cuda::bs, BPHDecayToResResBuilder::build(), BPHX3872ToJPsiPiPiBuilder::build(), BPHDecayToResFlyingBuilder::build(), BPHDecayToChargedXXbarBuilder::build(), BPHOniaToMuMuBuilder::build(), BPHDecayToResTrkBuilder::build(), BPHDecayToTkpTknSymChargeBuilder::build(), BPHDecayToV0Builder::build(), TwoTrackMinimumDistance::calculate(), edm::HandleBase::clear(), BPHPlusMinusCandidate::composite(), BPHRecoBuilder::createCollection(), BPHDecayToResFlyingBuilder::daughMap(), reco::CompositeCandidate::daughter(), BPHDecayMomentum::daughters(), alignCSCRings::e, Puppi_cff::etaMax, edm::EventSetup::get(), BPHOniaToMuMuBuilder::getConstrMass(), BPHOniaToMuMuBuilder::getConstrSigma(), BPHDecayMomentum::getDaug(), BPHOniaToMuMuBuilder::getList(), mps_fire::i, dqmiolumiharvest::j, visualization-live-secondInstance_cfg::m, HLT_FULL_cff::magneticField, dqmiodumpmetadata::n, HLT_FULL_cff::nPhi, reco::CompositeCandidate::numberOfDaughters(), BPHDecayMomentum::originalReco(), BPHOniaToMuMuBuilder::Phi, TwoTrackMinimumDistance::points(), reco::Vertex::position(), funct::pow(), edm::Handle< T >::product(), BPHOniaToMuMuBuilder::Psi1, BPHOniaToMuMuBuilder::Psi2, ptMin, MetAnalyzer::pv(), BPHRecoBuilder::sameTrack(), edm::second(), BPHDecayConstrainedBuilder::setConstr(), BPHX3872ToJPsiPiPiBuilder::setConstr(), BPHOniaToMuMuBuilder::setConstr(), BPHDecayToChargedXXbarBuilder::setEtaMax(), BPHOniaToMuMuBuilder::setEtaMax(), BPHKx0ToKPiBuilder::setEtaMax(), BPHDecayToV0Builder::setEtaMax(), BPHX3872ToJPsiPiPiBuilder::setJPsiMassMax(), BPHBsToJPsiPhiBuilder::setJPsiMassMax(), BPHBdToJPsiKxBuilder::setJPsiMassMax(), BPHBdToJPsiKsBuilder::setJPsiMassMax(), BPHLbToJPsiL0Builder::setJPsiMassMax(), BPHBuToJPsiKBuilder::setJPsiMassMax(), BPHBcToJPsiPiBuilder::setJPsiMassMax(), BPHX3872ToJPsiPiPiBuilder::setJPsiMassMin(), BPHBsToJPsiPhiBuilder::setJPsiMassMin(), BPHBdToJPsiKxBuilder::setJPsiMassMin(), BPHBdToJPsiKsBuilder::setJPsiMassMin(), BPHLbToJPsiL0Builder::setJPsiMassMin(), BPHBuToJPsiKBuilder::setJPsiMassMin(), BPHBcToJPsiPiBuilder::setJPsiMassMin(), BPHBdToJPsiKsBuilder::setK0MassMax(), BPHBdToJPsiKsBuilder::setK0MassMin(), BPHBuToJPsiKBuilder::setKEtaMax(), BPHBuToJPsiKBuilder::setKPtMin(), BPHBdToJPsiKxBuilder::setKxMassMax(), BPHBdToJPsiKxBuilder::setKxMassMin(), BPHLbToJPsiL0Builder::setLambda0MassMax(), BPHLbToJPsiL0Builder::setLambda0MassMin(), BPHDecayGenericBuilder::setMassFitMax(), BPHX3872ToJPsiPiPiBuilder::setMassFitMax(), BPHDecayGenericBuilder::setMassFitMin(), BPHX3872ToJPsiPiPiBuilder::setMassFitMin(), BPHDecayGenericBuilder::setMassMax(), BPHX3872ToJPsiPiPiBuilder::setMassMax(), BPHOniaToMuMuBuilder::setMassMax(), BPHDecayGenericBuilder::setMassMin(), BPHX3872ToJPsiPiPiBuilder::setMassMin(), BPHOniaToMuMuBuilder::setMassMin(), BPHBsToJPsiPhiBuilder::setPhiMassMax(), BPHBsToJPsiPhiBuilder::setPhiMassMin(), BPHX3872ToJPsiPiPiBuilder::setPiEtaMax(), BPHBcToJPsiPiBuilder::setPiEtaMax(), BPHX3872ToJPsiPiPiBuilder::setPiPtMin(), BPHBcToJPsiPiBuilder::setPiPtMin(), BPHDecayGenericBuilder::setProbMin(), BPHX3872ToJPsiPiPiBuilder::setProbMin(), BPHOniaToMuMuBuilder::setProbMin(), BPHDecayToChargedXXbarBuilder::setPtMin(), BPHKx0ToKPiBuilder::setPtMin(), BPHOniaToMuMuBuilder::setPtMin(), BPHDecayToV0Builder::setPtMin(), BPHOniaToMuMuBuilder::Ups, BPHOniaToMuMuBuilder::Ups1, BPHOniaToMuMuBuilder::Ups2, BPHOniaToMuMuBuilder::Ups3, pat::PATObject< ObjectType >::userData(), and BPHDecayVertex::vertex().

303  {
304  lFull.clear();
305  lJPsi.clear();
306  lSd.clear();
307  lSs.clear();
308  lBu.clear();
309  lBd.clear();
310  lBs.clear();
311  lK0.clear();
312  lL0.clear();
313  lB0.clear();
314  lLb.clear();
315  lBc.clear();
316  lX3872.clear();
317  jPsiOMap.clear();
318  daughMap.clear();
319  pvRefMap.clear();
320  ccRefMap.clear();
321 
322  // get magnetic field
324  es.get<IdealMagneticFieldRecord>().get(magneticField);
325 
326  // get object collections
327  // collections are got through "BPHTokenWrapper" interface to allow
328  // uniform access in different CMSSW versions
329 
331  pVertexToken.get(ev, pVertices);
332  int npv = pVertices->size();
333 
334  int nrc = 0;
335 
336  // get reco::PFCandidate collection (in full AOD )
338  if (usePF) {
339  pfCandsToken.get(ev, pfCands);
340  nrc = pfCands->size();
341  }
342 
343  // get pat::PackedCandidate collection (in MiniAOD)
344  // pat::PackedCandidate is not defined in CMSSW_5XY, so a
345  // typedef (BPHTrackReference::candidate) is used, actually referring
346  // to pat::PackedCandidate only for CMSSW versions where it's defined
348  if (usePC) {
349  pcCandsToken.get(ev, pcCands);
350  nrc = pcCands->size();
351  }
352 
353  // get pat::GenericParticle collection (in skimmed data)
355  if (useGP) {
356  gpCandsToken.get(ev, gpCands);
357  nrc = gpCands->size();
358  }
359 
360  // get pat::Muon collection (in full AOD and MiniAOD)
362  if (usePM) {
363  patMuonToken.get(ev, patMuon);
364  }
365 
366  // get K0 reco::VertexCompositeCandidate collection (in full AOD)
368  if (useK0) {
369  k0CandsToken.get(ev, k0Cand);
370  }
371 
372  // get Lambda0 reco::VertexCompositeCandidate collection (in full AOD)
374  if (useL0) {
375  l0CandsToken.get(ev, l0Cand);
376  }
377 
378  // get K0 reco::VertexCompositePtrCandidate collection (in MiniAOD)
380  if (useKS) {
381  kSCandsToken.get(ev, kSCand);
382  }
383 
384  // get Lambda0 reco::VertexCompositePtrCandidate collection (in MiniAOD)
386  if (useLS) {
387  lSCandsToken.get(ev, lSCand);
388  }
389 
390  // get muons from pat::CompositeCandidate objects describing onia;
391  // muons from all composite objects are copied to an unique std::vector
392  vector<const reco::Candidate*> muDaugs;
393  set<const pat::Muon*> muonSet;
394  typedef multimap<const reco::Candidate*, const pat::CompositeCandidate*> mu_cc_map;
395  mu_cc_map muCCMap;
396  if (useCC) {
398  ccCandsToken.get(ev, ccCands);
399  int n = ccCands->size();
400  muDaugs.clear();
401  muDaugs.reserve(n);
402  muonSet.clear();
403  set<const pat::Muon*>::const_iterator iter;
404  set<const pat::Muon*>::const_iterator iend;
405  int i;
406  for (i = 0; i < n; ++i) {
407  const pat::CompositeCandidate& cc = ccCands->at(i);
408  int j;
409  int m = cc.numberOfDaughters();
410  for (j = 0; j < m; ++j) {
411  const reco::Candidate* dp = cc.daughter(j);
412  const pat::Muon* mp = dynamic_cast<const pat::Muon*>(dp);
413  iter = muonSet.begin();
414  iend = muonSet.end();
415  bool add = (mp != nullptr) && (muonSet.find(mp) == iend);
416  while (add && (iter != iend)) {
417  if (BPHRecoBuilder::sameTrack(mp, *iter++, 1.0e-5))
418  add = false;
419  }
420  if (add)
421  muonSet.insert(mp);
422  // associate muon to the CompositeCandidate containing it
423  muCCMap.insert(pair<const reco::Candidate*, const pat::CompositeCandidate*>(dp, &cc));
424  }
425  }
426  iter = muonSet.begin();
427  iend = muonSet.end();
428  while (iter != iend)
429  muDaugs.push_back(*iter++);
430  }
431 
432  map<recoType, map<parType, double>>::const_iterator rIter = parMap.begin();
433  map<recoType, map<parType, double>>::const_iterator rIend = parMap.end();
434 
435  // reconstruct quarkonia
436 
437  BPHOniaToMuMuBuilder* onia = nullptr;
438  if (recoOnia) {
439  if (usePM)
440  onia = new BPHOniaToMuMuBuilder(
441  es, BPHRecoBuilder::createCollection(patMuon, "ingmcf"), BPHRecoBuilder::createCollection(patMuon, "ingmcf"));
442  else if (useCC)
443  onia = new BPHOniaToMuMuBuilder(
444  es, BPHRecoBuilder::createCollection(muDaugs, "ingmcf"), BPHRecoBuilder::createCollection(muDaugs, "ingmcf"));
445  }
446 
447  if (onia != nullptr) {
448  while (rIter != rIend) {
449  const map<recoType, map<parType, double>>::value_type& rEntry = *rIter++;
450  recoType rType = rEntry.first;
451  const map<parType, double>& pMap = rEntry.second;
453  switch (rType) {
454  case Pmm:
456  break;
457  case Psi1:
459  break;
460  case Psi2:
462  break;
463  case Ups:
465  break;
466  case Ups1:
468  break;
469  case Ups2:
471  break;
472  case Ups3:
474  break;
475  default:
476  continue;
477  }
478  map<parType, double>::const_iterator pIter = pMap.begin();
479  map<parType, double>::const_iterator pIend = pMap.end();
480  while (pIter != pIend) {
481  const map<parType, double>::value_type& pEntry = *pIter++;
482  parType id = pEntry.first;
483  double pv = pEntry.second;
484  switch (id) {
485  case ptMin:
486  onia->setPtMin(type, pv);
487  break;
488  case etaMax:
489  onia->setEtaMax(type, pv);
490  break;
491  case massMin:
492  onia->setMassMin(type, pv);
493  break;
494  case massMax:
495  onia->setMassMax(type, pv);
496  break;
497  case probMin:
498  onia->setProbMin(type, pv);
499  break;
500  case constrMass:
501  onia->setConstr(type, pv, onia->getConstrSigma(type));
502  break;
503  case constrSigma:
504  onia->setConstr(type, onia->getConstrMass(type), pv);
505  break;
506  default:
507  break;
508  }
509  }
510  }
511  lFull = onia->build();
512  }
513 
514  // associate onia to primary vertex
515 
516  int iFull;
517  int nFull = lFull.size();
518  map<const BPHRecoCandidate*, const reco::Vertex*> oniaVtxMap;
519 
520  typedef mu_cc_map::const_iterator mu_cc_iter;
521  for (iFull = 0; iFull < nFull; ++iFull) {
522  const reco::Vertex* pVtx = nullptr;
523  int pvId = 0;
524  const BPHPlusMinusCandidate* ptr = lFull[iFull].get();
525  const std::vector<const reco::Candidate*>& daugs = ptr->daughters();
526 
527  // try to recover primary vertex association in skim data:
528  // get the CompositeCandidate containing both muons
529  pair<mu_cc_iter, mu_cc_iter> cc0 = muCCMap.equal_range(ptr->originalReco(daugs[0]));
530  pair<mu_cc_iter, mu_cc_iter> cc1 = muCCMap.equal_range(ptr->originalReco(daugs[1]));
531  mu_cc_iter iter0 = cc0.first;
532  mu_cc_iter iend0 = cc0.second;
533  mu_cc_iter iter1 = cc1.first;
534  mu_cc_iter iend1 = cc1.second;
535  while ((iter0 != iend0) && (pVtx == nullptr)) {
536  const pat::CompositeCandidate* ccp = iter0++->second;
537  while (iter1 != iend1) {
538  if (ccp != iter1++->second)
539  continue;
540  pVtx = ccp->userData<reco::Vertex>("PVwithmuons");
541  const reco::Vertex* sVtx = nullptr;
542  const reco::Vertex::Point& pPos = pVtx->position();
543  float dMin = 999999.;
544  int ipv;
545  for (ipv = 0; ipv < npv; ++ipv) {
546  const reco::Vertex* tVtx = &pVertices->at(ipv);
547  const reco::Vertex::Point& tPos = tVtx->position();
548  float dist = pow(pPos.x() - tPos.x(), 2) + pow(pPos.y() - tPos.y(), 2) + pow(pPos.z() - tPos.z(), 2);
549  if (dist < dMin) {
550  dMin = dist;
551  sVtx = tVtx;
552  pvId = ipv;
553  }
554  }
555  pVtx = sVtx;
556  break;
557  }
558  }
559 
560  // if not found, as for other type of input data,
561  // try to get the nearest primary vertex in z direction
562  if (pVtx == nullptr) {
563  const reco::Vertex::Point& sVtp = ptr->vertex().position();
564  GlobalPoint cPos(sVtp.x(), sVtp.y(), sVtp.z());
565  const pat::CompositeCandidate& sCC = ptr->composite();
566  GlobalVector cDir(sCC.px(), sCC.py(), sCC.pz());
567  GlobalPoint bPos(0.0, 0.0, 0.0);
568  GlobalVector bDir(0.0, 0.0, 1.0);
570  bool state = ttmd.calculate(GlobalTrajectoryParameters(cPos, cDir, TrackCharge(0), &(*magneticField)),
571  GlobalTrajectoryParameters(bPos, bDir, TrackCharge(0), &(*magneticField)));
572  float minDz = 999999.;
573  float extrapZ = (state ? ttmd.points().first.z() : -9e20);
574  int ipv;
575  for (ipv = 0; ipv < npv; ++ipv) {
576  const reco::Vertex& tVtx = pVertices->at(ipv);
577  float deltaZ = fabs(extrapZ - tVtx.position().z());
578  if (deltaZ < minDz) {
579  minDz = deltaZ;
580  pVtx = &tVtx;
581  pvId = ipv;
582  }
583  }
584  }
585 
586  oniaVtxMap[ptr] = pVtx;
587  pvRefMap[ptr] = vertex_ref(pVertices, pvId);
588  }
589  pVertexToken.get(ev, pVertices);
590 
591  // get JPsi subsample and associate JPsi candidate to original
592  // generic onia candidate
593  if (nFull)
595 
596  int nJPsi = lJPsi.size();
597  delete onia;
598 
599  if (!nJPsi)
600  return;
601  if (!nrc)
602  return;
603 
604  int ij;
605  int io;
606  int nj = lJPsi.size();
607  int no = lFull.size();
608  for (ij = 0; ij < nj; ++ij) {
609  const BPHRecoCandidate* jp = lJPsi[ij].get();
610  for (io = 0; io < no; ++io) {
611  const BPHRecoCandidate* oc = lFull[io].get();
612  if ((jp->originalReco(jp->getDaug("MuPos")) == oc->originalReco(oc->getDaug("MuPos"))) &&
613  (jp->originalReco(jp->getDaug("MuNeg")) == oc->originalReco(oc->getDaug("MuNeg")))) {
614  jPsiOMap[jp] = oc;
615  break;
616  }
617  }
618  }
619 
620  // build and dump Bu
621 
622  BPHBuToJPsiKBuilder* bu = nullptr;
623  if (recoBu) {
624  if (usePF)
626  else if (usePC)
628  else if (useGP)
630  }
631 
632  if (bu != nullptr) {
633  rIter = parMap.find(Bu);
634  if (rIter != rIend) {
635  const map<parType, double>& pMap = rIter->second;
636  map<parType, double>::const_iterator pIter = pMap.begin();
637  map<parType, double>::const_iterator pIend = pMap.end();
638  while (pIter != pIend) {
639  const map<parType, double>::value_type& pEntry = *pIter++;
640  parType id = pEntry.first;
641  double pv = pEntry.second;
642  switch (id) {
643  case ptMin:
644  bu->setKPtMin(pv);
645  break;
646  case etaMax:
647  bu->setKEtaMax(pv);
648  break;
649  case mPsiMin:
650  bu->setJPsiMassMin(pv);
651  break;
652  case mPsiMax:
653  bu->setJPsiMassMax(pv);
654  break;
655  case massMin:
656  bu->setMassMin(pv);
657  break;
658  case massMax:
659  bu->setMassMax(pv);
660  break;
661  case probMin:
662  bu->setProbMin(pv);
663  break;
664  case mFitMin:
665  bu->setMassFitMin(pv);
666  break;
667  case mFitMax:
668  bu->setMassFitMax(pv);
669  break;
670  case constrMJPsi:
671  bu->setConstr(pv > 0);
672  break;
673  case writeCandidate:
674  writeBu = (pv > 0);
675  break;
676  default:
677  break;
678  }
679  }
680  }
681  lBu = bu->build();
682  delete bu;
683  }
684 
685  // build and dump Kx0
686 
687  vector<BPHPlusMinusConstCandPtr> lKx0;
688  BPHKx0ToKPiBuilder* kx0 = nullptr;
689  if (recoKx0) {
690  if (usePF)
691  kx0 = new BPHKx0ToKPiBuilder(
693  else if (usePC)
694  kx0 = new BPHKx0ToKPiBuilder(
696  else if (useGP)
697  kx0 = new BPHKx0ToKPiBuilder(
699  }
700 
701  if (kx0 != nullptr) {
702  rIter = parMap.find(Kx0);
703  if (rIter != rIend) {
704  const map<parType, double>& pMap = rIter->second;
705  map<parType, double>::const_iterator pIter = pMap.begin();
706  map<parType, double>::const_iterator pIend = pMap.end();
707  while (pIter != pIend) {
708  const map<parType, double>::value_type& pEntry = *pIter++;
709  parType id = pEntry.first;
710  double pv = pEntry.second;
711  switch (id) {
712  case ptMin:
713  kx0->setPtMin(pv);
714  break;
715  case etaMax:
716  kx0->setEtaMax(pv);
717  break;
718  case massMin:
719  kx0->setMassMin(pv);
720  break;
721  case massMax:
722  kx0->setMassMax(pv);
723  break;
724  case probMin:
725  kx0->setProbMin(pv);
726  break;
727  case writeCandidate:
728  writeKx0 = (pv > 0);
729  break;
730  default:
731  break;
732  }
733  }
734  }
735  lKx0 = kx0->build();
736  delete kx0;
737  }
738 
739  int nKx0 = lKx0.size();
740 
741  // build and dump Bd -> JPsi Kx0
742 
743  if (recoBd && nKx0) {
745  rIter = parMap.find(Bd);
746  if (rIter != rIend) {
747  const map<parType, double>& pMap = rIter->second;
748  map<parType, double>::const_iterator pIter = pMap.begin();
749  map<parType, double>::const_iterator pIend = pMap.end();
750  while (pIter != pIend) {
751  const map<parType, double>::value_type& pEntry = *pIter++;
752  parType id = pEntry.first;
753  double pv = pEntry.second;
754  switch (id) {
755  case mPsiMin:
756  bd->setJPsiMassMin(pv);
757  break;
758  case mPsiMax:
759  bd->setJPsiMassMax(pv);
760  break;
761  case mKx0Min:
762  bd->setKxMassMin(pv);
763  break;
764  case mKx0Max:
765  bd->setKxMassMax(pv);
766  break;
767  case massMin:
768  bd->setMassMin(pv);
769  break;
770  case massMax:
771  bd->setMassMax(pv);
772  break;
773  case probMin:
774  bd->setProbMin(pv);
775  break;
776  case mFitMin:
777  bd->setMassFitMin(pv);
778  break;
779  case mFitMax:
780  bd->setMassFitMax(pv);
781  break;
782  case constrMJPsi:
783  bd->setConstr(pv > 0);
784  break;
785  case writeCandidate:
786  writeBd = (pv > 0);
787  break;
788  default:
789  break;
790  }
791  }
792  }
793 
794  lBd = bd->build();
795  delete bd;
796 
797  set<BPHRecoConstCandPtr> sKx0;
798  int iBd;
799  int nBd = lBd.size();
800  for (iBd = 0; iBd < nBd; ++iBd)
801  sKx0.insert(lBd[iBd]->getComp("Kx0"));
802  set<BPHRecoConstCandPtr>::const_iterator iter = sKx0.begin();
803  set<BPHRecoConstCandPtr>::const_iterator iend = sKx0.end();
804  while (iter != iend)
805  lSd.push_back(*iter++);
806  }
807 
808  // build and dump Phi
809 
810  vector<BPHPlusMinusConstCandPtr> lPhi;
811  BPHPhiToKKBuilder* phi = nullptr;
812  if (recoPkk) {
813  if (usePF)
814  phi = new BPHPhiToKKBuilder(
816  else if (usePC)
817  phi = new BPHPhiToKKBuilder(
819  else if (useGP)
820  phi = new BPHPhiToKKBuilder(
822  }
823 
824  if (phi != nullptr) {
825  rIter = parMap.find(Pkk);
826  if (rIter != rIend) {
827  const map<parType, double>& pMap = rIter->second;
828  map<parType, double>::const_iterator pIter = pMap.begin();
829  map<parType, double>::const_iterator pIend = pMap.end();
830  while (pIter != pIend) {
831  const map<parType, double>::value_type& pEntry = *pIter++;
832  parType id = pEntry.first;
833  double pv = pEntry.second;
834  switch (id) {
835  case ptMin:
836  phi->setPtMin(pv);
837  break;
838  case etaMax:
839  phi->setEtaMax(pv);
840  break;
841  case massMin:
842  phi->setMassMin(pv);
843  break;
844  case massMax:
845  phi->setMassMax(pv);
846  break;
847  case probMin:
848  phi->setProbMin(pv);
849  break;
850  case writeCandidate:
851  writePkk = (pv > 0);
852  break;
853  default:
854  break;
855  }
856  }
857  }
858  lPhi = phi->build();
859  delete phi;
860  }
861 
862  int nPhi = lPhi.size();
863 
864  // build and dump Bs
865 
866  if (recoBs && nPhi) {
868  rIter = parMap.find(Bs);
869  if (rIter != rIend) {
870  const map<parType, double>& pMap = rIter->second;
871  map<parType, double>::const_iterator pIter = pMap.begin();
872  map<parType, double>::const_iterator pIend = pMap.end();
873  while (pIter != pIend) {
874  const map<parType, double>::value_type& pEntry = *pIter++;
875  parType id = pEntry.first;
876  double pv = pEntry.second;
877  switch (id) {
878  case mPsiMin:
879  bs->setJPsiMassMin(pv);
880  break;
881  case mPsiMax:
882  bs->setJPsiMassMax(pv);
883  break;
884  case mPhiMin:
885  bs->setPhiMassMin(pv);
886  break;
887  case mPhiMax:
888  bs->setPhiMassMax(pv);
889  break;
890  case massMin:
891  bs->setMassMin(pv);
892  break;
893  case massMax:
894  bs->setMassMax(pv);
895  break;
896  case probMin:
897  bs->setProbMin(pv);
898  break;
899  case mFitMin:
900  bs->setMassFitMin(pv);
901  break;
902  case mFitMax:
903  bs->setMassFitMax(pv);
904  break;
905  case constrMJPsi:
906  bs->setConstr(pv > 0);
907  break;
908  case writeCandidate:
909  writeBs = (pv > 0);
910  break;
911  default:
912  break;
913  }
914  }
915  }
916 
917  lBs = bs->build();
918  delete bs;
919 
920  set<BPHRecoConstCandPtr> sPhi;
921  int iBs;
922  int nBs = lBs.size();
923  for (iBs = 0; iBs < nBs; ++iBs)
924  sPhi.insert(lBs[iBs]->getComp("Phi"));
925  set<BPHRecoConstCandPtr>::const_iterator iter = sPhi.begin();
926  set<BPHRecoConstCandPtr>::const_iterator iend = sPhi.end();
927  while (iter != iend)
928  lSs.push_back(*iter++);
929  }
930 
931  // build K0
932 
933  BPHK0sToPiPiBuilder* k0s = nullptr;
934  if (recoK0s) {
935  if (useK0)
936  k0s = new BPHK0sToPiPiBuilder(es, k0Cand.product(), "cfp");
937  else if (useKS)
938  k0s = new BPHK0sToPiPiBuilder(es, kSCand.product(), "cfp");
939  }
940  if (k0s != nullptr) {
941  rIter = parMap.find(K0s);
942  if (rIter != rIend) {
943  const map<parType, double>& pMap = rIter->second;
944  map<parType, double>::const_iterator pIter = pMap.begin();
945  map<parType, double>::const_iterator pIend = pMap.end();
946  while (pIter != pIend) {
947  const map<parType, double>::value_type& pEntry = *pIter++;
948  parType id = pEntry.first;
949  double pv = pEntry.second;
950  switch (id) {
951  case ptMin:
952  k0s->setPtMin(pv);
953  break;
954  case etaMax:
955  k0s->setEtaMax(pv);
956  break;
957  case massMin:
958  k0s->setMassMin(pv);
959  break;
960  case massMax:
961  k0s->setMassMax(pv);
962  break;
963  case probMin:
964  k0s->setProbMin(pv);
965  break;
966  case writeCandidate:
967  writeK0s = (pv > 0);
968  break;
969  default:
970  break;
971  }
972  }
973  }
974  lK0 = k0s->build();
975  delete k0s;
976  }
977 
978  int nK0 = lK0.size();
979 
980  // build Lambda0
981 
982  BPHLambda0ToPPiBuilder* l0s = nullptr;
983  if (recoLambda0) {
984  if (useL0)
985  l0s = new BPHLambda0ToPPiBuilder(es, l0Cand.product(), "cfp");
986  else if (useLS)
987  l0s = new BPHLambda0ToPPiBuilder(es, lSCand.product(), "cfp");
988  }
989  if (l0s != nullptr) {
990  rIter = parMap.find(Lambda0);
991  if (rIter != rIend) {
992  const map<parType, double>& pMap = rIter->second;
993  map<parType, double>::const_iterator pIter = pMap.begin();
994  map<parType, double>::const_iterator pIend = pMap.end();
995  while (pIter != pIend) {
996  const map<parType, double>::value_type& pEntry = *pIter++;
997  parType id = pEntry.first;
998  double pv = pEntry.second;
999  switch (id) {
1000  case ptMin:
1001  l0s->setPtMin(pv);
1002  break;
1003  case etaMax:
1004  l0s->setEtaMax(pv);
1005  break;
1006  case massMin:
1007  l0s->setMassMin(pv);
1008  break;
1009  case massMax:
1010  l0s->setMassMax(pv);
1011  break;
1012  case probMin:
1013  l0s->setProbMin(pv);
1014  break;
1015  case writeCandidate:
1016  writeLambda0 = (pv > 0);
1017  break;
1018  default:
1019  break;
1020  }
1021  }
1022  }
1023  lL0 = l0s->build();
1024  delete l0s;
1025  }
1026 
1027  int nL0 = lL0.size();
1028 
1029  // build and dump Bd -> JPsi K0s
1030 
1031  if (recoB0 && nK0) {
1033  rIter = parMap.find(B0);
1034  if (rIter != rIend) {
1035  const map<parType, double>& pMap = rIter->second;
1036  map<parType, double>::const_iterator pIter = pMap.begin();
1037  map<parType, double>::const_iterator pIend = pMap.end();
1038  while (pIter != pIend) {
1039  const map<parType, double>::value_type& pEntry = *pIter++;
1040  parType id = pEntry.first;
1041  double pv = pEntry.second;
1042  switch (id) {
1043  case mPsiMin:
1044  b0->setJPsiMassMin(pv);
1045  break;
1046  case mPsiMax:
1047  b0->setJPsiMassMax(pv);
1048  break;
1049  case mK0sMin:
1050  b0->setK0MassMin(pv);
1051  break;
1052  case mK0sMax:
1053  b0->setK0MassMax(pv);
1054  break;
1055  case massMin:
1056  b0->setMassMin(pv);
1057  break;
1058  case massMax:
1059  b0->setMassMax(pv);
1060  break;
1061  case probMin:
1062  b0->setProbMin(pv);
1063  break;
1064  case mFitMin:
1065  b0->setMassFitMin(pv);
1066  break;
1067  case mFitMax:
1068  b0->setMassFitMax(pv);
1069  break;
1070  case constrMJPsi:
1071  b0->setConstr(pv > 0);
1072  break;
1073  case writeCandidate:
1074  writeB0 = (pv > 0);
1075  break;
1076  default:
1077  break;
1078  }
1079  }
1080  }
1081 
1082  lB0 = b0->build();
1083  const map<const BPHRecoCandidate*, const BPHRecoCandidate*>& b0Map = b0->daughMap();
1084  daughMap.insert(b0Map.begin(), b0Map.end());
1085  delete b0;
1086  }
1087 
1088  // build and dump Lambdab -> JPsi Lambda0
1089 
1090  if (recoLambdab && nL0) {
1092  rIter = parMap.find(Lambdab);
1093  if (rIter != rIend) {
1094  const map<parType, double>& pMap = rIter->second;
1095  map<parType, double>::const_iterator pIter = pMap.begin();
1096  map<parType, double>::const_iterator pIend = pMap.end();
1097  while (pIter != pIend) {
1098  const map<parType, double>::value_type& pEntry = *pIter++;
1099  parType id = pEntry.first;
1100  double pv = pEntry.second;
1101  switch (id) {
1102  case mPsiMin:
1103  lb->setJPsiMassMin(pv);
1104  break;
1105  case mPsiMax:
1106  lb->setJPsiMassMax(pv);
1107  break;
1108  case mLambda0Min:
1109  lb->setLambda0MassMin(pv);
1110  break;
1111  case mLambda0Max:
1112  lb->setLambda0MassMax(pv);
1113  break;
1114  case massMin:
1115  lb->setMassMin(pv);
1116  break;
1117  case massMax:
1118  lb->setMassMax(pv);
1119  break;
1120  case probMin:
1121  lb->setProbMin(pv);
1122  break;
1123  case mFitMin:
1124  lb->setMassFitMin(pv);
1125  break;
1126  case mFitMax:
1127  lb->setMassFitMax(pv);
1128  break;
1129  case constrMJPsi:
1130  lb->setConstr(pv > 0);
1131  break;
1132  case writeCandidate:
1133  writeLambdab = (pv > 0);
1134  break;
1135  default:
1136  break;
1137  }
1138  }
1139  }
1140 
1141  lLb = lb->build();
1142  const map<const BPHRecoCandidate*, const BPHRecoCandidate*>& ldMap = lb->daughMap();
1143  daughMap.insert(ldMap.begin(), ldMap.end());
1144  delete lb;
1145  }
1146 
1147  // build and dump Bc
1148 
1149  BPHBcToJPsiPiBuilder* bc = nullptr;
1150  if (recoBc) {
1151  if (usePF)
1153  else if (usePC)
1155  else if (useGP)
1157  }
1158 
1159  if (bc != nullptr) {
1160  rIter = parMap.find(Bc);
1161  if (rIter != rIend) {
1162  const map<parType, double>& pMap = rIter->second;
1163  map<parType, double>::const_iterator pIter = pMap.begin();
1164  map<parType, double>::const_iterator pIend = pMap.end();
1165  while (pIter != pIend) {
1166  const map<parType, double>::value_type& pEntry = *pIter++;
1167  parType id = pEntry.first;
1168  double pv = pEntry.second;
1169  switch (id) {
1170  case ptMin:
1171  bc->setPiPtMin(pv);
1172  break;
1173  case etaMax:
1174  bc->setPiEtaMax(pv);
1175  break;
1176  case mPsiMin:
1177  bc->setJPsiMassMin(pv);
1178  break;
1179  case mPsiMax:
1180  bc->setJPsiMassMax(pv);
1181  break;
1182  case massMin:
1183  bc->setMassMin(pv);
1184  break;
1185  case massMax:
1186  bc->setMassMax(pv);
1187  break;
1188  case probMin:
1189  bc->setProbMin(pv);
1190  break;
1191  case mFitMin:
1192  bc->setMassFitMin(pv);
1193  break;
1194  case mFitMax:
1195  bc->setMassFitMax(pv);
1196  break;
1197  case constrMJPsi:
1198  bc->setConstr(pv > 0);
1199  break;
1200  case writeCandidate:
1201  writeBc = (pv > 0);
1202  break;
1203  default:
1204  break;
1205  }
1206  }
1207  }
1208  lBc = bc->build();
1209  delete bc;
1210  }
1211 
1212  // build and dump X3872
1213 
1214  BPHX3872ToJPsiPiPiBuilder* x3872 = nullptr;
1215  if (recoX3872) {
1216  if (usePF)
1217  x3872 = new BPHX3872ToJPsiPiPiBuilder(
1219  else if (usePC)
1220  x3872 = new BPHX3872ToJPsiPiPiBuilder(
1222  else if (useGP)
1223  x3872 = new BPHX3872ToJPsiPiPiBuilder(
1225  }
1226 
1227  if (x3872 != nullptr) {
1228  rIter = parMap.find(X3872);
1229  if (rIter != rIend) {
1230  const map<parType, double>& pMap = rIter->second;
1231  map<parType, double>::const_iterator pIter = pMap.begin();
1232  map<parType, double>::const_iterator pIend = pMap.end();
1233  while (pIter != pIend) {
1234  const map<parType, double>::value_type& pEntry = *pIter++;
1235  parType id = pEntry.first;
1236  double pv = pEntry.second;
1237  switch (id) {
1238  case ptMin:
1239  x3872->setPiPtMin(pv);
1240  break;
1241  case etaMax:
1242  x3872->setPiEtaMax(pv);
1243  break;
1244  case mPsiMin:
1245  x3872->setJPsiMassMin(pv);
1246  break;
1247  case mPsiMax:
1248  x3872->setJPsiMassMax(pv);
1249  break;
1250  case massMin:
1251  x3872->setMassMin(pv);
1252  break;
1253  case massMax:
1254  x3872->setMassMax(pv);
1255  break;
1256  case probMin:
1257  x3872->setProbMin(pv);
1258  break;
1259  case mFitMin:
1260  x3872->setMassFitMin(pv);
1261  break;
1262  case mFitMax:
1263  x3872->setMassFitMax(pv);
1264  break;
1265  case constrMJPsi:
1266  x3872->setConstr(pv > 0);
1267  break;
1268  case writeCandidate:
1269  writeX3872 = (pv > 0);
1270  break;
1271  default:
1272  break;
1273  }
1274  }
1275  }
1276  lX3872 = x3872->build();
1277  delete x3872;
1278  }
1279 
1280  return;
1281 }
BPHTokenWrapper< std::vector< reco::VertexCompositeCandidate > > l0CandsToken
Analysis-level particle class.
static bool sameTrack(const reco::Candidate *lCand, const reco::Candidate *rCand, double minPDifference)
void setJPsiMassMax(double m)
const Candidate * daughter(size_type) const override
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
BPHTokenWrapper< std::vector< reco::Vertex > > pVertexToken
std::vector< BPHRecoConstCandPtr > build()
build candidates
std::vector< BPHPlusMinusConstCandPtr > lJPsi
void setJPsiMassMin(double m)
set cuts
void setJPsiMassMin(double m)
set cuts
BPHTokenWrapper< std::vector< reco::VertexCompositePtrCandidate > > kSCandsToken
std::vector< BPHPlusMinusConstCandPtr > lK0
std::vector< BPHPlusMinusConstCandPtr > build()
build candidates
edm::Ref< std::vector< reco::Vertex > > vertex_ref
void setMassMax(oniaType type, double m)
std::map< const BPHRecoCandidate *, compcc_ref > ccRefMap
std::vector< BPHPlusMinusConstCandPtr > lFull
double getConstrMass(oniaType type) const
void setJPsiMassMin(double m)
void setKPtMin(double pt)
set cuts
void setPiPtMin(double pt)
set cuts
std::vector< BPHRecoConstCandPtr > lBs
bool get(const edm::Event &ev, edm::Handle< Obj > &obj)
BPHTokenWrapper< std::vector< BPHTrackReference::candidate > > pcCandsToken
void setKEtaMax(double eta)
tuple magneticField
void setPtMin(double pt)
const Point & position() const
position
Definition: Vertex.h:127
virtual const std::vector< const reco::Candidate * > & daughters() const
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
static BPHGenericCollection * createCollection(const edm::Handle< T > &collection, const std::string &list="cfhpmig")
std::vector< BPHRecoConstCandPtr > lSs
std::vector< BPHRecoConstCandPtr > build()
build candidates
U second(std::pair< T, U > const &p)
deep_tau::DeepTauBase::BasicDiscriminator bd
Definition: DeepTauId.cc:1082
int TrackCharge
Definition: TrackCharge.h:4
void setJPsiMassMax(double m)
BPHTokenWrapper< std::vector< reco::VertexCompositeCandidate > > k0CandsToken
std::vector< BPHPlusMinusConstCandPtr > build()
build Phi candidates
void setEtaMax(double eta)
std::vector< BPHPlusMinusConstCandPtr > getList(oniaType type, BPHRecoSelect *dSel=nullptr, BPHMomentumSelect *mSel=nullptr, BPHVertexSelect *vSel=nullptr, BPHFitSelect *kSel=nullptr)
double getConstrSigma(oniaType type) const
BPHTokenWrapper< pat::MuonCollection > patMuonToken
void setEtaMax(double eta)
BPHTokenWrapper< std::vector< reco::PFCandidate > > pfCandsToken
void setJPsiMassMin(double m)
set cuts
BPHTokenWrapper< std::vector< reco::VertexCompositePtrCandidate > > lSCandsToken
void setMassMin(double m)
set cuts
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
BPHTokenWrapper< std::vector< pat::CompositeCandidate > > ccCandsToken
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
std::pair< GlobalPoint, GlobalPoint > points() const override
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
std::vector< BPHRecoConstCandPtr > lSd
virtual const reco::Candidate * originalReco(const reco::Candidate *daug) const
get the original particle from the clone
void add(std::map< std::string, TH1 * > &h, TH1 *hist)
void setLambda0MassMax(double m)
const std::map< const BPHRecoCandidate *, const BPHRecoCandidate * > & daughMap() const
get original daughters map
BPHTokenWrapper< std::vector< pat::GenericParticle > > gpCandsToken
T const * product() const
Definition: Handle.h:70
void setJPsiMassMin(double m)
set cuts
size_type numberOfDaughters() const override
number of daughters
void setPiEtaMax(double eta)
std::map< const BPHRecoCandidate *, vertex_ref > pvRefMap
std::vector< BPHPlusMinusConstCandPtr > lL0
std::vector< BPHRecoConstCandPtr > build()
build X3872 candidates
std::vector< BPHRecoConstCandPtr > lBu
std::vector< BPHRecoConstCandPtr > lLb
std::vector< BPHPlusMinusConstCandPtr > build()
build resonance candidates
void setJPsiMassMin(double m)
void setConstr(oniaType type, double mass, double sigma)
void setEtaMax(oniaType type, double eta)
T get() const
Definition: EventSetup.h:82
void setPtMin(double pt)
set cuts
std::map< recoType, std::map< parType, double > > parMap
const_iterator begin() const
first daughter const_iterator
Definition: Candidate.h:143
std::vector< BPHRecoConstCandPtr > build()
build candidates
std::map< const BPHRecoCandidate *, const BPHRecoCandidate * > jPsiOMap
static constexpr float b0
void setJPsiMassMin(double m)
set cuts
std::vector< BPHRecoConstCandPtr > lB0
virtual const reco::Candidate * getDaug(const std::string &name) const
void setMassMin(oniaType type, double m)
std::vector< BPHRecoConstCandPtr > lBc
std::vector< BPHRecoConstCandPtr > lX3872
std::vector< BPHPlusMinusConstCandPtr > build()
build candidates
std::vector< BPHRecoConstCandPtr > lBd
Analysis-level muon class.
Definition: Muon.h:51
void setJPsiMassMax(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)
void setPiPtMin(double pt)
set cuts
void BPHWriteSpecificDecay::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 217 of file BPHWriteSpecificDecay.cc.

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

217  {
219  desc.add<string>("pVertexLabel", "");
220  desc.add<string>("patMuonLabel", "");
221  desc.add<string>("ccCandsLabel", "");
222  desc.add<string>("pfCandsLabel", "");
223  desc.add<string>("pcCandsLabel", "");
224  desc.add<string>("gpCandsLabel", "");
225  desc.add<string>("k0CandsLabel", "");
226  desc.add<string>("l0CandsLabel", "");
227  desc.add<string>("kSCandsLabel", "");
228  desc.add<string>("lSCandsLabel", "");
229  desc.add<string>("oniaName", "oniaCand");
230  desc.add<string>("sdName", "kx0Cand");
231  desc.add<string>("ssName", "phiCand");
232  desc.add<string>("buName", "buFitted");
233  desc.add<string>("bdName", "bdFitted");
234  desc.add<string>("bsName", "bsFitted");
235  desc.add<string>("k0Name", "k0Fitted");
236  desc.add<string>("l0Name", "l0Fitted");
237  desc.add<string>("b0Name", "b0Fitted");
238  desc.add<string>("lbName", "lbFitted");
239  desc.add<string>("bcName", "bcFitted");
240  desc.add<string>("x3872Name", "x3872Fitted");
241  desc.add<bool>("writeVertex", true);
242  desc.add<bool>("writeMomentum", true);
244  dpar.add<string>("name");
245  dpar.add<double>("ptMin", -2.0e35);
246  dpar.add<double>("etaMax", -2.0e35);
247  dpar.add<double>("mJPsiMin", -2.0e35);
248  dpar.add<double>("mJPsiMax", -2.0e35);
249  dpar.add<double>("mKx0Min", -2.0e35);
250  dpar.add<double>("mKx0Max", -2.0e35);
251  dpar.add<double>("mPhiMin", -2.0e35);
252  dpar.add<double>("mPhiMax", -2.0e35);
253  dpar.add<double>("mK0sMin", -2.0e35);
254  dpar.add<double>("mK0sMax", -2.0e35);
255  dpar.add<double>("mLambda0Min", -2.0e35);
256  dpar.add<double>("mLambda0Max", -2.0e35);
257  dpar.add<double>("massMin", -2.0e35);
258  dpar.add<double>("massMax", -2.0e35);
259  dpar.add<double>("probMin", -2.0e35);
260  dpar.add<double>("massFitMin", -2.0e35);
261  dpar.add<double>("massFitMax", -2.0e35);
262  dpar.add<double>("constrMass", -2.0e35);
263  dpar.add<double>("constrSigma", -2.0e35);
264  dpar.add<bool>("constrMJPsi", true);
265  dpar.add<bool>("writeCandidate", true);
266  vector<edm::ParameterSet> rpar;
267  desc.addVPSet("recoSelect", dpar, rpar);
268  descriptions.add("bphWriteSpecificDecay", desc);
269  return;
270 }
ParameterDescriptionBase * addVPSet(U const &iLabel, ParameterSetDescription const &validator, std::vector< ParameterSet > const &defaults)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void BPHWriteSpecificDecay::produce ( edm::Event ev,
const edm::EventSetup es 
)
overridevirtual

Implements edm::one::EDProducerBase.

Definition at line 274 of file BPHWriteSpecificDecay.cc.

References PVValHelper::fill(), and TablePrint::write.

274  {
275  fill(ev, es);
276  if (writeOnia)
277  write(ev, lFull, oniaName);
278  if (writeKx0)
279  write(ev, lSd, sdName);
280  if (writePkk)
281  write(ev, lSs, ssName);
282  if (writeBu)
283  write(ev, lBu, buName);
284  if (writeBd)
285  write(ev, lBd, bdName);
286  if (writeBs)
287  write(ev, lBs, bsName);
288  if (writeK0s)
289  write(ev, lK0, k0Name);
290  if (writeLambda0)
291  write(ev, lL0, l0Name);
292  if (writeB0)
293  write(ev, lB0, b0Name);
294  if (writeLambdab)
295  write(ev, lLb, lbName);
296  if (writeBc)
297  write(ev, lBc, bcName);
298  if (writeX3872)
299  write(ev, lX3872, x3872Name);
300  return;
301 }
std::vector< BPHPlusMinusConstCandPtr > lK0
std::vector< BPHPlusMinusConstCandPtr > lFull
std::vector< BPHRecoConstCandPtr > lBs
virtual void fill(edm::Event &ev, const edm::EventSetup &es)
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)
std::vector< BPHPlusMinusConstCandPtr > lL0
std::vector< BPHRecoConstCandPtr > lBu
std::vector< BPHRecoConstCandPtr > lLb
std::vector< BPHRecoConstCandPtr > lB0
std::vector< BPHRecoConstCandPtr > lBc
std::vector< BPHRecoConstCandPtr > lX3872
std::vector< BPHRecoConstCandPtr > lBd
void BPHWriteSpecificDecay::setRecoParameters ( const edm::ParameterSet ps)
private

Definition at line 1285 of file BPHWriteSpecificDecay.cc.

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

1285  {
1286  const string& name = ps.getParameter<string>("name");
1287  bool writeCandidate = ps.getParameter<bool>("writeCandidate");
1288  switch (rMap[name]) {
1289  case Onia:
1290  recoOnia = true;
1292  break;
1293  case Pmm:
1294  case Psi1:
1295  case Psi2:
1296  case Ups:
1297  case Ups1:
1298  case Ups2:
1299  case Ups3:
1300  recoOnia = true;
1301  break;
1302  case Kx0:
1303  recoKx0 = true;
1305  break;
1306  case Pkk:
1307  recoPkk = true;
1309  break;
1310  case Bu:
1311  recoBu = true;
1313  break;
1314  case Bd:
1315  recoBd = true;
1317  break;
1318  case Bs:
1319  recoBs = true;
1321  break;
1322  case K0s:
1323  recoK0s = true;
1325  break;
1326  case Lambda0:
1327  recoLambda0 = true;
1329  break;
1330  case B0:
1331  recoB0 = true;
1333  break;
1334  case Lambdab:
1335  recoLambdab = true;
1337  break;
1338  case Bc:
1339  recoBc = true;
1341  break;
1342  case X3872:
1343  recoX3872 = true;
1345  break;
1346  }
1347 
1348  map<string, parType>::const_iterator pIter = pMap.begin();
1349  map<string, parType>::const_iterator pIend = pMap.end();
1350  while (pIter != pIend) {
1351  const map<string, parType>::value_type& entry = *pIter++;
1352  const string& pn = entry.first;
1353  parType id = entry.second;
1354  double pv = ps.getParameter<double>(pn);
1355  if (pv > -1.0e35)
1356  edm::LogVerbatim("Configuration") << "BPHWriteSpecificDecay::setRecoParameters: set " << pn << " for " << name
1357  << " : " << (parMap[rMap[name]][id] = pv);
1358  }
1359 
1360  map<string, parType>::const_iterator fIter = fMap.begin();
1361  map<string, parType>::const_iterator fIend = fMap.end();
1362  while (fIter != fIend) {
1363  const map<string, parType>::value_type& entry = *fIter++;
1364  const string& fn = entry.first;
1365  parType id = entry.second;
1366  edm::LogVerbatim("Configuration") << "BPHWriteSpecificDecay::setRecoParameters: set " << fn << " for " << name
1367  << " : " << (parMap[rMap[name]][id] = (ps.getParameter<bool>(fn) ? 1 : -1));
1368  }
1369 }
Log< level::Info, true > LogVerbatim
uint16_t *__restrict__ id
std::map< std::string, recoType > rMap
std::map< std::string, parType > pMap
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::map< recoType, std::map< parType, double > > parMap
list entry
Definition: mps_splice.py:68
std::map< std::string, parType > fMap
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 196 of file BPHWriteSpecificDecay.h.

References pat::PATObject< ObjectType >::addUserData(), pat::PATObject< ObjectType >::addUserFloat(), ccRefMap, daughMap, 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, edm::Event::put(), pvRefMap, AlCaHLTBitMon_QueryRunRegistry::string, writeMomentum, and writeVertex.

198  {
200  int i;
201  int n = list.size();
202  std::map<const BPHRecoCandidate*, const BPHRecoCandidate*>::const_iterator dauIter;
203  std::map<const BPHRecoCandidate*, const BPHRecoCandidate*>::const_iterator dauIend = daughMap.end();
204  std::map<const BPHRecoCandidate*, const BPHRecoCandidate*>::const_iterator jpoIter;
205  std::map<const BPHRecoCandidate*, const BPHRecoCandidate*>::const_iterator jpoIend = jPsiOMap.end();
206  std::map<const BPHRecoCandidate*, vertex_ref>::const_iterator pvrIter;
207  std::map<const BPHRecoCandidate*, vertex_ref>::const_iterator pvrIend = pvRefMap.end();
208  std::map<const BPHRecoCandidate*, compcc_ref>::const_iterator ccrIter;
209  std::map<const BPHRecoCandidate*, compcc_ref>::const_iterator ccrIend = ccRefMap.end();
210  for (i = 0; i < n; ++i) {
211  const T& ptr = list[i];
212  ccList->push_back(ptr->composite());
213  pat::CompositeCandidate& cc = ccList->back();
214  if ((pvrIter = pvRefMap.find(ptr.get())) != pvrIend)
215  cc.addUserData("primaryVertex", pvrIter->second);
216  const std::vector<std::string>& cNames = ptr->compNames();
217  int j = 0;
218  int m = cNames.size();
219  while (j < m) {
220  const std::string& compName = cNames[j++];
221  const BPHRecoCandidate* cptr = ptr->getComp(compName).get();
222  if ((ccrIter = ccRefMap.find(cptr)) == ccrIend) {
223  if ((dauIter = daughMap.find(cptr)) != dauIend)
224  cptr = dauIter->second;
225  if ((jpoIter = jPsiOMap.find(cptr)) != jpoIend)
226  cptr = jpoIter->second;
227  }
228  if ((ccrIter = ccRefMap.find(cptr)) != ccrIend) {
229  compcc_ref cref = ccrIter->second;
230  if (cref.isNonnull())
231  cc.addUserData("refTo" + compName, cref);
232  }
233  }
234  const BPHPlusMinusCandidate* pmp = dynamic_cast<const BPHPlusMinusCandidate*>(ptr.get());
235  if (pmp != nullptr) {
236  cc.addUserData("cowboy", pmp->isCowboy());
237  // cc.addUserFloat( "dca", pmp->cAppInRPhi().distance() );
238  }
239  if (writeVertex)
240  cc.addUserData("vertex", ptr->vertex());
241  if (ptr->isEmpty())
242  continue;
243  if (writeVertex)
244  cc.addUserData("fitVertex", reco::Vertex(*ptr->topDecayVertex()));
245  if (ptr->isValidFit()) {
246  const RefCountedKinematicParticle kinPart = ptr->topParticle();
247  const KinematicState kinStat = kinPart->currentState();
248  cc.addUserFloat("fitMass", kinStat.mass());
249  if (writeMomentum)
250  cc.addUserData("fitMomentum", kinStat.kinematicParameters().momentum());
251  }
252  }
253  typedef std::unique_ptr<pat::CompositeCandidateCollection> ccc_pointer;
254  edm::OrphanHandle<pat::CompositeCandidateCollection> ccHandle = ev.put(ccc_pointer(ccList), name);
255  for (i = 0; i < n; ++i) {
256  const BPHRecoCandidate* ptr = list[i].get();
258  ccRefMap[ptr] = ccRef;
259  }
260  return ccHandle;
261  }
Analysis-level particle class.
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
std::map< const BPHRecoCandidate *, compcc_ref > ccRefMap
edm::Ref< pat::CompositeCandidateCollection > compcc_ref
void addUserFloat(const std::string &label, float data, const bool overwrite=false)
Set user-defined float.
Definition: PATObject.h:892
std::map< const BPHRecoCandidate *, const BPHRecoCandidate * > daughMap
ParticleMass mass() const
bool isCowboy() const
get cowboy/sailor classification
KinematicParameters const & kinematicParameters() const
std::vector< CompositeCandidate > CompositeCandidateCollection
std::map< const BPHRecoCandidate *, vertex_ref > pvRefMap
virtual BPHRecoConstCandPtr getComp(const std::string &name) const
std::map< const BPHRecoCandidate *, const BPHRecoCandidate * > jPsiOMap
long double T
void addUserData(const std::string &label, const T &data, bool transientOnly=false, bool overwrite=false)
Definition: PATObject.h:348
GlobalVector momentum() const

Member Data Documentation

std::string BPHWriteSpecificDecay::b0Name
private

Definition at line 89 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::bcName
private

Definition at line 91 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::bdName
private

Definition at line 85 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::bsName
private

Definition at line 86 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::buName
private

Definition at line 84 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::ccCandsLabel
private

Definition at line 49 of file BPHWriteSpecificDecay.h.

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

Definition at line 61 of file BPHWriteSpecificDecay.h.

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

Definition at line 191 of file BPHWriteSpecificDecay.h.

Referenced by write().

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

Definition at line 187 of file BPHWriteSpecificDecay.h.

Referenced by write().

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

Definition at line 140 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::gpCandsLabel
private

Definition at line 52 of file BPHWriteSpecificDecay.h.

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

Definition at line 64 of file BPHWriteSpecificDecay.h.

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

Definition at line 186 of file BPHWriteSpecificDecay.h.

Referenced by write().

std::string BPHWriteSpecificDecay::k0CandsLabel
private

Definition at line 53 of file BPHWriteSpecificDecay.h.

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

Definition at line 65 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::k0Name
private

Definition at line 87 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::kSCandsLabel
private

Definition at line 55 of file BPHWriteSpecificDecay.h.

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

Definition at line 67 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::l0CandsLabel
private

Definition at line 54 of file BPHWriteSpecificDecay.h.

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

Definition at line 66 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::l0Name
private

Definition at line 88 of file BPHWriteSpecificDecay.h.

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

Definition at line 181 of file BPHWriteSpecificDecay.h.

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

Definition at line 183 of file BPHWriteSpecificDecay.h.

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

Definition at line 177 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::lbName
private

Definition at line 90 of file BPHWriteSpecificDecay.h.

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

Definition at line 178 of file BPHWriteSpecificDecay.h.

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

Definition at line 176 of file BPHWriteSpecificDecay.h.

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

Definition at line 172 of file BPHWriteSpecificDecay.h.

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

Definition at line 173 of file BPHWriteSpecificDecay.h.

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

Definition at line 179 of file BPHWriteSpecificDecay.h.

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

Definition at line 180 of file BPHWriteSpecificDecay.h.

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

Definition at line 182 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::lSCandsLabel
private

Definition at line 56 of file BPHWriteSpecificDecay.h.

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

Definition at line 68 of file BPHWriteSpecificDecay.h.

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

Definition at line 174 of file BPHWriteSpecificDecay.h.

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

Definition at line 175 of file BPHWriteSpecificDecay.h.

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

Definition at line 184 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::oniaName
private

Definition at line 81 of file BPHWriteSpecificDecay.h.

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

Definition at line 141 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::patMuonLabel
private

Definition at line 48 of file BPHWriteSpecificDecay.h.

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

Definition at line 60 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::pcCandsLabel
private

Definition at line 51 of file BPHWriteSpecificDecay.h.

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

Definition at line 63 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::pfCandsLabel
private

Definition at line 50 of file BPHWriteSpecificDecay.h.

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

Definition at line 62 of file BPHWriteSpecificDecay.h.

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

Definition at line 139 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::pVertexLabel
private

Definition at line 47 of file BPHWriteSpecificDecay.h.

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

Definition at line 59 of file BPHWriteSpecificDecay.h.

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

Definition at line 189 of file BPHWriteSpecificDecay.h.

Referenced by write().

bool BPHWriteSpecificDecay::recoB0
private

Definition at line 151 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::recoBc
private

Definition at line 153 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::recoBd
private

Definition at line 147 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::recoBs
private

Definition at line 148 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::recoBu
private

Definition at line 146 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::recoK0s
private

Definition at line 149 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::recoKx0
private

Definition at line 144 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::recoLambda0
private

Definition at line 150 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::recoLambdab
private

Definition at line 152 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::recoOnia
private

Definition at line 143 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::recoPkk
private

Definition at line 145 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::recoX3872
private

Definition at line 154 of file BPHWriteSpecificDecay.h.

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

Definition at line 138 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::sdName
private

Definition at line 82 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::ssName
private

Definition at line 83 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::useCC
private

Definition at line 72 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::useGP
private

Definition at line 75 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::useK0
private

Definition at line 76 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::useKS
private

Definition at line 78 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::useL0
private

Definition at line 77 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::useLS
private

Definition at line 79 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::usePC
private

Definition at line 74 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::usePF
private

Definition at line 73 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::usePM
private

Definition at line 71 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::usePV
private

Definition at line 70 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::writeB0
private

Definition at line 164 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::writeBc
private

Definition at line 166 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::writeBd
private

Definition at line 160 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::writeBs
private

Definition at line 161 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::writeBu
private

Definition at line 159 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::writeK0s
private

Definition at line 162 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::writeKx0
private

Definition at line 157 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::writeLambda0
private

Definition at line 163 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::writeLambdab
private

Definition at line 165 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::writeMomentum
private

Definition at line 170 of file BPHWriteSpecificDecay.h.

Referenced by write().

bool BPHWriteSpecificDecay::writeOnia
private

Definition at line 156 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::writePkk
private

Definition at line 158 of file BPHWriteSpecificDecay.h.

bool BPHWriteSpecificDecay::writeVertex
private

Definition at line 169 of file BPHWriteSpecificDecay.h.

Referenced by write().

bool BPHWriteSpecificDecay::writeX3872
private

Definition at line 167 of file BPHWriteSpecificDecay.h.

std::string BPHWriteSpecificDecay::x3872Name
private

Definition at line 92 of file BPHWriteSpecificDecay.h.