CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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

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

Private Member Functions

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

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 45 of file IntegralOverPhiFunction.cc.

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

Referenced by Clone().

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

Definition at line 64 of file IntegralOverPhiFunction.cc.

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

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

Member Function Documentation

virtual ROOT::Math::IGenFunction* IntegralOverPhiFunction::Clone ( ) const
inlinevirtual

Definition at line 40 of file IntegralOverPhiFunction.h.

References IntegralOverPhiFunction().

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

Definition at line 222 of file IntegralOverPhiFunction.cc.

223 {
224 //--- virtual function inherited from ROOT::Math::ParamFunction base class;
225 // not implemented, because not neccessary, but needs to be defined to make code compile...
226  edm::LogWarning("") << "Function not implemented yet !" << std::endl;
227 
228  return 0.;
229 }
double IntegralOverPhiFunction::DoEval ( double  x) const
private

Definition at line 118 of file IntegralOverPhiFunction.cc.

References Abs(), alpha_, ecal_dqm_sourceclient-live_cfg::cerr, checkSolutions(), debugLevel_, geometryDiff::epsilon, i, j, Max(), Min(), normalizedPhi(), numSolutionMax1_, numSolutionMax2_, numSolutionMax3_, numSolutionMax4_, numSolutionMin1_, numSolutionMin2_, numSolutionMin3_, numSolutionMin4_, phi, phi0_, Pi, alignCSCRings::r, alignCSCRings::s, lumiQTWidget::t, and theta0_.

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

Definition at line 107 of file IntegralOverPhiFunction.cc.

References x().

108 {
109  theta0_ = param[0];
110  phi0_ = param[1];
111  alpha_ = param[2];
112 
113  return DoEval(x);
114 }
double DoEval(double x) const
double IntegralOverPhiFunction::DoParameterDerivative ( double  ,
const double *  ,
unsigned  int 
) const
privatevirtual

Definition at line 231 of file IntegralOverPhiFunction.cc.

232 {
233 //--- virtual function inherited from ROOT::Math::ParamFunction base class;
234 // not implemented, because not neccessary, but needs to be defined to make code compile...
235  edm::LogWarning("") << "Function not implemented yet !" << std::endl;
236 
237  return 0.;
238 }
void IntegralOverPhiFunction::DoParameterGradient ( double  x,
double *  paramGradient 
) const
private

Definition at line 240 of file IntegralOverPhiFunction.cc.

241 {
242 //--- virtual function inherited from ROOT::Math::ParamFunction base class;
243 // not implemented, because not neccessary, but needs to be defined to make code compile...
244  edm::LogWarning("") << "Function not implemented yet !" << std::endl;
245 }
void IntegralOverPhiFunction::SetParameterAlpha ( double  alpha)

Definition at line 95 of file IntegralOverPhiFunction.cc.

References alpha, and alpha_.

Referenced by IntegrandThetaFunction::DoEval().

96 {
97  alpha_ = alpha;
98 }
float alpha
Definition: AMPTWrapper.h:95
void IntegralOverPhiFunction::SetParameterPhi0 ( double  phi0)

Definition at line 90 of file IntegralOverPhiFunction.cc.

References normalizedPhi(), and phi0_.

Referenced by IntegrandThetaFunction::DoEval().

91 {
92  phi0_ = normalizedPhi(phi0); // map azimuth angle into interval [-pi,+pi]
93 }
double normalizedPhi(double phi)
Definition: normalizedPhi.cc:5
void IntegralOverPhiFunction::SetParameters ( double *  param)
private

Definition at line 100 of file IntegralOverPhiFunction.cc.

References alpha_, phi0_, and theta0_.

101 {
102  theta0_ = param[0];
103  phi0_ = param[1];
104  alpha_ = param[2];
105 }
void IntegralOverPhiFunction::SetParameterTheta0 ( double  theta0)

Definition at line 85 of file IntegralOverPhiFunction.cc.

References theta0_.

Referenced by IntegrandThetaFunction::DoEval().

86 {
87  theta0_ = theta0;
88 }

Member Data Documentation

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

Definition at line 66 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