CMS 3D CMS Logo

CastorCluster.h
Go to the documentation of this file.
1 #ifndef CastorReco_CastorCluster_h
2 #define CastorReco_CastorCluster_h
3 
12 #include <vector>
13 #include <memory>
15 
19 
21 
22 namespace reco {
23 
24  class CastorCluster {
25  public:
28  : energy_(0.),
29  position_(ROOT::Math::XYZPoint(0., 0., 0.)),
30  emEnergy_(0.),
31  hadEnergy_(0.),
32  fem_(0.),
33  width_(0.),
34  depth_(0.),
35  fhot_(0.),
36  sigmaz_(0.) {}
37 
39  CastorCluster(const double energy,
41  const double emEnergy,
42  const double hadEnergy,
43  const double fem,
44  const double width,
45  const double depth,
46  const double fhot,
47  const double sigmaz,
48  const CastorTowerRefVector& usedTowers);
49 
51  virtual ~CastorCluster();
52 
54  double energy() const { return energy_; }
55 
58 
60  double emEnergy() const { return emEnergy_; }
61 
63  double hadEnergy() const { return hadEnergy_; }
64 
66  double fem() const { return fem_; }
67 
69  double width() const { return width_; }
70 
72  double depth() const { return depth_; }
73 
75  double fhot() const { return fhot_; }
76 
78  double sigmaz() const { return sigmaz_; }
79 
82 
85 
88 
90  size_t towersSize() const { return usedTowers_.size(); }
91 
94 
96  bool operator>=(const CastorCluster& rhs) const { return (energy_ >= rhs.energy_); }
97 
99  bool operator>(const CastorCluster& rhs) const { return (energy_ > rhs.energy_); }
100 
102  bool operator<=(const CastorCluster& rhs) const { return (energy_ <= rhs.energy_); }
103 
105  bool operator<(const CastorCluster& rhs) const { return (energy_ < rhs.energy_); }
106 
108  double eta() const { return position_.eta(); }
109 
111  double phi() const { return position_.phi(); }
112 
114  double x() const { return position_.x(); }
115 
117  double y() const { return position_.y(); }
118 
120  double rho() const { return position_.rho(); }
121 
122  private:
124  double energy_;
125 
128 
130  double emEnergy_;
131 
133  double hadEnergy_;
134 
136  double fem_;
137 
139  double width_;
140 
142  double depth_;
143 
145  double fhot_;
146 
148  double sigmaz_;
149 
152  };
153 
155  typedef std::vector<CastorCluster> CastorClusterCollection;
156 
157  // persistent reference to CastorCluster objects
159 
162 
165 } // namespace reco
166 
167 #endif
virtual ~CastorCluster()
destructor
void add(const CastorTowerRef &tower)
add reference to constituent CastorTower
Definition: CastorCluster.h:93
double hadEnergy() const
cluster had energy
Definition: CastorCluster.h:63
CastorTowerRefVector usedTowers_
references to CastorTower constituents
double width_
cluster width
bool operator<(const CastorCluster &rhs) const
comparison <= operator
double rho() const
rho of cluster centroid
double energy() const
cluster energy
Definition: CastorCluster.h:54
edm::RefVector< CastorClusterCollection > CastorClusterRefVector
vector of references to CastorCluster objects all in the same collection
double fhot() const
cluster hotcell/tot ratio
Definition: CastorCluster.h:75
CastorCluster()
default constructor. Sets energy and position to zero
Definition: CastorCluster.h:27
double emEnergy() const
cluster em energy
Definition: CastorCluster.h:60
std::vector< CastorCluster > CastorClusterCollection
collection of CastorCluster objects
double sigmaz_
cluster sigma z
CastorTowerRefVector getUsedTowers() const
vector of used Towers
Definition: CastorCluster.h:81
double energy_
cluster energy
bool operator<=(const CastorCluster &rhs) const
comparison <= operator
CastorTower_iterator towersEnd() const
last iterator over CastorTower constituents
Definition: CastorCluster.h:87
edm::Ref< CastorClusterCollection > CastorClusterRef
double eta() const
pseudorapidity of cluster centroid
CastorClusterRefVector::iterator CastorCluster_iterator
iterator over a vector of references to CastorCluster objects all in the same collection ...
double depth_
cluster depth
double fem_
cluster em/tot Ratio
ROOT::Math::XYZPoint position_
cluster centroid position
double depth() const
cluster depth in z
Definition: CastorCluster.h:72
double emEnergy_
cluster em energy
double sigmaz() const
cluster sigma z
Definition: CastorCluster.h:78
double fem() const
cluster em/tot ratio
Definition: CastorCluster.h:66
CastorTower_iterator towersBegin() const
fist iterator over CastorTower constituents
Definition: CastorCluster.h:84
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
double hadEnergy_
cluster had energy
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228
ROOT::Math::XYZPoint position() const
cluster centroid position
Definition: CastorCluster.h:57
bool operator>=(const CastorCluster &rhs) const
comparison >= operator
Definition: CastorCluster.h:96
bool operator>(const CastorCluster &rhs) const
comparison > operator
Definition: CastorCluster.h:99
size_t towersSize() const
number of CastorTower constituents
Definition: CastorCluster.h:90
double width() const
cluster width in phi
Definition: CastorCluster.h:69
fixed size matrix
Transform3DPJ::Point XYZPoint
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67
double y() const
y of cluster centroid
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223
double fhot_
cluster hotcell/tot ratio
double x() const
x of cluster centroid
double phi() const
azimuthal angle of cluster centroid