CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
METAlgo.cc
Go to the documentation of this file.
1 // File: METAlgo.cc
2 // Description: see METAlgo.h
3 // Author: Michael Schmitt, Richard Cavanaugh The University of Florida
4 // Creation Date: MHS May 31, 2005 Initial version.
5 //
6 //------------------------------------------------------------------------
7 
12 #include <iostream>
13 
14 using namespace std;
15 using namespace reco;
16 
17 //------------------------------------------------------------------------
18 // Default Constructer
19 //----------------------------------
21 //------------------------------------------------------------------------
22 
23 //------------------------------------------------------------------------
24 // Default Destructor
25 //----------------------------------
27 //------------------------------------------------------------------------
28 
29 //------------------------------------------------------------------------
30 // This method represents "the" implementation of the MET algorithm and is
31 // very simple:
32 // (1) It takes as input a collection of candidates (which can be
33 // calorimeter towers, HEPMC generator-level particles, etc).
34 // (2) It returns as output, a pointer to a struct of CommonMETData which
35 // contains the following members: MET, MEx, MEy, SumET, and MEz
36 // (The inclusion of MEz deserves some justification ; it is included here
37 // since it _may_ be useful for Data Quality Monitering as it should be
38 // symmetrically distributed about the origin.)
39 //----------------------------------
40 //void METAlgo::run(const CandidateCollection *input, CommonMETData *met, double globalThreshold)
42 {
43  double sum_et = 0.0;
44  double sum_ex = 0.0;
45  double sum_ey = 0.0;
46  double sum_ez = 0.0;
47  // Loop over Candidate Objects and calculate MET and related quantities
48  /*
49  CandidateCollection::const_iterator candidate;
50  for( candidate = input->begin(); candidate != input->end(); candidate++ )
51  */
52  for (unsigned int candidate_i = 0; candidate_i < input->size(); candidate_i++)
53  {
54  const Candidate *candidate = &((*input)[candidate_i]);
55  if( candidate->et() > globalThreshold )
56  {
57  double phi = candidate->phi();
58  double theta = candidate->theta();
59  double e = candidate->energy();
60  double et = e*sin(theta);
61  sum_ez += e*cos(theta);
62  sum_et += et;
63  sum_ex += et*cos(phi);
64  sum_ey += et*sin(phi);
65  }
66  }
67  met->mex = -sum_ex;
68  met->mey = -sum_ey;
69  met->mez = -sum_ez;
70  met->met = sqrt( sum_ex*sum_ex + sum_ey*sum_ey );
71  // cout << "MET = " << met->met << endl;
72  met->sumet = sum_et;
73  met->phi = atan2( -sum_ey, -sum_ex ); // since MET is now a candidate,
74 } // this is no longer needed
75 //------------------------------------------------------------------------
76 
virtual double energy() const =0
energy
virtual double et() const =0
transverse energy
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
virtual void run(edm::Handle< edm::View< reco::Candidate > >, CommonMETData *, double)
Definition: METAlgo.cc:41
virtual double theta() const =0
momentum polar angle
T sqrt(T t)
Definition: SSEVec.h:46
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Structure containing data common to all types of MET.
Definition: CommonMETData.h:22
virtual ~METAlgo()
Definition: METAlgo.cc:26
METAlgo()
Definition: METAlgo.cc:20
virtual double phi() const =0
momentum azimuthal angle
Definition: DDAxes.h:10