CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/RecoVZero/VZeroFinding/src/VZeroFinder.cc

Go to the documentation of this file.
00001 #include "RecoVZero/VZeroFinding/interface/VZeroFinder.h"
00002 
00003 #include "DataFormats/Math/interface/Vector3D.h"
00004 
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 
00007 #include "TrackingTools/PatternTools/interface/TwoTrackMinimumDistance.h"
00008 
00009 #include <utility>
00010 /*****************************************************************************/
00011 VZeroFinder::VZeroFinder
00012   (const edm::EventSetup& es,
00013    const edm::ParameterSet& pset)
00014 {
00015   // Get track-pair level cuts
00016   maxDca = pset.getParameter<double>("maxDca");
00017 
00018   // Get mother cuts
00019   minCrossingRadius = pset.getParameter<double>("minCrossingRadius");
00020   maxCrossingRadius = pset.getParameter<double>("maxCrossingRadius");
00021   maxImpactMother   = pset.getParameter<double>("maxImpactMother");
00022 
00023   // Get magnetic field
00024   edm::ESHandle<MagneticField> magField;
00025   es.get<IdealMagneticFieldRecord>().get(magField);
00026   theMagField = magField.product();
00027 }
00028 
00029 /*****************************************************************************/
00030 VZeroFinder::~VZeroFinder()
00031 {
00032 }
00033 
00034 /*****************************************************************************/
00035 GlobalTrajectoryParameters VZeroFinder::getGlobalTrajectoryParameters
00036   (const reco::Track& track)
00037 {
00038   GlobalPoint position(track.vertex().x(),
00039                        track.vertex().y(),
00040                        track.vertex().z());
00041 
00042   GlobalVector momentum(track.momentum().x(),
00043                         track.momentum().y(),
00044                         track.momentum().z());
00045 
00046   GlobalTrajectoryParameters gtp(position,momentum,
00047                                  track.charge(),theMagField);
00048 
00049   return gtp;
00050 }
00051 
00052 /*****************************************************************************/
00053 GlobalVector VZeroFinder::rotate(const GlobalVector & p, double a)
00054 {
00055   double pt = p.perp();
00056   return GlobalVector( -pt*cos(a), -pt*sin(a), p.z());
00057 }
00058 
00059 /*****************************************************************************/
00060 bool VZeroFinder::checkTrackPair(const reco::Track& posTrack,
00061                                  const reco::Track& negTrack,
00062                                  const reco::VertexCollection* vertices,
00063                                        reco::VZeroData& data)
00064 {
00065 /*
00066   LogTrace("MinBiasTracking") << " [VZeroFinder] tracks"
00067    << " +" << posTrack.algoName() << " " << posTrack.d0()
00068    << " -" << negTrack.algoName() << " " << negTrack.d0();
00069 */
00070 
00071   // Get trajectories
00072   GlobalTrajectoryParameters posGtp = getGlobalTrajectoryParameters(posTrack);
00073   GlobalTrajectoryParameters negGtp = getGlobalTrajectoryParameters(negTrack);
00074 
00075   // Two track minimum distance
00076   TwoTrackMinimumDistance theMinimum(TwoTrackMinimumDistance::SlowMode);
00077 
00078   // Closest approach
00079   theMinimum.calculate(posGtp,negGtp);
00080 
00081   // Closest points
00082   // std::pair<GlobalPoint, GlobalPoint> points = theMinimum.points();
00083 
00084   // Momenta at closest point
00085   std::pair<GlobalVector,GlobalVector> momenta;
00086   momenta.first  = rotate(posGtp.momentum(), theMinimum.firstAngle() );
00087   momenta.second = rotate(negGtp.momentum(), theMinimum.secondAngle());
00088 
00089   // dcaR
00090   float dca            = theMinimum.distance();
00091   GlobalPoint crossing = theMinimum.crossingPoint();
00092 
00093 /*
00094   LogTrace("MinBiasTracking") << fixed << setprecision(2)
00095     << "  [VZeroFinder] dca    = "<<dca<<" (<"<<maxDca<<")";
00096   LogTrace("MinBiasTracking") << fixed << setprecision(2)
00097     << "  [VZeroFinder] crossR = "<<crossing.perp()<<" ("
00098     <<minCrossingRadius<<"< <"<<maxCrossingRadius<<")";
00099 */
00100 
00101   if(dca < maxDca &&
00102      crossing.perp() > minCrossingRadius &&
00103      crossing.perp() < maxCrossingRadius)
00104   {
00105     // Momentum of the mother
00106     GlobalVector momentum = momenta.first + momenta.second;
00107     float impact = -1.;
00108 
00109     if(vertices->size() > 0)
00110     {
00111       // Impact parameter of the mother wrt vertices, choose smallest
00112       for(reco::VertexCollection::const_iterator
00113           vertex = vertices->begin(); vertex!= vertices->end(); vertex++)
00114       {
00115         GlobalVector r(crossing.x(),
00116                        crossing.y(),
00117                        crossing.z() - vertex->position().z());
00118         GlobalVector p(momentum.x(),momentum.y(),momentum.z());
00119   
00120         GlobalVector b = r - (r*p)*p / p.mag2();
00121   
00122         float im = b.mag();
00123         if(im < impact || vertex == vertices->begin())
00124           impact = im; 
00125       }
00126     }
00127     else
00128     {
00129       // Impact parameter of the mother in the plane
00130       GlobalVector r_(crossing.x(),crossing.y(),0);
00131       GlobalVector p_(momentum.x(),momentum.y(),0);
00132 
00133       GlobalVector b_ = r_ - (r_*p_)*p_ / p_.mag2();
00134       impact = b_.mag(); 
00135     }
00136 
00137 /*
00138     LogTrace("MinBiasTracking") << fixed << setprecision(2)
00139       << "  [VZeroFinder] impact = "<<impact<<" (<"<<maxImpactMother<<")";
00140 */
00141 
00142     if(impact < maxImpactMother)
00143     {
00144       // Armenteros variables
00145       float pt =
00146         (momenta.first.cross(momenta.second)).mag()/momentum.mag();
00147       float alpha =
00148         (momenta.first.mag2() - momenta.second.mag2())/momentum.mag2();
00149 
00150       // Fill data
00151       data.dca             = dca;
00152       data.crossingPoint   = crossing;
00153       data.impactMother    = impact;
00154       data.armenterosPt    = pt;
00155       data.armenterosAlpha = alpha;
00156       data.momenta         = momenta;
00157 
00158 /*
00159       LogTrace("MinBiasTracking") << fixed << setprecision(2)
00160       << "  [VZeroFinder] success {alpha = "<<alpha<<", pt = "<<pt<<"}";
00161 */
00162 
00163       return true;
00164     }
00165   }
00166 
00167   return false;
00168 }
00169