CMS 3D CMS Logo

IntegrandThetaFunction.cc
Go to the documentation of this file.
2 
3 // -*- C++ -*-
4 //
5 // Package: IntegrandThetaFunction
6 // Class: IntegrandThetaFunction
7 //
16 //
17 // Original Author: Christian Veelken, UC Davis
18 // Created: Thu Nov 2 13:47:40 CST 2006
19 //
20 //
21 
22 // system include files
23 #include <iostream>
24 #include <iomanip>
25 #include <vector>
26 
27 // ROOT include files
28 #include <TMath.h>
29 
30 // CMSSW include files
32 
35 
36 //
37 // constructors and destructor
38 //
39 
41  : ROOT::Math::ParamFunction<ROOT::Math::IParametricGradFunctionOneDim>(3) {
42  theta0_ = 0.;
43  phi0_ = 0.;
44  alpha_ = 0.;
45 
47 }
48 
50  : ROOT::Math::ParamFunction<ROOT::Math::IParametricGradFunctionOneDim>(bluePrint) {
51  theta0_ = bluePrint.theta0_;
52  phi0_ = bluePrint.phi0_;
53  alpha_ = bluePrint.alpha_;
54 
55  fPhi_ = new IntegralOverPhiFunction(*bluePrint.fPhi_);
56 }
57 
59 
60 //
61 // assignment operator
62 //
63 
65  theta0_ = bluePrint.theta0_;
66  phi0_ = bluePrint.phi0_;
67  alpha_ = bluePrint.alpha_;
68 
69  (*fPhi_) = (*bluePrint.fPhi_);
70 
71  return (*this);
72 }
73 
74 //
75 // member functions
76 //
77 
78 void IntegrandThetaFunction::SetParameterTheta0(double theta0) { theta0_ = theta0; }
79 
81  phi0_ = normalizedPhi(phi0); // map azimuth angle into interval [-pi,+pi]
82 }
83 
85 
86 void IntegrandThetaFunction::SetParameters(const double* param) {
87  theta0_ = param[0];
88  phi0_ = param[1];
89  alpha_ = param[2];
90 }
91 
92 double IntegrandThetaFunction::DoEvalPar(double x, const double* param) const //FIXME: constness
93 {
94  theta0_ = param[0];
95  phi0_ = param[1];
96  alpha_ = param[2];
97 
98  return DoEval(x);
99 }
100 
101 double IntegrandThetaFunction::DoEval(double x) const {
102  //--- return zero if theta either close to zero or close to Pi
103  // (numerical expressions might become "NaN"s)
104  const double epsilon = 1.e-3;
105  if (x < epsilon || x > (TMath::Pi() - epsilon))
106  return 0.;
107 
108  //--- calculate trigonometric expressions
109  // (dependend on angle theta;
110  // polar angle of point within cone)
111  double sinTheta = TMath::Sin(x);
112  double cscTheta = 1. / sinTheta;
113 
114  double detJacobi =
115  -cscTheta; // partial derrivative dEta/dTheta (for constant particle density in tau id. isolation cone)
116  //double detJacobi = 1.; // ordinary solid angle (FOR TESTING ONLY)
117 
118  //--- evaluate integral over angle phi
119  // (azimuth angle of point within cone)
123 
124  double integralOverPhi = (*fPhi_)(x);
125 
126  if (debugLevel_ > 0) {
127  edm::LogVerbatim("") << "integralOverPhi = " << integralOverPhi << std::endl
128  << " theta0 = " << theta0_ << std::endl
129  << " phi0 = " << phi0_ << std::endl
130  << " alpha = " << alpha_ << std::endl
131  << " theta = " << x << std::endl
132  << std::endl;
133  }
134 
135  //--- integrand for integration over theta
136  // equals
137  // |dEta/dTheta| * integral over phi * sin(theta)
138  //
139  // (o the factor dEta/dTheta represents the particle density as function of theta,
140  // assuming that the particle density is flat in eta;
141  // o the factor sin(theta) originates from the solid angle surface element
142  // expressed in spherical polar coordinates)
143  //
144  //return TMath::Abs(detJacobi)*integralOverPhi*sinTheta;
145  return TMath::Abs(detJacobi) * integralOverPhi;
146 }
147 
149  //--- virtual function inherited from ROOT::Math::ParamFunction base class;
150  // not implemented, because not neccessary, but needs to be defined to make code compile...
151  edm::LogWarning("") << "Function not implemented yet !" << std::endl;
152 
153  return 0.;
154 }
155 
156 double IntegrandThetaFunction::DoParameterDerivative(double, const double*, unsigned int) const {
157  //--- virtual function inherited from ROOT::Math::ParamFunction base class;
158  // not implemented, because not neccessary, but needs to be defined to make code compile...
159  edm::LogWarning("") << "Function not implemented yet !" << std::endl;
160 
161  return 0.;
162 }
163 
164 void IntegrandThetaFunction::DoParameterGradient(double x, double* paramGradient) const {
165  //--- virtual function inherited from ROOT::Math::ParamFunction base class;
166  // not implemented, because not neccessary, but needs to be defined to make code compile...
167  edm::LogWarning("") << "Function not implemented yet !" << std::endl;
168 }
IntegrandThetaFunction::SetParameters
void SetParameters(double const *param) override
Definition: IntegrandThetaFunction.cc:86
IntegrandThetaFunction::alpha_
double alpha_
Definition: IntegrandThetaFunction.h:58
MessageLogger.h
IntegralOverPhiFunction::SetParameterTheta0
void SetParameterTheta0(double theta0)
Definition: IntegralOverPhiFunction.cc:87
IntegrandThetaFunction::SetParameterAlpha
void SetParameterAlpha(double alpha)
Definition: IntegrandThetaFunction.cc:84
IntegrandThetaFunction::DoEval
double DoEval(double x) const override
Definition: IntegrandThetaFunction.cc:101
IntegralOverPhiFunction::SetParameterPhi0
void SetParameterPhi0(double phi0)
Definition: IntegralOverPhiFunction.cc:89
IntegrandThetaFunction::IntegrandThetaFunction
IntegrandThetaFunction()
Definition: IntegrandThetaFunction.cc:40
IntegrandThetaFunction::SetParameterTheta0
void SetParameterTheta0(double theta0)
Definition: IntegrandThetaFunction.cc:78
IntegrandThetaFunction::SetParameterPhi0
void SetParameterPhi0(double phi0)
Definition: IntegrandThetaFunction.cc:80
IntegrandThetaFunction::~IntegrandThetaFunction
~IntegrandThetaFunction() override
Definition: IntegrandThetaFunction.cc:58
alpha
float alpha
Definition: AMPTWrapper.h:105
DDAxes::x
IntegrandThetaFunction.h
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
geometryDiff.epsilon
int epsilon
Definition: geometryDiff.py:26
Abs
T Abs(T a)
Definition: MathUtil.h:49
normalizedPhi
constexpr T normalizedPhi(T phi)
Definition: normalizedPhi.h:8
IntegrandThetaFunction::DoParameterGradient
void DoParameterGradient(double x, double *paramGradient) const
Definition: IntegrandThetaFunction.cc:164
IntegrandThetaFunction
Definition: IntegrandThetaFunction.h:33
IntegrandThetaFunction::fPhi_
IntegralOverPhiFunction * fPhi_
Definition: IntegrandThetaFunction.h:60
IntegrandThetaFunction::DoParameterDerivative
double DoParameterDerivative(double, const double *, unsigned int) const override
Definition: IntegrandThetaFunction.cc:156
IntegralOverPhiFunction::SetParameterAlpha
void SetParameterAlpha(double alpha)
Definition: IntegralOverPhiFunction.cc:93
normalizedPhi.h
IntegrandThetaFunction::debugLevel_
static const unsigned int debugLevel_
Definition: IntegrandThetaFunction.h:62
IntegrandThetaFunction::operator=
IntegrandThetaFunction & operator=(const IntegrandThetaFunction &bluePrint)
Definition: IntegrandThetaFunction.cc:64
IntegrandThetaFunction::theta0_
double theta0_
Definition: IntegrandThetaFunction.h:56
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
genVertex_cff.x
x
Definition: genVertex_cff.py:12
IntegrandThetaFunction::DoDerivative
double DoDerivative(double x) const
Definition: IntegrandThetaFunction.cc:148
IntegrandThetaFunction::phi0_
double phi0_
Definition: IntegrandThetaFunction.h:57
Pi
const double Pi
Definition: CosmicMuonParameters.h:18
IntegrandThetaFunction::DoEvalPar
double DoEvalPar(double, const double *) const override
Definition: IntegrandThetaFunction.cc:92
IntegralOverPhiFunction
Definition: IntegralOverPhiFunction.h:30
IntegralOverPhiFunction.h
ROOT
Definition: Transform3DPJ.h:35