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();
104 LogDebug(
"Alignment") <<
"> cands size : "<< cands.size();
106 if (cands.size()<2)
return result;
108 TLorentzVector track0;
109 TLorentzVector track1;
110 TLorentzVector mother;
111 typedef pair<const reco::Track*,const reco::Track*> constTrackPair;
112 typedef pair<double,constTrackPair> candCollectionItem;
113 vector<candCollectionItem> candCollection;
115 for (
unsigned int iCand = 0; iCand < cands.size(); iCand++) {
117 track0.SetXYZT(cands.at(iCand)->px(),
118 cands.at(iCand)->py(),
119 cands.at(iCand)->pz(),
120 sqrt( cands.at(iCand)->p()*cands.at(iCand)->p() + theDaughterMass*theDaughterMass ));
122 for (
unsigned int jCand = iCand+1; jCand < cands.size(); jCand++) {
124 track1.SetXYZT(cands.at(jCand)->px(),
125 cands.at(jCand)->py(),
126 cands.at(jCand)->pz(),
127 sqrt( cands.at(jCand)->p()*cands.at(jCand)->p() + theDaughterMass*theDaughterMass ));
128 if (secThrBool==
true && track1.Pt() < thesecThr && track0.Pt()< thesecThr)
continue;
129 mother = track0 + track1;
134 bool correctCharge =
true;
135 if (theChargeSwitch) correctCharge = this->
checkCharge(trk1, trk2);
137 bool acoplanarTracks =
true;
138 if (theAcoplanarityFilterSwitch) acoplanarTracks = this->checkAcoplanarity(trk1, trk2);
140 if (mother.M() > theMinMass &&
141 mother.M() < theMaxMass &&
144 candCollection.push_back(candCollectionItem(mother.Pt(),
145 constTrackPair(trk1, trk2)));
150 if (candCollection.empty())
return result;
152 sort(candCollection.begin(), candCollection.end(), [](
auto&
a,
auto&
b){
return a.first >
b.first ;});
154 std::map<const reco::Track*,unsigned int> uniqueTrackIndex;
155 std::map<const reco::Track*,unsigned int>::iterator it;
156 for (
unsigned int i=0;
157 i<candCollection.size() &&
i<theCandNumber;
159 constTrackPair & trackPair = candCollection[
i].second;
161 it = uniqueTrackIndex.find(trackPair.first);
162 if (it==uniqueTrackIndex.end()) {
163 result.push_back(trackPair.first);
164 uniqueTrackIndex[trackPair.first] =
i;
167 it = uniqueTrackIndex.find(trackPair.second);
168 if (it==uniqueTrackIndex.end()) {
169 result.push_back(trackPair.second);
170 uniqueTrackIndex[trackPair.second] =
i;
183 LogDebug(
"Alignment") <<
"> cands size : "<< cands.size();
185 if (cands.empty())
return result;
187 TLorentzVector
track;
189 TLorentzVector mother;
192 iEvent.
getByToken(theMissingETToken ,missingET);
194 LogError(
"Alignment")<<
"@SUB=AlignmentTwoBodyDecayTrackSelector::checkMETMass" 195 <<
"> could not optain missingET Collection!";
199 typedef pair<double,const reco::Track*> candCollectionItem;
200 vector<candCollectionItem> candCollection;
202 for (reco::CaloMETCollection::const_iterator itMET = missingET->begin();
203 itMET != missingET->end();
206 met4.SetXYZT((*itMET).px(),
211 for (
unsigned int iCand = 0; iCand < cands.size(); iCand++) {
213 track.SetXYZT(cands.at(iCand)->px(),
214 cands.at(iCand)->py(),
215 cands.at(iCand)->pz(),
216 sqrt( cands.at(iCand)->p()*cands.at(iCand)->p() + theDaughterMass*theDaughterMass ));
218 mother = track + met4;
223 bool correctCharge =
true;
224 if (theChargeSwitch) correctCharge = this->
checkCharge(trk);
226 bool acoplanarTracks =
true;
227 if (theAcoplanarityFilterSwitch) acoplanarTracks = this->checkMETAcoplanarity(trk, met);
229 if (mother.M() > theMinMass &&
230 mother.M() < theMaxMass &&
233 candCollection.push_back(candCollectionItem(mother.Pt(), trk));
238 if (candCollection.empty())
return result;
240 sort(candCollection.begin(), candCollection.end(), [](
auto&
a,
auto&
b){
return a.first >
b.first ;});
242 std::map<const reco::Track*,unsigned int> uniqueTrackIndex;
243 std::map<const reco::Track*,unsigned int>::iterator it;
244 for (
unsigned int i=0;
245 i<candCollection.size() &&
i<theCandNumber;
247 it = uniqueTrackIndex.find(candCollection[
i].
second);
248 if (it==uniqueTrackIndex.end()) {
249 result.push_back(candCollection[
i].second);
250 uniqueTrackIndex[candCollection[
i].second] =
i;
261 int sumCharge = trk1->
charge();
262 if (trk2) sumCharge += trk2->
charge();
263 if (theUnsignedSwitch) sumCharge =
std::abs(sumCharge);
290 LogDebug(
"Alignment") <<
">......................................";
291 for(Tracks::const_iterator it = col.begin();it < col.end();++it,++
count){
293 <<
"> Track No. "<< count <<
": p = ("<<(*it)->px()<<
","<<(*it)->py()<<
","<<(*it)->pz()<<
")\n" 294 <<
"> pT = "<<(*it)->pt()<<
" eta = "<<(*it)->eta()<<
" charge = "<<(*it)->charge();
296 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
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
std::vector< reco::CaloMET > CaloMETCollection
collection of CaloMET objects
int charge() const
track electric charge
double phi() const final
momentum azimuthal angle
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