#include <Alignment/CommonAlignmentProducer/interface/AlignmentTwoBodyDecayTrackSelector.h>
Public Types | |
typedef std::vector< const reco::Track * > | Tracks |
Public Member Functions | |
AlignmentTwoBodyDecayTrackSelector (const edm::ParameterSet &cfg) | |
constructor | |
Tracks | select (const Tracks &tracks, const edm::Event &iEvent) |
select tracks | |
bool | useThisFilter () |
returns if any of the Filters is used. | |
~AlignmentTwoBodyDecayTrackSelector () | |
destructor | |
Private Member Functions | |
Tracks | checkAcoplanarity (const Tracks &cands) const |
checks if the [cands] are acoplanar (returns empty set if not) | |
Tracks | checkCharge (const Tracks &cands) const |
checks if the mother has charge = [theCharge] | |
Tracks | checkMass (const Tracks &cands) const |
checks if the mass of the mother is in the mass region | |
Tracks | checkMETAcoplanarity (const Tracks &cands, const edm::Event &iEvent) const |
checks if [cands] contains a acoplanar track w.r.t missing ET (returns empty set if not) | |
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 | |
void | printTracks (const Tracks &col) const |
print Information on Track-Collection | |
Private Attributes | |
double | theAcoplanarDistance |
bool | theAcoplanarityFilterSwitch |
int | theCharge |
bool | theChargeSwitch |
double | theDaughterMass |
bool | theMassrangeSwitch |
private data members | |
double | theMaxMass |
double | theMinMass |
edm::InputTag | theMissingETSource |
bool | theMissingETSwitch |
bool | theUnsignedSwitch |
Definition at line 14 of file AlignmentTwoBodyDecayTrackSelector.h.
typedef std::vector<const reco::Track*> AlignmentTwoBodyDecayTrackSelector::Tracks |
Definition at line 18 of file AlignmentTwoBodyDecayTrackSelector.h.
AlignmentTwoBodyDecayTrackSelector::AlignmentTwoBodyDecayTrackSelector | ( | const edm::ParameterSet & | cfg | ) |
constructor
Definition at line 24 of file AlignmentTwoBoyDecayTrackSelector.cc.
References funct::abs(), edm::ParameterSet::getParameter(), LogDebug, theAcoplanarDistance, theAcoplanarityFilterSwitch, theCharge, theChargeSwitch, theDaughterMass, theMassrangeSwitch, theMaxMass, theMinMass, theMissingETSource, theMissingETSwitch, and theUnsignedSwitch.
00024 : 00025 theMissingETSource("met") 00026 { 00027 LogDebug("Alignment") << "> applying two body decay Trackfilter ..."; 00028 theMassrangeSwitch = cfg.getParameter<bool>( "applyMassrangeFilter" ); 00029 if (theMassrangeSwitch){ 00030 theMinMass = cfg.getParameter<double>( "minXMass" ); 00031 theMaxMass = cfg.getParameter<double>( "maxXMass" ); 00032 theDaughterMass = cfg.getParameter<double>( "daughterMass" ); 00033 LogDebug("Alignment") << "> Massrange min,max : " << theMinMass << "," << theMaxMass 00034 << "\n> Mass of daughter Particle : " << theDaughterMass; 00035 00036 }else{ 00037 theMinMass = 0; 00038 theMaxMass = 0; 00039 theDaughterMass = 0; 00040 } 00041 theChargeSwitch = cfg.getParameter<bool>( "applyChargeFilter" ); 00042 if(theChargeSwitch){ 00043 theCharge = cfg.getParameter<int>( "charge" ); 00044 theUnsignedSwitch = cfg.getParameter<bool>( "useUnsignedCharge" ); 00045 if(theUnsignedSwitch) 00046 theCharge=abs(theCharge); 00047 LogDebug("Alignment") << "> Desired Charge, unsigned: "<<theCharge<<" , "<<theUnsignedSwitch; 00048 }else{ 00049 theCharge =0; 00050 theUnsignedSwitch = true; 00051 } 00052 theMissingETSwitch = cfg.getParameter<bool>( "applyMissingETFilter" ); 00053 if(theMissingETSwitch){ 00054 theMissingETSource = cfg.getParameter<InputTag>( "missingETSource" ); 00055 LogDebug("Alignment") << "> missing Et Source: "<< theMissingETSource; 00056 } 00057 theAcoplanarityFilterSwitch = cfg.getParameter<bool>( "applyAcoplanarityFilter" ); 00058 if(theAcoplanarityFilterSwitch){ 00059 theAcoplanarDistance = cfg.getParameter<double>( "acoplanarDistance" ); 00060 LogDebug("Alignment") << "> Acoplanar Distance: "<<theAcoplanarDistance; 00061 } 00062 00063 }
AlignmentTwoBodyDecayTrackSelector::~AlignmentTwoBodyDecayTrackSelector | ( | ) |
AlignmentTwoBodyDecayTrackSelector::Tracks AlignmentTwoBodyDecayTrackSelector::checkAcoplanarity | ( | const Tracks & | cands | ) | const [private] |
checks if the [cands] are acoplanar (returns empty set if not)
Definition at line 171 of file AlignmentTwoBoyDecayTrackSelector.cc.
References deltaPhi(), lat::endl(), LogDebug, HLT_VtxMuL3::result, and theAcoplanarDistance.
Referenced by select().
00172 { 00173 Tracks result; result.clear(); 00174 //TODO return the biggest set of acoplanar tracks or two tracks with smallest distance? 00175 if(cands.size() == 2){ 00176 LogDebug("Alignment") <<"> Acoplanarity: "<<fabs(fabs(deltaPhi(cands.at(0)->phi(),cands.at(1)->phi()))-M_PI)<<endl; 00177 if(fabs(fabs(deltaPhi(cands.at(0)->phi(),cands.at(1)->phi()))-M_PI)<theAcoplanarDistance) 00178 result = cands; 00179 } 00180 return result; 00181 }
AlignmentTwoBodyDecayTrackSelector::Tracks AlignmentTwoBodyDecayTrackSelector::checkCharge | ( | const Tracks & | cands | ) | const [private] |
checks if the mother has charge = [theCharge]
Definition at line 155 of file AlignmentTwoBoyDecayTrackSelector.cc.
References funct::abs(), it, HLT_VtxMuL3::result, theCharge, and theUnsignedSwitch.
Referenced by select().
00156 { 00157 Tracks result; result.clear(); 00158 int sumCharge = 0; 00159 for(Tracks::const_iterator it = cands.begin();it < cands.end();++it) 00160 sumCharge += (*it)->charge(); 00161 if(theUnsignedSwitch) 00162 sumCharge = abs(sumCharge); 00163 if(sumCharge == theCharge) 00164 result = cands; 00165 00166 return result; 00167 }
AlignmentTwoBodyDecayTrackSelector::Tracks AlignmentTwoBodyDecayTrackSelector::checkMass | ( | const Tracks & | cands | ) | const [private] |
checks if the mass of the mother is in the mass region
checks if the mass of the X is in the mass region
Definition at line 107 of file AlignmentTwoBoyDecayTrackSelector.cc.
References LogDebug, HLT_VtxMuL3::result, funct::sqrt(), theDaughterMass, theMaxMass, and theMinMass.
Referenced by select().
00108 { 00109 Tracks result; result.clear(); 00110 //TODO perhaps try combinations if there are more than 2 tracks .... 00111 if(cands.size() == 2){ 00112 //TODO use other vectors here 00113 TLorentzVector track0(cands.at(0)->px(),cands.at(0)->py(),cands.at(0)->pz(), 00114 sqrt((cands.at(0)->p()*cands.at(0)->p())+theDaughterMass*theDaughterMass)); 00115 TLorentzVector track1(cands.at(1)->px(),cands.at(1)->py(),cands.at(1)->pz(), 00116 sqrt((cands.at(1)->p()*cands.at(1)->p())+theDaughterMass*theDaughterMass)); 00117 TLorentzVector mother = track0+track1; 00118 if(mother.M() > theMinMass && mother.M() < theMaxMass) 00119 result = cands; 00120 LogDebug("Alignment") <<"> mass of mother: "<<mother.M()<<"GeV"; 00121 } 00122 return result; 00123 }
AlignmentTwoBodyDecayTrackSelector::Tracks AlignmentTwoBodyDecayTrackSelector::checkMETAcoplanarity | ( | const Tracks & | cands, | |
const edm::Event & | iEvent | |||
) | const [private] |
checks if [cands] contains a acoplanar track w.r.t missing ET (returns empty set if not)
Definition at line 184 of file AlignmentTwoBoyDecayTrackSelector.cc.
References deltaPhi(), lat::endl(), edm::Event::getByLabel(), edm::Handle< T >::isValid(), LogDebug, HLT_VtxMuL3::result, theAcoplanarDistance, and theMissingETSource.
Referenced by select().
00185 { 00186 Tracks result; result.clear(); 00187 if(cands.size() == 1){ 00188 Handle<reco::CaloMETCollection> missingET; 00189 iEvent.getByLabel(theMissingETSource ,missingET); 00190 if(missingET.isValid()){ 00191 //TODO return the biggest set of acoplanar tracks or the one with smallest distance? 00192 LogDebug("Alignment") <<"> METAcoplanarity: "<<fabs(fabs(deltaPhi(cands.at(0)->phi(),(*missingET).at(0).phi()))-M_PI)<<endl; 00193 if(fabs(fabs(deltaPhi(cands.at(0)->phi(),(*missingET).at(0).phi()))-M_PI)<theAcoplanarDistance) 00194 result = cands; 00195 00196 }else 00197 LogError("Alignment")<<"@SUB=AlignmentTwoBodyDecayTrackSelector::checkMETAcoplanarity" 00198 <<"> could not optain missingET Collection!"; 00199 } 00200 return result; 00201 }
AlignmentTwoBodyDecayTrackSelector::Tracks AlignmentTwoBodyDecayTrackSelector::checkMETMass | ( | const Tracks & | cands, | |
const edm::Event & | iEvent | |||
) | const [private] |
checks if the mass of the mother is in the mass region adding missing E_T
checks if the mass of the X is in the mass region adding missing E_T
Definition at line 127 of file AlignmentTwoBoyDecayTrackSelector.cc.
References edm::Event::getByLabel(), edm::Handle< T >::isValid(), LogDebug, CaloMET_cfi::met, HLT_VtxMuL3::result, funct::sqrt(), theDaughterMass, theMaxMass, theMinMass, theMissingETSource, and track.
Referenced by select().
00128 { 00129 Tracks result; result.clear(); 00130 if(cands.size() == 1){ 00131 Handle<reco::CaloMETCollection> missingET; 00132 iEvent.getByLabel(theMissingETSource ,missingET); 00133 if(missingET.isValid()){ 00134 //TODO use the one with highest pt instead of the first one? 00135 // for(reco::CaloMETCollection::const_iterator itMET = missingET->begin(); itMET != missingET->end() ; ++itMET){ 00136 // cout <<"missingET p = ("<<(*itMET).px()<<","<<(*itMET).py()<<","<<(*itMET).pz()<<")"<<endl; 00137 //} 00138 TLorentzVector track(cands.at(0)->px(),cands.at(0)->py(),cands.at(0)->pz(), 00139 sqrt((cands.at(0)->p()*cands.at(0)->p())+theDaughterMass*theDaughterMass)); 00140 TLorentzVector met((*missingET).at(0).px(),(*missingET).at(0).py(),(*missingET).at(0).pz(), 00141 (*missingET).at(0).p());//ignoring nuetralino masses for now ;) 00142 TLorentzVector motherSystem = track + met; 00143 if(motherSystem.M() > theMinMass && motherSystem.M() < theMaxMass) 00144 result = cands; 00145 LogDebug("Alignment") <<"> mass of motherSystem: "<<motherSystem.M()<<"GeV"; 00146 }else 00147 LogError("Alignment")<<"@SUB=AlignmentTwoBodyDecayTrackSelector::checkMETMass" 00148 <<"> could not optain missingET Collection!"; 00149 } 00150 return cands; 00151 }
print Information on Track-Collection
Definition at line 206 of file AlignmentTwoBoyDecayTrackSelector.cc.
References count, it, and LogDebug.
00207 { 00208 int count = 0; 00209 LogDebug("Alignment") << ">......................................"; 00210 for(Tracks::const_iterator it = col.begin();it < col.end();++it,++count){ 00211 LogDebug("Alignment") 00212 <<"> Track No. "<< count <<": p = ("<<(*it)->px()<<","<<(*it)->py()<<","<<(*it)->pz()<<")\n" 00213 <<"> pT = "<<(*it)->pt()<<" eta = "<<(*it)->eta()<<" charge = "<<(*it)->charge(); 00214 } 00215 LogDebug("Alignment") << ">......................................"; 00216 }
AlignmentTwoBodyDecayTrackSelector::Tracks AlignmentTwoBodyDecayTrackSelector::select | ( | const Tracks & | tracks, | |
const edm::Event & | iEvent | |||
) |
select tracks
Definition at line 80 of file AlignmentTwoBoyDecayTrackSelector.cc.
References checkAcoplanarity(), checkCharge(), checkMass(), checkMETAcoplanarity(), checkMETMass(), LogDebug, HLT_VtxMuL3::result, theAcoplanarityFilterSwitch, theChargeSwitch, theMassrangeSwitch, and theMissingETSwitch.
Referenced by TrackConfigSelector::select().
00081 { 00082 Tracks result=tracks; 00083 00084 if(theMassrangeSwitch){ 00085 if(theMissingETSwitch) 00086 result = checkMETMass(result,iEvent); 00087 else 00088 result = checkMass(result); 00089 } 00090 if(theChargeSwitch) 00091 result = checkCharge(result); 00092 if(theAcoplanarityFilterSwitch){ 00093 if(theMissingETSwitch) 00094 result = checkMETAcoplanarity(result,iEvent); 00095 else 00096 result = checkAcoplanarity(result); 00097 } 00098 LogDebug("Alignment") << "> TwoBodyDecay tracks all,kept: " << tracks.size() << "," << result.size(); 00099 // LogDebug("AlignmentTwoBodyDecayTrackSelector")<<"> o kept:"; 00100 //printTracks(result); 00101 return result; 00102 00103 }
bool AlignmentTwoBodyDecayTrackSelector::useThisFilter | ( | ) |
returns if any of the Filters is used.
Definition at line 72 of file AlignmentTwoBoyDecayTrackSelector.cc.
References theAcoplanarityFilterSwitch, theChargeSwitch, and theMassrangeSwitch.
Referenced by TrackConfigSelector::TrackConfigSelector().
00073 { 00074 return theMassrangeSwitch || theChargeSwitch || theAcoplanarityFilterSwitch; 00075 }
double AlignmentTwoBodyDecayTrackSelector::theAcoplanarDistance [private] |
Definition at line 58 of file AlignmentTwoBodyDecayTrackSelector.h.
Referenced by AlignmentTwoBodyDecayTrackSelector(), checkAcoplanarity(), and checkMETAcoplanarity().
Definition at line 47 of file AlignmentTwoBodyDecayTrackSelector.h.
Referenced by AlignmentTwoBodyDecayTrackSelector(), select(), and useThisFilter().
Definition at line 53 of file AlignmentTwoBodyDecayTrackSelector.h.
Referenced by AlignmentTwoBodyDecayTrackSelector(), and checkCharge().
Definition at line 45 of file AlignmentTwoBodyDecayTrackSelector.h.
Referenced by AlignmentTwoBodyDecayTrackSelector(), select(), and useThisFilter().
double AlignmentTwoBodyDecayTrackSelector::theDaughterMass [private] |
Definition at line 51 of file AlignmentTwoBodyDecayTrackSelector.h.
Referenced by AlignmentTwoBodyDecayTrackSelector(), checkMass(), and checkMETMass().
private data members
Definition at line 44 of file AlignmentTwoBodyDecayTrackSelector.h.
Referenced by AlignmentTwoBodyDecayTrackSelector(), select(), and useThisFilter().
double AlignmentTwoBodyDecayTrackSelector::theMaxMass [private] |
Definition at line 50 of file AlignmentTwoBodyDecayTrackSelector.h.
Referenced by AlignmentTwoBodyDecayTrackSelector(), checkMass(), and checkMETMass().
double AlignmentTwoBodyDecayTrackSelector::theMinMass [private] |
Definition at line 49 of file AlignmentTwoBodyDecayTrackSelector.h.
Referenced by AlignmentTwoBodyDecayTrackSelector(), checkMass(), and checkMETMass().
Definition at line 56 of file AlignmentTwoBodyDecayTrackSelector.h.
Referenced by AlignmentTwoBodyDecayTrackSelector(), checkMETAcoplanarity(), and checkMETMass().
Definition at line 46 of file AlignmentTwoBodyDecayTrackSelector.h.
Referenced by AlignmentTwoBodyDecayTrackSelector(), and select().
Definition at line 54 of file AlignmentTwoBodyDecayTrackSelector.h.
Referenced by AlignmentTwoBodyDecayTrackSelector(), and checkCharge().