CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/DataFormats/ParticleFlowReco/src/PFDisplacedVertexSeed.cc

Go to the documentation of this file.
00001 #include "DataFormats/ParticleFlowReco/interface/PFDisplacedVertexSeed.h"
00002 
00003 #include "DataFormats/TrackReco/interface/Track.h"
00004 #include "DataFormats/Math/interface/Point3D.h"
00005 
00006 using namespace std;
00007 using namespace reco;
00008 
00009 
00010 PFDisplacedVertexSeed::PFDisplacedVertexSeed() :
00011   seedPoint_(GlobalPoint(0,0,0)),
00012   totalWeight_(0)
00013 {}
00014 
00015 
00016 void PFDisplacedVertexSeed::addElement(TrackBaseRef element) {
00017   elements_.insert( element ); 
00018 }
00019 
00020     
00021 void PFDisplacedVertexSeed::updateSeedPoint(const GlobalPoint& dcaPoint, TrackBaseRef r1, TrackBaseRef r2, double weight){
00022     
00023   
00024   if ( isEmpty() ) {
00025     seedPoint_ = dcaPoint;
00026     totalWeight_ = weight;
00027   }
00028   else {
00029     Basic3DVector<double>vertexSeedVector(seedPoint_);
00030     Basic3DVector<double>dcaVector(dcaPoint);
00031 
00032  
00033     dcaVector = (dcaVector*weight + vertexSeedVector*totalWeight_)/(totalWeight_+weight);
00034     GlobalPoint P(dcaVector.x(), dcaVector.y(), dcaVector.z());
00035     totalWeight_ += weight;
00036     seedPoint_ = P;
00037 
00038   }
00039 
00040   addElement(r1); 
00041   addElement(r2);
00042 
00043 }
00044 
00045 
00046 void PFDisplacedVertexSeed::mergeWith(const PFDisplacedVertexSeed& displacedVertex){
00047 
00048   
00049   double weight = displacedVertex.totalWeight();
00050   set<TrackBaseRef, Compare> newElements= displacedVertex.elements();
00051   GlobalPoint dcaPoint = displacedVertex.seedPoint();
00052 
00053   Basic3DVector<double>vertexSeedVector(seedPoint_);
00054   Basic3DVector<double>dcaVector(dcaPoint);
00055 
00056   dcaVector = (dcaVector*weight + vertexSeedVector*totalWeight_)/(totalWeight_+weight);
00057   GlobalPoint P(dcaVector.x(), dcaVector.y(), dcaVector.z());
00058   totalWeight_ += weight;
00059   seedPoint_ = P;
00060 
00061 
00062 
00063   for (  set<TrackBaseRef, Compare>::const_iterator il = newElements.begin(); il != newElements.end(); il++)
00064     addElement(*il);
00065 
00066 
00067 
00068 }
00069 
00070 
00071 void PFDisplacedVertexSeed::Dump( ostream& out ) const {
00072   if(! out ) return;
00073 
00074   out<<"\t--- DisplacedVertexSeed ---  "<<endl;
00075   out<<"\tnumber of elements: "<<elements_.size()<<endl;
00076   
00077   out<<"\t Seed Point x = " << seedPoint().x() 
00078      <<"\t Seed Point y = " << seedPoint().y()
00079      <<"\t Seed Point z = " << seedPoint().z() << endl;
00080 
00081   // Build element label (string) : elid from type, layer and occurence number
00082   // use stringstream instead of sprintf to concatenate string and integer into string
00083   for(IEset ie = elements_.begin(); ie !=  elements_.end(); ie++){
00084 
00085     math::XYZPoint Pi((*ie).get()->innerPosition());
00086     math::XYZPoint Po((*ie).get()->outerPosition());
00087 
00088     float innermost_radius = sqrt(Pi.x()*Pi.x() + Pi.y()*Pi.y() + Pi.z()*Pi.z());
00089     float outermost_radius = sqrt(Po.x()*Po.x() + Po.y()*Po.y() + Po.z()*Po.z());
00090     float innermost_rho = sqrt(Pi.x()*Pi.x() + Pi.y()*Pi.y());
00091     float outermost_rho = sqrt(Po.x()*Po.x() + Po.y()*Po.y());
00092     
00093     double pt = (*ie)->pt();
00094 
00095 
00096     out<<"ie = " << (*ie).key() << " pt = " << pt
00097        <<" innermost hit radius = " << innermost_radius << " rho = " << innermost_rho
00098        <<" outermost hit radius = " << outermost_radius << " rho = " << outermost_rho
00099        <<endl;
00100 
00101     out<<"ie = " << (*ie).key() << " pt = " << pt
00102       //       <<" inn hit pos x = " << Pi.x() << " y = " << Pi.y() << " z = " << Pi.z() 
00103        <<" out hit pos x = " << Po.x() << " y = " << Po.y() << " z = " << Po.z() 
00104        <<endl;
00105 
00106   }
00107    
00108   out<<endl;
00109 
00110 
00111 }
00112   
00113