16 #include "TLorentzVector.h"
26 LogDebug(
"Alignment") <<
"> applying two body decay Trackfilter ...";
27 theMassrangeSwitch = cfg.
getParameter<
bool>(
"applyMassrangeFilter" );
28 if (theMassrangeSwitch){
31 theDaughterMass = cfg.
getParameter<
double>(
"daughterMass" );
32 theCandNumber = cfg.
getParameter<
unsigned int>(
"numberOfCandidates" );
33 secThrBool = cfg.
getParameter<
bool> (
"applySecThreshold" );
34 thesecThr = cfg.
getParameter<
double>(
"secondThreshold" );
35 LogDebug(
"Alignment") <<
"> Massrange min,max : " << theMinMass <<
"," << theMaxMass
36 <<
"\n> Mass of daughter Particle : " << theDaughterMass;
43 theChargeSwitch = cfg.
getParameter<
bool>(
"applyChargeFilter" );
46 theUnsignedSwitch = cfg.
getParameter<
bool>(
"useUnsignedCharge" );
49 LogDebug(
"Alignment") <<
"> Desired Charge, unsigned: "<<theCharge<<
" , "<<theUnsignedSwitch;
52 theUnsignedSwitch =
true;
54 theMissingETSwitch = cfg.
getParameter<
bool>(
"applyMissingETFilter" );
55 if(theMissingETSwitch){
58 LogDebug(
"Alignment") <<
"> missing Et Source: "<< theMissingETSource;
60 theAcoplanarityFilterSwitch = cfg.
getParameter<
bool>(
"applyAcoplanarityFilter" );
61 if(theAcoplanarityFilterSwitch){
62 theAcoplanarDistance = cfg.
getParameter<
double>(
"acoplanarDistance" );
63 LogDebug(
"Alignment") <<
"> Acoplanar Distance: "<<theAcoplanarDistance;
77 return theMassrangeSwitch || theChargeSwitch || theAcoplanarityFilterSwitch;
87 if (theMassrangeSwitch) {
88 if (theMissingETSwitch)
89 result = checkMETMass(result,iEvent);
91 result = checkMass(result);
94 LogDebug(
"Alignment") <<
"> TwoBodyDecay tracks all,kept: " << tracks.size() <<
"," << result.size();
103 return a.first > b.first ;
113 LogDebug(
"Alignment") <<
"> cands size : "<< cands.size();
115 if (cands.size()<2)
return result;
117 TLorentzVector track0;
118 TLorentzVector track1;
119 TLorentzVector mother;
120 typedef pair<const reco::Track*,const reco::Track*> constTrackPair;
121 typedef pair<double,constTrackPair> candCollectionItem;
122 vector<candCollectionItem> candCollection;
124 for (
unsigned int iCand = 0; iCand < cands.size(); iCand++) {
126 track0.SetXYZT(cands.at(iCand)->px(),
127 cands.at(iCand)->py(),
128 cands.at(iCand)->pz(),
129 sqrt( cands.at(iCand)->p()*cands.at(iCand)->p() + theDaughterMass*theDaughterMass ));
131 for (
unsigned int jCand = iCand+1; jCand < cands.size(); jCand++) {
133 track1.SetXYZT(cands.at(jCand)->px(),
134 cands.at(jCand)->py(),
135 cands.at(jCand)->pz(),
136 sqrt( cands.at(jCand)->p()*cands.at(jCand)->p() + theDaughterMass*theDaughterMass ));
137 if (secThrBool==
true && track1.Pt() < thesecThr && track0.Pt()< thesecThr)
continue;
138 mother = track0 + track1;
143 bool correctCharge =
true;
144 if (theChargeSwitch) correctCharge = this->
checkCharge(trk1, trk2);
146 bool acoplanarTracks =
true;
147 if (theAcoplanarityFilterSwitch) acoplanarTracks = this->checkAcoplanarity(trk1, trk2);
149 if (mother.M() > theMinMass &&
150 mother.M() < theMaxMass &&
153 candCollection.push_back(candCollectionItem(mother.Pt(),
154 constTrackPair(trk1, trk2)));
159 if (candCollection.size()==0)
return result;
161 sort(candCollection.begin(), candCollection.end(),
164 std::map<const reco::Track*,unsigned int> uniqueTrackIndex;
165 std::map<const reco::Track*,unsigned int>::iterator it;
166 for (
unsigned int i=0;
167 i<candCollection.size() &&
i<theCandNumber;
169 constTrackPair & trackPair = candCollection[
i].second;
171 it = uniqueTrackIndex.find(trackPair.first);
172 if (it==uniqueTrackIndex.end()) {
173 result.push_back(trackPair.first);
174 uniqueTrackIndex[trackPair.first] =
i;
177 it = uniqueTrackIndex.find(trackPair.second);
178 if (it==uniqueTrackIndex.end()) {
179 result.push_back(trackPair.second);
180 uniqueTrackIndex[trackPair.second] =
i;
193 LogDebug(
"Alignment") <<
"> cands size : "<< cands.size();
195 if (cands.size()==0)
return result;
197 TLorentzVector track;
199 TLorentzVector mother;
202 iEvent.
getByToken(theMissingETToken ,missingET);
204 LogError(
"Alignment")<<
"@SUB=AlignmentTwoBodyDecayTrackSelector::checkMETMass"
205 <<
"> could not optain missingET Collection!";
209 typedef pair<double,const reco::Track*> candCollectionItem;
210 vector<candCollectionItem> candCollection;
212 for (reco::CaloMETCollection::const_iterator itMET = missingET->begin();
213 itMET != missingET->end();
216 met4.SetXYZT((*itMET).px(),
221 for (
unsigned int iCand = 0; iCand < cands.size(); iCand++) {
223 track.SetXYZT(cands.at(iCand)->px(),
224 cands.at(iCand)->py(),
225 cands.at(iCand)->pz(),
226 sqrt( cands.at(iCand)->p()*cands.at(iCand)->p() + theDaughterMass*theDaughterMass ));
228 mother = track + met4;
233 bool correctCharge =
true;
234 if (theChargeSwitch) correctCharge = this->
checkCharge(trk);
236 bool acoplanarTracks =
true;
237 if (theAcoplanarityFilterSwitch) acoplanarTracks = this->checkMETAcoplanarity(trk, met);
239 if (mother.M() > theMinMass &&
240 mother.M() < theMaxMass &&
243 candCollection.push_back(candCollectionItem(mother.Pt(), trk));
248 if (candCollection.size()==0)
return result;
250 sort(candCollection.begin(), candCollection.end(),
253 std::map<const reco::Track*,unsigned int> uniqueTrackIndex;
254 std::map<const reco::Track*,unsigned int>::iterator it;
255 for (
unsigned int i=0;
256 i<candCollection.size() &&
i<theCandNumber;
258 it = uniqueTrackIndex.find(candCollection[
i].
second);
259 if (it==uniqueTrackIndex.end()) {
260 result.push_back(candCollection[
i].second);
261 uniqueTrackIndex[candCollection[
i].second] =
i;
272 int sumCharge = trk1->
charge();
273 if (trk2) sumCharge += trk2->
charge();
274 if (theUnsignedSwitch) sumCharge =
std::abs(sumCharge);
275 if (sumCharge == theCharge)
return true;
301 LogDebug(
"Alignment") <<
">......................................";
302 for(Tracks::const_iterator it = col.begin();it < col.end();++it,++
count){
304 <<
"> Track No. "<< count <<
": p = ("<<(*it)->px()<<
","<<(*it)->py()<<
","<<(*it)->pz()<<
")\n"
305 <<
"> pT = "<<(*it)->pt()<<
" eta = "<<(*it)->eta()<<
" charge = "<<(*it)->charge();
307 LogDebug(
"Alignment") <<
">......................................";
std::vector< const reco::Track * > Tracks
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
AlignmentTwoBodyDecayTrackSelector(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
constructor
bool getByToken(EDGetToken token, Handle< PROD > &result) const
virtual double phi() const final
momentum azimuthal angle
double phi() const
azimuthal angle of momentum vector
bool checkMETAcoplanarity(const reco::Track *trk, const reco::CaloMET *met) const
checks if [cands] contains a acoplanar track w.r.t missing ET (returns empty set if not) ...
bool checkAcoplanarity(const reco::Track *trk1, const reco::Track *trk2) const
checks if the [cands] are acoplanar (returns empty set if not)
U second(std::pair< T, U > const &p)
~AlignmentTwoBodyDecayTrackSelector()
destructor
bool checkCharge(const reco::Track *trk1, const reco::Track *trk2=0) const
checks if the mother has charge = [theCharge]
bool useThisFilter()
returns if any of the Filters is used.
void printTracks(const Tracks &col) const
print Information on Track-Collection
Abs< T >::type abs(const T &t)
Tracks select(const Tracks &tracks, const edm::Event &iEvent, const edm::EventSetup &iSetup)
select tracks
bool operator()(const T &a, const T &b)
std::vector< reco::CaloMET > CaloMETCollection
collection of CaloMET objects
int charge() const
track electric charge
Tracks checkMETMass(const Tracks &cands, const edm::Event &iEvent) const
checks if the mass of the mother is in the mass region adding missing E_T
Tracks checkMass(const Tracks &cands) const
checks if the mass of the mother is in the mass region