CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/DataFormats/ParticleFlowReco/interface/PFCluster.h

Go to the documentation of this file.
00001 #ifndef DataFormats_ParticleFlowReco_PFCluster_h
00002 #define DataFormats_ParticleFlowReco_PFCluster_h
00003 
00004 #include "DataFormats/CaloRecHit/interface/CaloCluster.h"
00005 
00006 #include "Math/GenVector/PositionVector3D.h"
00007 #include "DataFormats/Math/interface/Point3D.h"
00008 #include "Rtypes.h" 
00009 
00010 #include "DataFormats/ParticleFlowReco/interface/PFRecHitFraction.h"
00011 #include "DataFormats/ParticleFlowReco/interface/PFLayer.h"
00012 
00013 #include <iostream>
00014 #include <vector>
00015 
00016 
00017 
00018 class PFClusterAlgo;
00019 
00020 namespace reco {
00021 
00042   class PFCluster : public CaloCluster {
00043   public:
00044 
00045 
00046     typedef ROOT::Math::PositionVector3D<ROOT::Math::CylindricalEta3D<Double32_t> > REPPoint;
00047   
00048     PFCluster() : CaloCluster(CaloCluster::particleFlow), color_(1) {}
00049 
00051     PFCluster(PFLayer::Layer layer, double energy,
00052               double x, double y, double z );
00053 
00054 
00055 
00057     void reset();
00058     
00060     void addRecHitFraction( const reco::PFRecHitFraction& frac);
00061     
00063     const std::vector< reco::PFRecHitFraction >& recHitFractions() const 
00064       { return rechits_; }
00065     
00067     void setLayer( PFLayer::Layer layer);
00068     
00070     PFLayer::Layer  layer() const;     
00071     
00073     double        energy() const {return energy_;}
00074     
00076     const REPPoint&       positionREP() const {return posrep_;}
00077     
00079     void calculatePositionREP() {
00080       posrep_.SetCoordinates( position_.Rho(), 
00081                               position_.Eta(), 
00082                               position_.Phi() ); 
00083     }
00084     
00086     static double getDepthCorrection(double energy, bool isBelowPS = false,
00087                                      bool isHadron = false);
00088     
00090     void         setColor(int color) {color_ = color;}
00091     
00093     int          color() const {return color_;}
00094     
00095     
00096     PFCluster& operator=(const PFCluster&);
00097     
00098     friend    std::ostream& operator<<(std::ostream& out, 
00099                                        const PFCluster& cluster);
00101     static unsigned     instanceCounter_;
00102     
00104     static void setDepthCorParameters(int mode, 
00105                                       double a, double b, 
00106                                       double ap, double bp ) {
00107       depthCorMode_ = mode;
00108       depthCorA_ = a; 
00109       depthCorB_ = b; 
00110       depthCorAp_ = ap; 
00111       depthCorBp_ = bp; 
00112     } 
00113     
00114 
00118     
00120     double charge() const { return 0;}
00121 
00123     double pt() const { 
00124       return (energy() * sin(position_.theta()));
00125     }
00126 
00128     double theta() const { 
00129       return position_.theta();
00130     }
00131     
00133     math::XYZPoint const & vertex() const { 
00134       static math::XYZPoint dummyVtx(0,0,0);
00135       return dummyVtx;      
00136     }
00137     double vx() const { return vertex().x(); }
00138     double vy() const { return vertex().y(); }
00139     double vz() const { return vertex().z(); }    
00140 
00141   private:
00142     
00144     std::vector< reco::PFRecHitFraction >  rechits_;
00145     
00147     REPPoint            posrep_;
00148     
00149     
00151     static int    depthCorMode_;
00152     
00154     static double depthCorA_;
00155     
00157     static double depthCorB_ ;
00158     
00160     static double depthCorAp_;
00161     
00163     static double depthCorBp_;
00164     
00165     
00167     int                 color_;
00168     
00169     friend class ::PFClusterAlgo;
00170   };
00171 }
00172 
00173 #endif