ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitEdm/Producers/interface/VertexZProducer.h
Revision: 1.1
Committed: Mon Nov 30 10:11:40 2009 UTC (15 years, 5 months ago) by loizides
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_032, Mit_031, Mit_025c_branch2, Mit_025c_branch1, Mit_030, Mit_029c, Mit_029b, Mit_030_pre1, Mit_029a, Mit_029, Mit_029_pre1, Mit_028a, Mit_025c_branch0, Mit_028, Mit_027a, Mit_027, Mit_026, Mit_025e, Mit_025d, Mit_025c, Mit_025b, Mit_025a, Mit_025, Mit_025pre2, Mit_024b, Mit_025pre1, Mit_024a, Mit_024, Mit_023, Mit_022a, Mit_022, Mit_020d, TMit_020d, Mit_020c, Mit_021, Mit_021pre2, Mit_021pre1, Mit_020b, Mit_020a, Mit_020, Mit_020pre1, Mit_018, Mit_017, Mit_017pre3, Mit_017pre2, Mit_017pre1, V07-05-00, Mit_016, Mit_015b, Mit_015a, Mit_015, Mit_014e, Mit_014d, Mit_014c, Mit_014b, ConvRejection-10-06-09, Mit_014a, Mit_014, Mit_014pre3, Mit_014pre2, Mit_014pre1, Mit_013d, Mit_013c, Mit_013b, Mit_013a, Mit_013, Mit_013pre1, Mit_012i, Mit_012h, Mit_012g, Mit_012f, Mit_012e, Mit_012d, Mit_012c, HEAD
Branch point for: Mit_025c_branch
Log Message:
Added vertex z producer.

File Contents

# Content
1 // $Id: QcdLowPtDQM.h,v 1.9 2009/11/19 22:33:15 loizides Exp $
2
3 #ifndef MITEDM_PRODUCERS_VERTEXZPRODUCER_H
4 #define MITEDM_PRODUCERS_VERTEXZPRODUCER_H
5
6 #include "DataFormats/GeometryVector/interface/VectorUtil.h"
7 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
8 #include "FWCore/Framework/interface/Frameworkfwd.h"
9 #include "FWCore/Framework/interface/EDProducer.h"
10 #include "FWCore/Framework/interface/Event.h"
11 #include "FWCore/MessageLogger/interface/MessageLogger.h"
12 #include <TMath.h>
13 #include <vector>
14
15 class TrackerGeometry;
16
17 class VertexZProducer : public edm::EDProducer
18 {
19 public:
20 class Pixel {
21 public:
22 Pixel(double x=0, double y=0, double z=0, double eta=0, double phi=0,
23 double adc=0, double sx=0, double sy=0) :
24 x_(x), y_(y), z_(z), rho_(TMath::Sqrt(x_*x_+y_*y_)),
25 eta_(eta), phi_(phi), adc_(adc), sizex_(sx), sizey_(sy) {}
26 Pixel(const GlobalPoint &p, double adc=0, double sx=0, double sy=0) :
27 x_(p.x()), y_(p.y()), z_(p.z()), rho_(TMath::Sqrt(x_*x_+y_*y_)),
28 eta_(p.eta()), phi_(p.phi()), adc_(adc), sizex_(sx), sizey_(sy) {}
29 double adc() const { return adc_; }
30 double eta() const { return eta_; }
31 double rho() const { return rho_; }
32 double phi() const { return phi_; }
33 double sizex() const { return sizex_; }
34 double sizey() const { return sizey_; }
35 double x() const { return x_; }
36 double y() const { return y_; }
37 double z() const { return z_; }
38 protected:
39 double x_,y_,z_,rho_,eta_,phi_;
40 double adc_,sizex_,sizey_;
41 };
42 class Tracklet {
43 public:
44 Tracklet() : i1_(-1), i2_(-2), deta_(0), dphi_(0) {}
45 Tracklet(const Pixel &p1, const Pixel &p2) :
46 p1_(p1), p2_(p2), i1_(-1), i2_(-1),
47 deta_(p1.eta()-p2.eta()), dphi_(Geom::deltaPhi(p1.phi(),p2.phi())) {}
48 double deta() const { return deta_; }
49 double dphi() const { return dphi_; }
50 int i1() const { return i1_; }
51 int i2() const { return i2_; }
52 double eta() const { return p1_.eta(); }
53 const Pixel &p1() const { return p1_; }
54 const Pixel &p2() const { return p2_; }
55 void seti1(int i1) { i1_ = i1; }
56 void seti2(int i2) { i2_ = i2; }
57 protected:
58 Pixel p1_,p2_;
59 int i1_, i2_;
60 double deta_,dphi_;
61 };
62 class Vertex {
63 public:
64 Vertex(double x=0, double y=0, double z=0,
65 double xs=0, double ys=0, double zs=0,
66 double chi2=0, int n=0) :
67 x_(x), y_(y), z_(z), xs_(xs), ys_(ys), zs_(zs), chi2_(chi2), n_(n) {}
68 double chi2() const { return chi2_; }
69 int n() const { return n_; }
70 double x() const { return x_; }
71 double y() const { return y_; }
72 double z() const { return z_; }
73 double xs() const { return xs_; }
74 double ys() const { return ys_; }
75 double zs() const { return zs_; }
76 void set(int n, double z, double zs) { n_= n; z_ = z; zs_ = zs; }
77 void set(int n, double x,double y,double z, double xs,double ys,double zs)
78 { n_= n; x_ = x; xs_ = xs; y_ = y; ys_ = ys; z_ = z; zs_ = zs; }
79 void setchi2(double chi2) { chi2_ = chi2; }
80 void setinvalid() { n_ = 0; z_=-9999; zs_ = 0; chi2_ = 0; }
81 protected:
82 double x_,y_,z_,xs_,ys_,zs_,chi2_;
83 int n_;
84 };
85
86 explicit VertexZProducer(const edm::ParameterSet &parameters);
87 ~VertexZProducer();
88
89 void produce(edm::Event &iEvent, const edm::EventSetup &iSetup);
90
91 private:
92 void fillPixels(const edm::Event &iEvent);
93 template <typename TYPE>
94 void getProduct(const std::string name, edm::Handle<TYPE> &prod,
95 const edm::Event &event) const;
96 template <typename TYPE>
97 bool getProductSafe(const std::string name, edm::Handle<TYPE> &prod,
98 const edm::Event &event) const;
99 void print(int level, const char *msg);
100 void print(int level, const std::string &msg)
101 { print(level,msg.c_str()); }
102 void reallyPrint(int level, const char *msg);
103 void trackletVertexUnbinned(const edm::Event &iEvent, int which=12);
104 void trackletVertexUnbinned(std::vector<Pixel> &pix1,
105 std::vector<Pixel> &pix2,
106 Vertex &vtx);
107 void vertexFromClusters(const edm::Event &iEvent, int which=12);
108 void vertexFromClusters(const std::vector<Pixel> &pix, Vertex &vtx);
109
110 std::string pixelName_; //pixel reconstructed hits name
111 double dPhiVc_; //dPhi vertex cut for tracklet based vertex
112 double dZVc_; //dZ vertex cut for tracklet based vertex
113 int verbose_; //verbosity (0=debug,1=warn,2=error,3=throw)
114 int pixLayers_; //either 12,13,23
115 bool doClusVertex_; //if true run cluster vertex finder
116 bool useRecHitQ_; //if true use rec hit quality word
117 bool usePixelQ_; //if true use pixel hit quality word
118 std::vector<Pixel> bpix1_; //barrel pixels layer 1
119 std::vector<Pixel> bpix2_; //barrel pixels layer 2
120 std::vector<Pixel> bpix3_; //barrel pixels layer 3
121 Vertex trackletV_; //reconstructed tracklet vertex 12
122 const TrackerGeometry *tgeo_; //tracker geometry
123 };
124
125 //--------------------------------------------------------------------------------------------------
126 template <typename TYPE>
127 inline void VertexZProducer::getProduct(const std::string name, edm::Handle<TYPE> &prod,
128 const edm::Event &event) const
129 {
130 // Try to access data collection from EDM file. We check if we really get just one
131 // product with the given name. If not we throw an exception.
132
133 event.getByLabel(edm::InputTag(name),prod);
134 if (!prod.isValid())
135 throw edm::Exception(edm::errors::Configuration, "VertexZProducer::GetProduct()\n")
136 << "Collection with label " << name << " is not valid" << std::endl;
137 }
138
139 //--------------------------------------------------------------------------------------------------
140 template <typename TYPE>
141 inline bool VertexZProducer::getProductSafe(const std::string name, edm::Handle<TYPE> &prod,
142 const edm::Event &event) const
143 {
144 // Try to safely access data collection from EDM file. We check if we really get just one
145 // product with the given name. If not, we return false.
146
147 if (name.size()==0)
148 return false;
149
150 try {
151 event.getByLabel(edm::InputTag(name),prod);
152 if (!prod.isValid())
153 return false;
154 } catch (...) {
155 return false;
156 }
157 return true;
158 }
159
160 //--------------------------------------------------------------------------------------------------
161 inline void VertexZProducer::print(int level, const char *msg)
162 {
163 // Print out message if it
164
165 if (level>=verbose_)
166 reallyPrint(level,msg);
167 }
168 #endif