CMS 3D CMS Logo

PFTrackTransformer Class Reference

#include <RecoParticleFlow/PFTracking/interface/PFTrackTransformer.h>

List of all members.

Public Member Functions

bool addPoints (reco::PFRecTrack &pftrack, const reco::Track &track, const Trajectory &traj) const
 Add points to a PFTrack. return false if a TSOS is invalid.
bool addPointsAndBrems (reco::GsfPFRecTrack &pftrack, const reco::Track &track, const Trajectory &traj, const bool &GetMode) const
void OnlyProp ()
 PFTrackTransformer (math::XYZVector)
 ~PFTrackTransformer ()

Public Attributes

bool onlyprop_

Private Attributes

math::XYZVector B_
 B field.


Detailed Description

Definition at line 36 of file PFTrackTransformer.h.


Constructor & Destructor Documentation

PFTrackTransformer::PFTrackTransformer ( math::XYZVector  B  ) 

Definition at line 34 of file PFTrackTransformer.cc.

References onlyprop_.

00034                                                      :B_(B){
00035   LogInfo("PFTrackTransformer")<<"PFTrackTransformer built";
00036 
00037   PFGeometry pfGeometry;
00038   onlyprop_=false;
00039 }

PFTrackTransformer::~PFTrackTransformer (  ) 

Definition at line 41 of file PFTrackTransformer.cc.

00041                                        {
00042   
00043 }


Member Function Documentation

bool PFTrackTransformer::addPoints ( reco::PFRecTrack pftrack,
const reco::Track track,
const Trajectory traj 
) const

Add points to a PFTrack. return false if a TSOS is invalid.

Definition at line 47 of file PFTrackTransformer.cc.

References reco::PFTrack::addPoint(), reco::PFRecTrack::algoType(), alongMomentum, B_, PFGeometry::BeamPipe, reco::TrackBase::charge(), direction, Trajectory::direction(), BaseParticlePropagator::getSuccess(), reco::PFRecTrack::KF_ELCAND, LogDebug, PV3DBase< T, PVType, FrameType >::mag(), Trajectory::measurements(), RawParticle::momentum(), onlyprop_, reco::Track::outerMomentum(), reco::Track::outerPosition(), PFGeometry::outerRadius(), PFGeometry::outerZ(), reco::TrackBase::p(), p, BaseParticlePropagator::propagate(), reco::TrackBase::pt(), PT, reco::TrackBase::px(), reco::TrackBase::py(), reco::TrackBase::pz(), HcalSimpleRecAlgoImpl::reco(), RawParticle::setCharge(), BaseParticlePropagator::setPropagationConditions(), funct::sqrt(), std, v, RawParticle::vertex(), reco::TrackBase::vertex(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by PFConversionsProducer::fillPFConversions(), LightPFTrackProducer::produce(), GoodSeedProducer::produce(), and PFNuclearProducer::produce().

00049                                                              {
00050   
00051   LogDebug("PFTrackTransformer")<<"Trajectory propagation started";
00052   using namespace reco;
00053   using namespace std;
00054 
00055   float PT= track.pt();
00056   float pfmass= (pftrack.algoType()==reco::PFRecTrack::KF_ELCAND) ? 0.0005 : 0.139; 
00057   float pfenergy=sqrt((pfmass*pfmass)+(track.p()*track.p()));
00058    // closest approach
00059   BaseParticlePropagator theParticle = 
00060     BaseParticlePropagator( 
00061                            RawParticle(XYZTLorentzVector(track.px(),
00062                                                          track.py(),
00063                                                          track.pz(),
00064                                                          pfenergy),
00065                                        XYZTLorentzVector(track.vertex().x(),
00066                                                          track.vertex().y(),
00067                                                          track.vertex().z(),
00068                                                          0.)),
00069                            0.,0.,B_.z());
00070 
00071   theParticle.setCharge(track.charge());
00072   float pfoutenergy=sqrt((pfmass*pfmass)+track.outerMomentum().Mag2());
00073   BaseParticlePropagator theOutParticle = 
00074     BaseParticlePropagator( 
00075                            RawParticle(XYZTLorentzVector(track.outerMomentum().x(),
00076                                                          track.outerMomentum().y(),
00077                                                          track.outerMomentum().z(),
00078                                                          pfoutenergy),
00079                                        XYZTLorentzVector(track.outerPosition().x(),
00080                                                          track.outerPosition().y(),
00081                                                          track.outerPosition().z(),
00082                                                          0.)),
00083                            0.,0.,B_.z());
00084   theOutParticle.setCharge(track.charge());
00085   
00086   
00087   math::XYZTLorentzVector momClosest 
00088     = math::XYZTLorentzVector(track.px(), track.py(), 
00089                               track.pz(), track.p());
00090   math::XYZPoint posClosest = track.vertex();
00091   
00092   pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ClosestApproach,
00093                                      posClosest,momClosest));
00094   
00095   
00096   //BEAMPIPE
00097   theParticle.setPropagationConditions(PFGeometry::outerRadius(PFGeometry::BeamPipe), 
00098                                        PFGeometry::outerZ(PFGeometry::BeamPipe), false);
00099   theParticle.propagate();
00100   if(theParticle.getSuccess()!=0)
00101     pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::BeamPipeOrEndVertex,
00102                                        math::XYZPoint(theParticle.vertex()),
00103                                        math::XYZTLorentzVector(theParticle.momentum())));
00104   else {
00105     PFTrajectoryPoint dummyMaxSh;
00106     pftrack.addPoint(dummyMaxSh); 
00107   }
00108   
00109 
00110 
00111   //trajectory points
00112 
00113   if (!onlyprop_){
00114     bool direction =(traj.direction() == alongMomentum);
00115     vector<TrajectoryMeasurement> measurements =traj.measurements();
00116     int iTrajFirst = (direction) ? 0 :  measurements.size() - 1;
00117     int increment = (direction) ? +1 : -1;
00118     int iTrajLast  =  (direction) ? int(measurements.size()) : -1;
00119     
00120 
00121     for (int iTraj = iTrajFirst; iTraj != iTrajLast; iTraj += increment) {
00122       GlobalPoint v=measurements[iTraj].updatedState().globalPosition();
00123       GlobalVector p=measurements[iTraj].updatedState().globalMomentum();
00124       uint iid=measurements[iTraj].recHit()->det()->geographicalId().rawId();
00125       pftrack.addPoint(PFTrajectoryPoint(iid,-1,
00126                                          math::XYZPoint(v.x(), v.y(), v.z()),
00127                                          math::XYZTLorentzVector(p.x(),p.y(),p.z(),p.mag())));
00128     }
00129   }
00130 
00131   bool isBelowPS=false; 
00132    theOutParticle.propagateToPreshowerLayer1(false);
00133    if(theOutParticle.getSuccess()!=0)
00134      pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS1,
00135                                         math::XYZPoint(theOutParticle.vertex()),
00136                                         math::XYZTLorentzVector(theOutParticle.momentum())));
00137    else {
00138      PFTrajectoryPoint dummyPS1;
00139      pftrack.addPoint(dummyPS1); 
00140    }
00141    
00142 
00143    theOutParticle.propagateToPreshowerLayer2(false);
00144    if(theOutParticle.getSuccess()!=0){
00145      pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS2,
00146                                         math::XYZPoint(theOutParticle.vertex()),
00147                                         math::XYZTLorentzVector(theOutParticle.momentum())));
00148      isBelowPS=true;
00149    }   else {
00150      PFTrajectoryPoint dummyPS2;
00151      pftrack.addPoint(dummyPS2); 
00152    }
00153 
00154    theOutParticle.propagateToEcalEntrance(false);
00155 
00156    if(theOutParticle.getSuccess()!=0){
00157      pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALEntrance,
00158                                         math::XYZPoint(theOutParticle.vertex()),
00159                                         math::XYZTLorentzVector(theOutParticle.momentum())));
00160    double ecalShowerDepth     
00161      = PFCluster::getDepthCorrection(theOutParticle.momentum().E(),
00162                                      isBelowPS, 
00163                                      false);
00164 
00165    math::XYZPoint meanShower=math::XYZPoint(theOutParticle.vertex())+
00166      math::XYZTLorentzVector(theOutParticle.momentum()).Vect().Unit()*ecalShowerDepth;
00167  
00168    pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALShowerMax,
00169                                       meanShower,
00170                                       math::XYZTLorentzVector(theOutParticle.momentum())));}
00171    else {
00172      if (PT>5.)
00173        LogWarning("PFTrackTransformer")<<"KF TRACK "<<pftrack<< " PROPAGATION TO THE ECAL HAS FAILED";
00174      PFTrajectoryPoint dummyECAL;
00175      pftrack.addPoint(dummyECAL); 
00176      PFTrajectoryPoint dummyMaxSh;
00177      pftrack.addPoint(dummyMaxSh); 
00178    }
00179 
00180 
00181  
00182    //HCAL entrance
00183    theOutParticle.propagateToHcalEntrance(false);
00184    if(theOutParticle.getSuccess()!=0)
00185      pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALEntrance,
00186                                         math::XYZPoint(theOutParticle.vertex()),
00187                                         math::XYZTLorentzVector(theOutParticle.momentum())));
00188    else{
00189      if (PT>5.)
00190        LogWarning("PFTrackTransformer")<<"KF TRACK "<<pftrack<< " PROPAGATION TO THE HCAL ENTRANCE HAS FAILED";
00191      PFTrajectoryPoint dummyHCALentrance;
00192      pftrack.addPoint(dummyHCALentrance); 
00193    }
00194 
00195    //HCAL exit
00196    theOutParticle.propagateToHcalExit(false);
00197    if(theOutParticle.getSuccess()!=0)
00198      pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALExit,
00199                                         math::XYZPoint(theOutParticle.vertex()),
00200                                         math::XYZTLorentzVector(theOutParticle.momentum())));
00201    else{
00202      if (PT>5.)
00203        LogWarning("PFTrackTransformer")<<"KF TRACK "<<pftrack<< " PROPAGATION TO THE HCAL EXIT HAS FAILED";
00204      PFTrajectoryPoint dummyHCALexit;
00205      pftrack.addPoint(dummyHCALexit); 
00206    }
00207    
00208    return true;
00209 }

bool PFTrackTransformer::addPointsAndBrems ( reco::GsfPFRecTrack pftrack,
const reco::Track track,
const Trajectory traj,
const bool GetMode 
) const

Definition at line 211 of file PFTrackTransformer.cc.

References funct::abs(), reco::GsfPFRecTrack::addBrem(), reco::PFTrack::addPoint(), reco::PFRecTrack::algoType(), alongMomentum, B_, PFGeometry::BeamPipe, reco::TrackBase::charge(), PFGsfHelper::computeP(), direction, Trajectory::direction(), PFGsfHelper::fittedDP(), BaseParticlePropagator::getSuccess(), reco::PFRecTrack::KF_ELCAND, PV3DBase< T, PVType, FrameType >::mag(), Trajectory::measurements(), RawParticle::momentum(), PFGeometry::outerRadius(), PFGeometry::outerZ(), p, PV3DBase< T, PVType, FrameType >::perp(), BaseParticlePropagator::propagate(), BaseParticlePropagator::propagateToEcalEntrance(), BaseParticlePropagator::propagateToHcalEntrance(), BaseParticlePropagator::propagateToHcalExit(), BaseParticlePropagator::propagateToPreshowerLayer1(), BaseParticlePropagator::propagateToPreshowerLayer2(), reco::TrackBase::pt(), PT, RawParticle::setCharge(), BaseParticlePropagator::setPropagationConditions(), PFGsfHelper::sigmafittedDP(), funct::sqrt(), v, RawParticle::vertex(), reco::TrackBase::vertex(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by PFElecTkProducer::produce().

00214                                                                   {
00215 
00216   float PT= track.pt();
00217   // Trajectory for each trajectory point
00218 
00219   bool direction =(traj.direction() == alongMomentum);
00220   vector<TrajectoryMeasurement> measurements =traj.measurements();
00221   int iTrajFirst = (direction) ? 0 :  measurements.size() - 1;
00222   int increment = (direction) ? +1 : -1;
00223   int iTrajLast  =  (direction) ? int(measurements.size()) : -1;
00224   
00225   
00226   uint iTrajPos = 0;
00227   for (int iTraj = iTrajFirst; iTraj != iTrajLast; iTraj += increment) {
00228     
00229     GlobalPoint v=measurements[iTraj].updatedState().globalPosition();
00230     PFGsfHelper* PFGsf = new PFGsfHelper(measurements[iTraj]);
00231     //if (PFGsf->isValid()){ 
00232     bool ComputeMODE = GetMode;
00233     GlobalVector p = PFGsf->computeP(ComputeMODE);
00234     double DP = PFGsf->fittedDP();
00235     double SigmaDP =  PFGsf->sigmafittedDP();   
00236     uint iid=measurements[iTraj].recHit()->det()->geographicalId().rawId();
00237     delete PFGsf;
00238 
00239     // --------------------------   Fill GSF Track ------------------------------------- 
00240     
00241 
00242     float pfmass= (pftrack.algoType()==reco::PFRecTrack::KF_ELCAND) ? 0.0005 : 0.139; 
00243     float ptot =  sqrt((p.x()*p.x())+(p.y()*p.y())+(p.z()*p.z()));
00244     float pfenergy=sqrt((pfmass*pfmass)+(ptot *ptot));
00245 
00246     if (iTraj == iTrajFirst) {
00247 
00248       math::XYZTLorentzVector momClosest 
00249         = math::XYZTLorentzVector(p.x(), p.y(), 
00250                                   p.z(), ptot);
00251       math::XYZPoint posClosest = track.vertex();
00252       pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ClosestApproach,
00253                                          posClosest,momClosest));
00254       
00255       BaseParticlePropagator theInnerParticle = 
00256         BaseParticlePropagator( 
00257                                RawParticle(XYZTLorentzVector(p.x(),
00258                                                              p.y(),
00259                                                              p.z(),
00260                                                              pfenergy),
00261                                            XYZTLorentzVector(track.vertex().x(),
00262                                                              track.vertex().y(),
00263                                                              track.vertex().z(),
00264                                                              0.)),  //DANIELE Same thing v.x(),v.y(),v.()? 
00265                                0.,0.,B_.z());
00266       theInnerParticle.setCharge(track.charge());  
00267 
00268       //BEAMPIPE
00269       theInnerParticle.setPropagationConditions(PFGeometry::outerRadius(PFGeometry::BeamPipe), 
00270                                            PFGeometry::outerZ(PFGeometry::BeamPipe), false);
00271       theInnerParticle.propagate();
00272       if(theInnerParticle.getSuccess()!=0)
00273         pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::BeamPipeOrEndVertex,
00274                                            math::XYZPoint(theInnerParticle.vertex()),
00275                                            math::XYZTLorentzVector(theInnerParticle.momentum())));
00276       else {
00277         PFTrajectoryPoint dummyMaxSh;
00278         pftrack.addPoint(dummyMaxSh); 
00279       }
00280       
00281       // First Point for the trajectory == Vertex ?? 
00282       pftrack.addPoint(PFTrajectoryPoint(iid,-1,
00283                                          math::XYZPoint(v.x(), v.y(), v.z()),
00284                                          math::XYZTLorentzVector(p.x(),p.y(),p.z(),p.mag())));
00285       
00286      
00287     }
00288     if (iTraj != iTrajFirst && iTraj != (abs(iTrajLast)-1)) {
00289       pftrack.addPoint(PFTrajectoryPoint(iid,-1,
00290                                          math::XYZPoint(v.x(), v.y(), v.z()),
00291                                          math::XYZTLorentzVector(p.x(),p.y(),p.z(),p.mag())));
00292       
00293  
00294     }
00295     if (iTraj == (abs(iTrajLast)-1)) {
00296       
00297       // Last Trajectory Meas
00298       pftrack.addPoint(PFTrajectoryPoint(iid,-1,
00299                                          math::XYZPoint(v.x(), v.y(), v.z()),
00300                                          math::XYZTLorentzVector(p.x(),p.y(),p.z(),p.mag())));
00301 
00302 
00303 
00304 
00305       BaseParticlePropagator theOutParticle = 
00306         BaseParticlePropagator( 
00307                                RawParticle(XYZTLorentzVector(p.x(),
00308                                                              p.y(),
00309                                                              p.z(),
00310                                                              pfenergy),
00311                                            XYZTLorentzVector(v.x(),
00312                                                              v.y(),
00313                                                              v.z(),
00314                                                              0.)), 
00315                                0.,0.,B_.z());
00316       theOutParticle.setCharge(track.charge());  
00317       bool isBelowPS=false; 
00318       theOutParticle.propagateToPreshowerLayer1(false);
00319       if(theOutParticle.getSuccess()!=0)
00320         pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS1,
00321                                            math::XYZPoint(theOutParticle.vertex()),
00322                                            math::XYZTLorentzVector(theOutParticle.momentum())));
00323       else {
00324         PFTrajectoryPoint dummyPS1;
00325         pftrack.addPoint(dummyPS1); 
00326       }
00327       
00328       
00329       theOutParticle.propagateToPreshowerLayer2(false);
00330       if(theOutParticle.getSuccess()!=0){
00331         pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS2,
00332                                            math::XYZPoint(theOutParticle.vertex()),
00333                                            math::XYZTLorentzVector(theOutParticle.momentum())));
00334         isBelowPS=true;
00335       }   else {
00336         PFTrajectoryPoint dummyPS2;
00337         pftrack.addPoint(dummyPS2); 
00338       }
00339       
00340       theOutParticle.propagateToEcalEntrance(false);
00341       
00342       if(theOutParticle.getSuccess()!=0){
00343         pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALEntrance,
00344                                            math::XYZPoint(theOutParticle.vertex()),
00345                                            math::XYZTLorentzVector(theOutParticle.momentum())));
00346         double ecalShowerDepth     
00347           = PFCluster::getDepthCorrection(theOutParticle.momentum().E(),
00348                                           isBelowPS, 
00349                                           false);
00350         
00351         math::XYZPoint meanShower=math::XYZPoint(theOutParticle.vertex())+
00352           math::XYZTLorentzVector(theOutParticle.momentum()).Vect().Unit()*ecalShowerDepth;
00353         
00354         pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALShowerMax,
00355                                            meanShower,
00356                                            math::XYZTLorentzVector(theOutParticle.momentum())));}
00357       else {
00358         if (PT>5.)
00359           LogWarning("PFTrackTransformer")<<"GSF TRACK "<<pftrack<< " PROPAGATION TO THE ECAL HAS FAILED";
00360         PFTrajectoryPoint dummyECAL;
00361         pftrack.addPoint(dummyECAL); 
00362         PFTrajectoryPoint dummyMaxSh;
00363         pftrack.addPoint(dummyMaxSh); 
00364       }
00365       
00366       
00367       
00368       //HCAL entrance
00369       theOutParticle.propagateToHcalEntrance(false);
00370       if(theOutParticle.getSuccess()!=0)
00371         pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALEntrance,
00372                                            math::XYZPoint(theOutParticle.vertex()),
00373                                            math::XYZTLorentzVector(theOutParticle.momentum())));
00374       else{
00375         if (PT>5.)
00376           LogWarning("PFTrackTransformer")<<"GSF TRACK "<<pftrack<< " PROPAGATION TO THE HCAL ENTRANCE HAS FAILED";
00377         PFTrajectoryPoint dummyHCALentrance;
00378         pftrack.addPoint(dummyHCALentrance); 
00379       }  
00380       
00381       //HCAL exit
00382       theOutParticle.propagateToHcalExit(false);
00383       if(theOutParticle.getSuccess()!=0)
00384         pftrack.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALExit,
00385                                            math::XYZPoint(theOutParticle.vertex()),
00386                                            math::XYZTLorentzVector(theOutParticle.momentum())));
00387       else{
00388         if (PT>5.)
00389           LogWarning("PFTrackTransformer")<<"GSF TRACK "<<pftrack<< " PROPAGATION TO THE HCAL EXIT HAS FAILED";
00390         PFTrajectoryPoint dummyHCALexit;
00391         pftrack.addPoint(dummyHCALexit); 
00392       } 
00393       
00394     }
00395 
00396     // --------------------------   END GSF Track ------------------------------------- 
00397   
00398     // --------------------------   Fill Brem "Track" --------------------------------- 
00399     // Fill the brem for each traj point
00400 
00401     //check that the vertex of the brem is in the tracker volume
00402     if ((v.perp()>110) ||(fabs(v.z())>280)) continue;    
00403     uint iTrajPoint =  iTrajPos + 2;
00404 
00405     PFBrem brem(DP,SigmaDP,iTrajPoint);
00406 //     cout << " DP " << DP <<  " Sigma " << SigmaDP << endl;   // Daniele: Remove
00407 
00408 
00409     GlobalVector p_gamma= p*(fabs(DP)/p.mag());   // Direction from the electron (tangent), DP without any sign!;
00410     float e_gamma = fabs(DP); // DP = pout-pin so could be negative
00411     BaseParticlePropagator theBremParticle = 
00412       BaseParticlePropagator( 
00413                              RawParticle(XYZTLorentzVector(p_gamma.x(),
00414                                                            p_gamma.y(),
00415                                                            p_gamma.z(),
00416                                                            e_gamma),
00417                                          XYZTLorentzVector(v.x(),
00418                                                            v.y(),
00419                                                            v.z(),
00420                                                            0.)),
00421                              0.,0.,B_.z());
00422     int gamma_charge = 0;
00423     theBremParticle.setCharge(gamma_charge);  
00424 
00425 
00426     // add TrajectoryPoint for Brem, PS, ECAL, ECALShowMax, HCAL
00427     // Brem Entrance PS Layer1
00428 
00429     PFTrajectoryPoint dummyClosest;   // Added just to have the right number order in PFTrack.cc
00430     brem.addPoint(dummyClosest); 
00431 
00432     
00433     PFTrajectoryPoint dummyBeamPipe;  // Added just to have the right number order in PFTrack.cc
00434     brem.addPoint(dummyBeamPipe); 
00435 
00436 
00437     
00438     bool isBelowPS=false; 
00439     theBremParticle.propagateToPreshowerLayer1(false);
00440     if(theBremParticle.getSuccess()!=0)
00441       brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS1,
00442                                          math::XYZPoint(theBremParticle.vertex()),
00443                                          math::XYZTLorentzVector(theBremParticle.momentum())));
00444     else {
00445       PFTrajectoryPoint dummyPS1;
00446       brem.addPoint(dummyPS1); 
00447     }
00448     
00449     // Brem Entrance PS Layer 2
00450 
00451     theBremParticle.propagateToPreshowerLayer2(false);
00452     if(theBremParticle.getSuccess()!=0){
00453       brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::PS2,
00454                                          math::XYZPoint(theBremParticle.vertex()),
00455                                          math::XYZTLorentzVector(theBremParticle.momentum())));
00456       isBelowPS=true;
00457     }   else {
00458       PFTrajectoryPoint dummyPS2;
00459       brem.addPoint(dummyPS2); 
00460     }
00461 
00462    theBremParticle.propagateToEcalEntrance(false);
00463 
00464    if(theBremParticle.getSuccess()!=0){
00465      brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALEntrance,
00466                                         math::XYZPoint(theBremParticle.vertex()),
00467                                         math::XYZTLorentzVector(theBremParticle.momentum())));
00468    double ecalShowerDepth     
00469      = PFCluster::getDepthCorrection(theBremParticle.momentum().E(),
00470                                      isBelowPS, 
00471                                      false);
00472 
00473    math::XYZPoint meanShower=math::XYZPoint(theBremParticle.vertex())+
00474      math::XYZTLorentzVector(theBremParticle.momentum()).Vect().Unit()*ecalShowerDepth;
00475  
00476    brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::ECALShowerMax,
00477                                       meanShower,
00478                                       math::XYZTLorentzVector(theBremParticle.momentum())));}
00479    else {
00480      if ((DP>5.) && ((DP/SigmaDP)>3))
00481        LogWarning("PFTrackTransformer")<<"BREM "<<brem<<" PROPAGATION TO THE ECAL HAS FAILED";
00482      PFTrajectoryPoint dummyECAL;
00483      brem.addPoint(dummyECAL); 
00484      PFTrajectoryPoint dummyMaxSh;
00485      brem.addPoint(dummyMaxSh); 
00486    }
00487 
00488 
00489  
00490    //HCAL entrance
00491    theBremParticle.propagateToHcalEntrance(false);
00492    if(theBremParticle.getSuccess()!=0)
00493      brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALEntrance,
00494                                         math::XYZPoint(theBremParticle.vertex()),
00495                                         math::XYZTLorentzVector(theBremParticle.momentum())));
00496    else{
00497      if ((DP>5.) && ((DP/SigmaDP)>3))
00498        LogWarning("PFTrackTransformer")<<"BREM "<<brem<<" PROPAGATION TO THE HCAL ENTRANCE HAS FAILED";
00499      PFTrajectoryPoint dummyHCALentrance;
00500      pftrack.addPoint(dummyHCALentrance); 
00501    }  
00502 
00503    //HCAL exit
00504    theBremParticle.propagateToHcalExit(false);
00505    if(theBremParticle.getSuccess()!=0)
00506      brem.addPoint(PFTrajectoryPoint(-1,PFTrajectoryPoint::HCALExit,
00507                                         math::XYZPoint(theBremParticle.vertex()),
00508                                         math::XYZTLorentzVector(theBremParticle.momentum())));
00509    else{  
00510      if ((DP>5.) && ((DP/SigmaDP)>3))
00511        LogWarning("PFTrackTransformer")<<"BREM "<<brem<<" PROPAGATION TO THE HCAL EXIT HAS FAILED";
00512      PFTrajectoryPoint dummyHCALexit;
00513      pftrack.addPoint(dummyHCALexit); 
00514    }
00515 
00516    pftrack.addBrem(brem);
00517    iTrajPos++;
00518   }
00519   return true;
00520 }

void PFTrackTransformer::OnlyProp (  )  [inline]

Definition at line 56 of file PFTrackTransformer.h.

References onlyprop_.

Referenced by LightPFTrackProducer::beginJob(), and PFNuclearProducer::beginJob().

00056                  {
00057     onlyprop_=true;
00058   }


Member Data Documentation

math::XYZVector PFTrackTransformer::B_ [private]

B field.

Definition at line 63 of file PFTrackTransformer.h.

Referenced by addPoints(), and addPointsAndBrems().

bool PFTrackTransformer::onlyprop_

Definition at line 59 of file PFTrackTransformer.h.

Referenced by addPoints(), OnlyProp(), and PFTrackTransformer().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:29:50 2009 for CMSSW by  doxygen 1.5.4