CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
IntegralOverPhiFunction Class Reference

#include <PhysicsTools/IsolationUtils/src/IntegralOverPhiFunction.cc>

Inheritance diagram for IntegralOverPhiFunction:

Public Member Functions

ROOT::Math::IGenFunction * Clone () const override
 
 IntegralOverPhiFunction ()
 
void SetParameterAlpha (double alpha)
 
void SetParameterPhi0 (double phi0)
 
void SetParameterTheta0 (double theta0)
 
 ~IntegralOverPhiFunction () override
 

Private Member Functions

double DoDerivative (double x) const
 
double DoEval (double x) const override
 
double DoEvalPar (double x, const double *param) const override
 
double DoParameterDerivative (double, const double *, unsigned int) const override
 
void DoParameterGradient (double x, double *paramGradient) const
 
void SetParameters (double const *param) override
 

Private Attributes

double alpha_
 
unsigned int numSolutionMax1_
 
unsigned int numSolutionMax2_
 
unsigned int numSolutionMax3_
 
unsigned int numSolutionMax4_
 
unsigned int numSolutionMin1_
 
unsigned int numSolutionMin2_
 
unsigned int numSolutionMin3_
 
unsigned int numSolutionMin4_
 
double phi0_
 
double theta0_
 

Static Private Attributes

static const unsigned int debugLevel_ = 0
 

Detailed Description

Description: auxialiary class for fixed area isolation cone computation (this class performs the integration over the azimuthal angle)

Implementation: imported into CMSSW on 05/18/2007

Definition at line 30 of file IntegralOverPhiFunction.h.

Constructor & Destructor Documentation

IntegralOverPhiFunction::IntegralOverPhiFunction ( )

Definition at line 49 of file IntegralOverPhiFunction.cc.

References alpha_, numSolutionMax1_, numSolutionMax2_, numSolutionMax3_, numSolutionMax4_, numSolutionMin1_, numSolutionMin2_, numSolutionMin3_, numSolutionMin4_, phi0_, and theta0_.

Referenced by Clone().

50  : ROOT::Math::ParamFunction<ROOT::Math::IParametricGradFunctionOneDim>(3) {
51  theta0_ = 0.;
52  phi0_ = 0.;
53  alpha_ = 0.;
54 
55  // !!! ONLY FOR TESTING
56  numSolutionMin1_ = 0;
57  numSolutionMax1_ = 0;
58  numSolutionMin2_ = 0;
59  numSolutionMax2_ = 0;
60  numSolutionMin3_ = 0;
61  numSolutionMax3_ = 0;
62  numSolutionMin4_ = 0;
63  numSolutionMax4_ = 0;
64  // FOR TESTING ONLY !!!
65 }
IntegralOverPhiFunction::~IntegralOverPhiFunction ( )
override

Definition at line 67 of file IntegralOverPhiFunction.cc.

References debugLevel_, numSolutionMax1_, numSolutionMax2_, numSolutionMax3_, numSolutionMax4_, numSolutionMin1_, numSolutionMin2_, numSolutionMin3_, and numSolutionMin4_.

67  {
68  // !!! ONLY FOR TESTING
69  if (debugLevel_ > 0) {
70  edm::LogVerbatim("") << "<IntegralOverPhiFunction::~IntegralOverPhiFunction>:" << std::endl
71  << " numSolutionMin1 = " << numSolutionMin1_ << std::endl
72  << " numSolutionMax1 = " << numSolutionMax1_ << std::endl
73  << " numSolutionMin2 = " << numSolutionMin2_ << std::endl
74  << " numSolutionMax2 = " << numSolutionMax2_ << std::endl
75  << " numSolutionMin3 = " << numSolutionMin3_ << std::endl
76  << " numSolutionMax3 = " << numSolutionMax3_ << std::endl
77  << " numSolutionMin4 = " << numSolutionMin4_ << std::endl
78  << " numSolutionMax4 = " << numSolutionMax4_ << std::endl;
79  }
80  // FOR TESTING ONLY !!!
81 }
Log< level::Info, true > LogVerbatim
static const unsigned int debugLevel_

Member Function Documentation

ROOT::Math::IGenFunction* IntegralOverPhiFunction::Clone ( ) const
inlineoverride

Definition at line 39 of file IntegralOverPhiFunction.h.

References IntegralOverPhiFunction().

39 { return new IntegralOverPhiFunction(*this); }
double IntegralOverPhiFunction::DoDerivative ( double  x) const
private

Definition at line 230 of file IntegralOverPhiFunction.cc.

230  {
231  //--- virtual function inherited from ROOT::Math::ParamFunction base class;
232  // not implemented, because not neccessary, but needs to be defined to make code compile...
233  edm::LogWarning("") << "Function not implemented yet !" << std::endl;
234 
235  return 0.;
236 }
Log< level::Warning, false > LogWarning
double IntegralOverPhiFunction::DoEval ( double  x) const
overrideprivate

Definition at line 111 of file IntegralOverPhiFunction.cc.

References alpha_, EcnaPython_AdcPeg12_S1_10_R170298_1_0_150_Dee0::cerr, checkSolutions(), debugLevel_, geometryDiff::epsilon, mps_fire::i, dqmiolumiharvest::j, siStripFEDMonitor_P5_cff::Max, siStripFEDMonitor_P5_cff::Min, normalizedPhi(), numSolutionMax1_, numSolutionMax2_, numSolutionMax3_, numSolutionMax4_, numSolutionMin1_, numSolutionMin2_, numSolutionMin3_, numSolutionMin4_, phi, phi0_, L1TMuonDQMOffline_cfi::phiMax, L1TMuonDQMOffline_cfi::phiMin, Pi, alignCSCRings::r, alignCSCRings::s, submitPVValidationJobs::t, and theta0_.

111  {
112  //--- return zero if theta either close to zero or close to Pi
113  // (phi not well-defined in that case;
114  // numerical expressions might become "NaN"s)
115  const double epsilon = 1.e-3;
116  if (x < epsilon || x > (TMath::Pi() - epsilon))
117  return 0.;
118 
119  //--- calculate trigonometric expressions
120  // (dependend on angle theta0;
121  // polar angle of cone axis);
122  // avoid theta0 exactly zero or exactly Pi
123  // (numerical expressions become "NaN"s)
124  double theta0 = theta0_;
125  if (theta0 < epsilon)
126  theta0 = epsilon;
127  if (theta0 > (TMath::Pi() - epsilon))
128  theta0 = (TMath::Pi() - epsilon);
129  double cosTheta0 = TMath::Cos(theta0);
130  double cos2Theta0 = TMath::Cos(2 * theta0);
131  double sinTheta0 = TMath::Sin(theta0);
132  double cotTheta0 = 1. / TMath::Tan(theta0);
133  double cscTheta0 = 1. / TMath::Sin(theta0);
134  // (dependend on angle phi0;
135  // azimuth angle of cone axis)
136  double cosPhi0 = TMath::Cos(phi0_);
137  double tanPhi0 = TMath::Tan(phi0_);
138  // (dependend on angle theta;
139  // polar angle of point within cone)
140  double cosTheta = TMath::Cos(x);
141  double cos2Theta = TMath::Cos(2 * x);
142  double sinTheta = TMath::Sin(x);
143  double cotTheta = 1. / TMath::Tan(x);
144  double cscTheta = 1. / TMath::Sin(x);
145  // (dependent on angle alpha;
146  // opening angle of cone, measured from cone axis)
147  double cosAlpha = TMath::Cos(alpha_);
148 
149  double s = -cosPhi0 * cosPhi0 *
150  (2 * cosAlpha * cosAlpha + cos2Theta - 4 * cosAlpha * cosTheta * cosTheta0 + cos2Theta0) * sinTheta *
151  sinTheta * sinTheta0 * sinTheta0;
152 
153  //--- return either zero or 2 Pi
154  // in case argument of square-root is zero or negative
155  // (negative values solely arise from rounding errors);
156  // check whether to return zero or 2 Pi:
157  // o return zero
158  // if |theta- theta0| > alpha,
159  // (theta outside cone of opening angle alpha, so vanishing integral over phi)
160  // o return 2 Pi
161  // if |theta- theta0| < alpha
162  // (theta within cone of opening angle alpha;
163  // actually theta0 < alpha in forward/backward direction,
164  // so that integral includes all phi values, hence yielding 2 Pi)
165  if (s <= 0) {
166  if (TMath::Abs(x - theta0) > alpha_)
167  return 0;
168  if (TMath::Abs(x - theta0) <= alpha_)
169  return 2 * TMath::Pi();
170  std::cerr << "Error in <IntegralOverPhiFunction::operator()>: failed to compute return value !" << std::endl;
171  }
172 
173  double r = (1. / TMath::Sqrt(2.)) * (cscTheta * cscTheta * cscTheta0 * cscTheta0 * TMath::Sqrt(s) * tanPhi0);
174  double t = cosPhi0 * (-cotTheta * cotTheta0 + cosAlpha * cscTheta * cscTheta0);
175 
176  double phi[4];
177  phi[0] = -TMath::ACos(t - r);
178  phi[1] = TMath::ACos(t - r);
179  phi[2] = -TMath::ACos(t + r);
180  phi[3] = TMath::ACos(t + r);
181 
182  if (debugLevel_ > 0) {
183  edm::LogVerbatim("") << "phi0 = " << phi0_ << std::endl
184  << "phi[0] = " << phi[0] << " (phi[0] - phi0 = " << (phi[0] - phi0_) * 180 / TMath::Pi() << ")"
185  << std::endl
186  << "phi[1] = " << phi[1] << " (phi[1] - phi0 = " << (phi[1] - phi0_) * 180 / TMath::Pi() << ")"
187  << std::endl
188  << "phi[2] = " << phi[2] << " (phi[2] - phi0 = " << (phi[2] - phi0_) * 180 / TMath::Pi() << ")"
189  << std::endl
190  << "phi[3] = " << phi[3] << " (phi[3] - phi0 = " << (phi[3] - phi0_) * 180 / TMath::Pi() << ")"
191  << std::endl;
192  }
193 
194  double phiMin = 0.;
195  double phiMax = 0.;
196  for (unsigned int i = 0; i < 4; ++i) {
197  for (unsigned int j = (i + 1); j < 4; ++j) {
198  //--- search for the two solutions for phi
199  // that have an equal distance to phi0
200  // and are on either side
201  double dPhi_i = phi[i] - phi0_;
202  double dPhi_j = phi0_ - phi[j];
203  if (TMath::Abs(normalizedPhi(dPhi_i - dPhi_j)) <
204  epsilon) { // map difference in azimuth angle into interval [-pi,+pi] and require it to be negligible
205  //if ( TMath::Abs((phi[i] - phi0_) - (phi0_ - phi[j])) < epsilon ) { // map difference in azimuth angle into interval [-pi,+pi] and require it to be negligible
206  phiMin = TMath::Min(phi[i], phi[j]);
207  phiMax = TMath::Max(phi[i], phi[j]);
208 
209  // !!! ONLY FOR TESTING
210  if (phi[i] == phiMin)
212  if (phi[i] == phiMax)
214  if (phi[j] == phiMin)
216  if (phi[j] == phiMax)
218  // FOR TESTING ONLY !!!
219  }
220  }
221  }
222 
223  if (phiMin == 0 && phiMax == 0) {
224  edm::LogError("") << "failed to compute Return Value !" << std::endl;
225  }
226 
227  return TMath::Abs(normalizedPhi(phi0_ - phiMin)) + TMath::Abs(normalizedPhi(phiMax - phi0_));
228 }
const double Pi
Log< level::Info, true > LogVerbatim
constexpr T normalizedPhi(T phi)
Definition: normalizedPhi.h:8
Log< level::Error, false > LogError
void checkSolutions(unsigned int i, unsigned int &numSolution1, unsigned int &numSolution2, unsigned int &numSolution3, unsigned int &numSolution4)
static const unsigned int debugLevel_
double IntegralOverPhiFunction::DoEvalPar ( double  x,
const double *  param 
) const
overrideprivate

Definition at line 101 of file IntegralOverPhiFunction.cc.

References x.

103 {
104  theta0_ = param[0];
105  phi0_ = param[1];
106  alpha_ = param[2];
107 
108  return DoEval(x);
109 }
double DoEval(double x) const override
double IntegralOverPhiFunction::DoParameterDerivative ( double  ,
const double *  ,
unsigned  int 
) const
overrideprivate

Definition at line 238 of file IntegralOverPhiFunction.cc.

238  {
239  //--- virtual function inherited from ROOT::Math::ParamFunction base class;
240  // not implemented, because not neccessary, but needs to be defined to make code compile...
241  edm::LogWarning("") << "Function not implemented yet !" << std::endl;
242 
243  return 0.;
244 }
Log< level::Warning, false > LogWarning
void IntegralOverPhiFunction::DoParameterGradient ( double  x,
double *  paramGradient 
) const
private

Definition at line 246 of file IntegralOverPhiFunction.cc.

246  {
247  //--- virtual function inherited from ROOT::Math::ParamFunction base class;
248  // not implemented, because not neccessary, but needs to be defined to make code compile...
249  edm::LogWarning("") << "Function not implemented yet !" << std::endl;
250 }
Log< level::Warning, false > LogWarning
void IntegralOverPhiFunction::SetParameterAlpha ( double  alpha)

Definition at line 93 of file IntegralOverPhiFunction.cc.

References alpha, and alpha_.

Referenced by IntegrandThetaFunction::DoEval().

93 { alpha_ = alpha; }
float alpha
Definition: AMPTWrapper.h:105
void IntegralOverPhiFunction::SetParameterPhi0 ( double  phi0)

Definition at line 89 of file IntegralOverPhiFunction.cc.

References normalizedPhi(), and phi0_.

Referenced by IntegrandThetaFunction::DoEval().

89  {
90  phi0_ = normalizedPhi(phi0); // map azimuth angle into interval [-pi,+pi]
91 }
constexpr T normalizedPhi(T phi)
Definition: normalizedPhi.h:8
void IntegralOverPhiFunction::SetParameters ( double const *  param)
overrideprivate

Definition at line 95 of file IntegralOverPhiFunction.cc.

References alpha_, phi0_, and theta0_.

95  {
96  theta0_ = param[0];
97  phi0_ = param[1];
98  alpha_ = param[2];
99 }
void IntegralOverPhiFunction::SetParameterTheta0 ( double  theta0)

Definition at line 87 of file IntegralOverPhiFunction.cc.

References theta0_.

Referenced by IntegrandThetaFunction::DoEval().

87 { theta0_ = theta0; }

Member Data Documentation

double IntegralOverPhiFunction::alpha_
mutableprivate
const unsigned int IntegralOverPhiFunction::debugLevel_ = 0
staticprivate

Definition at line 65 of file IntegralOverPhiFunction.h.

Referenced by DoEval(), and ~IntegralOverPhiFunction().

unsigned int IntegralOverPhiFunction::numSolutionMax1_
mutableprivate
unsigned int IntegralOverPhiFunction::numSolutionMax2_
mutableprivate
unsigned int IntegralOverPhiFunction::numSolutionMax3_
mutableprivate
unsigned int IntegralOverPhiFunction::numSolutionMax4_
mutableprivate
unsigned int IntegralOverPhiFunction::numSolutionMin1_
mutableprivate
unsigned int IntegralOverPhiFunction::numSolutionMin2_
mutableprivate
unsigned int IntegralOverPhiFunction::numSolutionMin3_
mutableprivate
unsigned int IntegralOverPhiFunction::numSolutionMin4_
mutableprivate
double IntegralOverPhiFunction::phi0_
mutableprivate
double IntegralOverPhiFunction::theta0_
mutableprivate