Go to the documentation of this file.00001 #include "DataFormats/ParticleFlowReco/interface/PFDisplacedVertex.h"
00002
00003 #include "TMath.h"
00004
00005 using namespace std;
00006 using namespace reco;
00007
00008
00009 PFDisplacedVertex::PFDisplacedVertex() : Vertex(),
00010 vertexType_(ANY),
00011 primaryDirection_(0,0,0)
00012 {}
00013
00014 PFDisplacedVertex::PFDisplacedVertex(Vertex& v) : Vertex(v),
00015 vertexType_(ANY),
00016 primaryDirection_(0,0,0)
00017 {}
00018
00019 void
00020 PFDisplacedVertex::addElement( const TrackBaseRef & r, const Track & refTrack,
00021 const PFTrackHitFullInfo& hitInfo ,
00022 VertexTrackType trackType, float w ) {
00023 add(r, refTrack, w );
00024 trackTypes_.push_back(trackType);
00025 trackHitFullInfos_.push_back(hitInfo);
00026 }
00027
00028 void
00029 PFDisplacedVertex::cleanTracks() {
00030
00031 removeTracks();
00032 trackTypes_.clear();
00033 trackHitFullInfos_.clear();
00034
00035 }
00036
00037
00038 const bool
00039 PFDisplacedVertex::isThereKindTracks(VertexTrackType T) const {
00040
00041 vector <VertexTrackType>::const_iterator iter =
00042 find (trackTypes_.begin(), trackTypes_.end(), T);
00043 return (iter != trackTypes_.end()) ;
00044
00045 }
00046
00047 const int
00048 PFDisplacedVertex::nKindTracks(VertexTrackType T) const {
00049
00050 return count ( trackTypes_.begin(), trackTypes_.end(), T);
00051
00052 }
00053
00054
00055 const size_t
00056 PFDisplacedVertex::trackPosition(const reco::TrackBaseRef& originalTrack) const {
00057
00058 size_t pos = -1;
00059
00060 const Track refittedTrack = PFDisplacedVertex::refittedTrack(originalTrack);
00061
00062 std::vector<Track> refitTrks = refittedTracks();
00063 for (size_t i = 0; i < refitTrks.size(); i++){
00064 if ( fabs(refitTrks[i].pt() - refittedTrack.pt()) < 1.e-5 ){
00065 pos = i;
00066 continue;
00067 }
00068
00069 }
00070
00071
00072 return pos;
00073
00074 }
00075
00076
00077 void
00078 PFDisplacedVertex::setPrimaryDirection(const math::XYZPoint& pvtx){
00079 primaryDirection_ = math::XYZVector(position().x(), position().y(), position().z());
00080 math::XYZVector vtx(pvtx.x(), pvtx.y(), pvtx.z());
00081
00082 primaryDirection_ = primaryDirection_ - vtx;
00083 primaryDirection_ /= (sqrt(primaryDirection_.Mag2())+1e-10);
00084 }
00085
00086
00087 std::string
00088 PFDisplacedVertex::nameVertexType() const {
00089 switch (vertexType_){
00090 case ANY: return "ANY";
00091 case FAKE: return "FAKE";
00092 case LOOPER: return "LOOPER";
00093 case NUCL: return "NUCL";
00094 case NUCL_LOOSE: return "NUCL_LOOSE";
00095 case NUCL_KINK: return "NUCL_KINK";
00096 case CONVERSION: return "CONVERSION";
00097 case CONVERSION_LOOSE: return "CONVERSION_LOOSE";
00098 case CONVERTED_BREMM: return "CONVERTED_BREMM";
00099 case K0_DECAY: return "K0_DECAY";
00100 case LAMBDA_DECAY: return "LAMBDA_DECAY";
00101 case LAMBDABAR_DECAY: return "LAMBDABAR_DECAY";
00102 case KPLUS_DECAY: return "KPLUS_DECAY";
00103 case KMINUS_DECAY: return "KMINUS_DECAY";
00104 case KPLUS_DECAY_LOOSE: return "KPLUS_DECAY_LOOSE";
00105 case KMINUS_DECAY_LOOSE: return "KMINUS_DECAY_LOOSE";
00106 case BSM_VERTEX: return "BSM_VERTEX";
00107 default: return "?";
00108 }
00109 return "?";
00110 }
00111
00112
00113 const math::XYZTLorentzVector
00114 PFDisplacedVertex::momentum(string massHypo, VertexTrackType T, bool useRefitted, double mass) const {
00115
00116 M_Hypo mHypo = M_CUSTOM;
00117
00118 if (massHypo.find("PI")!=string::npos) mHypo = M_PION;
00119 else if (massHypo.find("KAON")!=string::npos) mHypo = M_KAON;
00120 else if (massHypo.find("LAMBDA")!=string::npos) mHypo = M_LAMBDA;
00121 else if (massHypo.find("MASSLESS")!=string::npos) mHypo = M_MASSLESS;
00122 else if (massHypo.find("CUSTOM")!=string::npos) mHypo = M_CUSTOM;
00123
00124 return momentum(mHypo, T, useRefitted, mass);
00125
00126 }
00127
00128
00129 const math::XYZTLorentzVector
00130 PFDisplacedVertex::momentum(M_Hypo massHypo, VertexTrackType T, bool useRefitted, double mass) const {
00131
00132 const double m2 = getMass2(massHypo, mass);
00133
00134
00135
00136 math::XYZTLorentzVector P;
00137
00138 for (size_t i = 0; i< tracksSize(); i++){
00139 bool bType = (trackTypes_[i]== T);
00140 if (T == T_TO_VERTEX || T == T_MERGED)
00141 bType = (trackTypes_[i] == T_TO_VERTEX || trackTypes_[i] == T_MERGED);
00142
00143 if ( bType ) {
00144
00145 if (!useRefitted) {
00146
00147 TrackBaseRef trackRef = originalTrack(refittedTracks()[i]);
00148
00149 double p2 = trackRef->innerMomentum().Mag2();
00150 P += math::XYZTLorentzVector (trackRef->momentum().x(),
00151 trackRef->momentum().y(),
00152 trackRef->momentum().z(),
00153 sqrt(m2 + p2));
00154 } else {
00155
00156
00157
00158 double p2 = refittedTracks()[i].momentum().Mag2();
00159 P += math::XYZTLorentzVector (refittedTracks()[i].momentum().x(),
00160 refittedTracks()[i].momentum().y(),
00161 refittedTracks()[i].momentum().z(),
00162 sqrt(m2 + p2));
00163
00164
00165 }
00166 }
00167 }
00168
00169 return P;
00170
00171 }
00172
00173
00174 const int
00175 PFDisplacedVertex::totalCharge() const {
00176
00177 int charge = 0;
00178
00179 for (size_t i = 0; i< tracksSize(); i++){
00180 if(trackTypes_[i] == T_TO_VERTEX) charge += refittedTracks()[i].charge();
00181 else if(trackTypes_[i] == T_FROM_VERTEX) charge -= refittedTracks()[i].charge();
00182 }
00183
00184 return charge;
00185 }
00186
00187
00188 const double
00189 PFDisplacedVertex::angle_io() const {
00190 math::XYZTLorentzVector momentumSec = secondaryMomentum((string) "PI", true);
00191
00192 math::XYZVector p_out = momentumSec.Vect();
00193
00194 math::XYZVector p_in = primaryDirection();
00195
00196 if (p_in.Mag2() < 1e-10) return -1;
00197 return acos(p_in.Dot(p_out)/sqrt(p_in.Mag2()*p_out.Mag2()))/TMath::Pi()*180.0;
00198
00199 }
00200
00201 const math::XYZVector
00202 PFDisplacedVertex:: primaryDirection() const {
00203
00204 math::XYZTLorentzVector momentumPrim = primaryMomentum((string) "PI", true);
00205 math::XYZTLorentzVector momentumSec = secondaryMomentum((string) "PI", true);
00206
00207 math::XYZVector p_out = momentumSec.Vect();
00208
00209 math::XYZVector p_in;
00210
00211 if (( isThereKindTracks(T_TO_VERTEX) || isThereKindTracks(T_MERGED) ) &&
00212 momentumPrim.E() > momentumSec.E()){
00213 p_in = momentumPrim.Vect()/sqrt(momentumPrim.Vect().Mag2()+1e-10);
00214 } else {
00215 p_in = primaryDirection_;
00216 }
00217
00218 return p_in;
00219 }
00220
00221
00222 const double
00223 PFDisplacedVertex::getMass2(M_Hypo massHypo, double mass) const {
00224
00225
00226 double pion_mass2 = 0.0194;
00227
00228 double kaon_mass2 = 0.2476;
00229
00230 double lambda_mass2 = 1.267;
00231
00232 if (massHypo == M_PION) return pion_mass2;
00233 else if (massHypo == M_KAON) return kaon_mass2;
00234 else if (massHypo == M_LAMBDA) return lambda_mass2;
00235 else if (massHypo == M_MASSLESS) return 0;
00236 else if (massHypo == M_CUSTOM) return mass*mass;
00237
00238 cout << "Warning: undefined mass hypothesis" << endl;
00239 return 0;
00240
00241 }
00242
00243 void PFDisplacedVertex::Dump( ostream& out ) const {
00244 if(! out ) return;
00245
00246 out << "" << endl;
00247 out << "==================== This is a Displaced Vertex type " <<
00248 nameVertexType() << " ===============" << endl;
00249
00250 out << " Vertex chi2 = " << chi2() << " ndf = " << ndof()<< " normalised chi2 = " << normalizedChi2()<< endl;
00251
00252 out << " The vertex Fitted Position is: x = " << position().x()
00253 << " y = " << position().y()
00254 << " rho = " << position().rho()
00255 << " z = " << position().z()
00256 << endl;
00257
00258 out<< "\t--- Structure --- " << endl;
00259 out<< "Number of tracks: " << nTracks()
00260 << " nPrimary " << nPrimaryTracks()
00261 << " nMerged " << nMergedTracks()
00262 << " nSecondary " << nSecondaryTracks() << endl;
00263
00264 vector <PFDisplacedVertex::PFTrackHitFullInfo> pattern = trackHitFullInfos();
00265 vector <PFDisplacedVertex::VertexTrackType> trackType = trackTypes();
00266 for (unsigned i = 0; i < pattern.size(); i++){
00267 out << "track " << i
00268 << " type = " << trackType[i]
00269 << " nHit BeforeVtx = " << pattern[i].first.first
00270 << " AfterVtx = " << pattern[i].second.first
00271 << " MissHit BeforeVtx = " << pattern[i].first.second
00272 << " AfterVtx = " << pattern[i].second.second
00273 << endl;
00274 }
00275
00276 math::XYZTLorentzVector mom_prim = primaryMomentum((string) "PI", true);
00277 math::XYZTLorentzVector mom_sec = secondaryMomentum((string) "PI", true);
00278
00279
00280 out << "Primary P:\t E " << mom_prim.E()
00281 << "\tPt = " << mom_prim.Pt()
00282 << "\tPz = " << mom_prim.Pz()
00283 << "\tM = " << mom_prim.M()
00284 << "\tEta = " << mom_prim.Eta()
00285 << "\tPhi = " << mom_prim.Phi() << endl;
00286
00287 out << "Secondary P:\t E " << mom_sec.E()
00288 << "\tPt = " << mom_sec.Pt()
00289 << "\tPz = " << mom_sec.Pz()
00290 << "\tM = " << mom_sec.M()
00291 << "\tEta = " << mom_sec.Eta()
00292 << "\tPhi = " << mom_sec.Phi() << endl;
00293
00294 out << " The vertex Direction is x = " << primaryDirection().x()
00295 << " y = " << primaryDirection().y()
00296 << " z = " << primaryDirection().z()
00297 << " eta = " << primaryDirection().eta()
00298 << " phi = " << primaryDirection().phi() << endl;
00299
00300 out << " Angle_io = " << angle_io() << " deg" << endl << endl;
00301
00302 }
00303