50 template<
typename T1,
typename T2>
85 template<
typename T1,
typename T2>
87 : srcCandsToken_ (consumes<edm::
View<T1> >(iConfig.getParameter<edm::InputTag> (
"srcObject")))
88 , srcObjectsTokens_ (edm::
vector_transform(iConfig.getParameter<std::vector<edm::InputTag> >(
"srcObjectsToMatch"), [this](edm::InputTag
const &
tag){
return consumes<edm::View<T2> >(
tag);}))
89 ,
deltaRMax_ (iConfig.getParameter<
double> (
"deltaRMax"))
96 produces<std::vector<T1> >();
101 template<
typename T1,
typename T2>
111 template<
typename T1,
typename T2>
114 std::auto_ptr<std::vector<T1> > cleanObjects(
new std::vector<T1 >);
119 bool* isMatch =
new bool[candidates->size()];
120 for (
unsigned int iObject=0;iObject<candidates->size();iObject++) isMatch[iObject] =
false;
122 for (
unsigned int iSrc=0;iSrc<srcObjectsTokens_.size();iSrc++) {
124 iEvent.
getByToken(srcObjectsTokens_[iSrc],objects);
126 if(objects->size()==0)
continue;
128 for (
unsigned int iObject=0;iObject<candidates->size();iObject++) {
129 const T1& candidate = candidates->at(iObject);
130 if (!
objCut_(candidate))
continue;
133 for (
unsigned int iObj=0;iObj<objects->size();iObj++) {
134 const T2&
obj = objects->at(iObj);
137 if (deltaR<
deltaRMax_) isMatch[iObject] =
true;
146 for (tIt = candidates->begin(); tIt != endcands; ++tIt, ++
counter) {
147 if(isMatch[counter]) cleanObjects->push_back( *tIt );
154 iEvent.
put(cleanObjects);
159 template<
typename T1,
typename T2>
162 std::stringstream
ss;
165 std::cout<<
"++++++++++++++++++++++++++++++++++++++++++++++++++"
166 <<
"\n"<<
moduleLabel_<<
"(ObjectViewMatcher) SUMMARY:\n"<<ss.str()
167 <<
"++++++++++++++++++++++++++++++++++++++++++++++++++"
virtual ~ObjectViewMatcher()
unsigned int nObjectsMatch_
unsigned int nObjectsTot_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
#define DEFINE_FWK_MODULE(type)
StringCutObjectSelector< T2, true > objMatchCut_
std::vector< edm::EDGetTokenT< edm::View< T2 > > > srcObjectsTokens_
objCut_(iConfig.existsAs< std::string >("srcObjectSelection")?iConfig.getParameter< std::string >("srcObjectSelection"):"", true)
objMatchCut_(iConfig.existsAs< std::string >("srcObjectsToMatchSelection")?iConfig.getParameter< std::string >("srcObjectsToMatchSelection"):"", true)
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
double deltaR(const T1 &t1, const T2 &t2)
edm::EDGetTokenT< edm::View< T1 > > srcCandsToken_
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
ObjectViewMatcher< reco::Photon, reco::Track > TrackMatchedPhotonProducer
ObjectViewMatcher< reco::Jet, reco::Track > TrackMatchedJetProducer
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
double deltaR(double eta1, double eta2, double phi1, double phi2)
ObjectViewMatcher(const edm::ParameterSet &iConfig)
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
static std::atomic< unsigned int > counter
StringCutObjectSelector< T1, true > objCut_
deltaRMax_(iConfig.getParameter< double >("deltaRMax"))