ATLAS Offline Software
Loading...
Searching...
No Matches
eflowCellIntegrator.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
5/*
6 * eflowCellIntegration.h
7 *
8 * Created on: 12.08.2013
9 * Author: tlodd
10 */
11
12#ifndef EFLOWCELLINTEGRATION_H_
13#define EFLOWCELLINTEGRATION_H_
14
15#include <stdexcept>
16#include <iostream>
17
18#include <math.h>
19
20#include "eflowRec/eflowUtil.h"
23
26
33template <class IntegrandType>
35public:
36 eflowRecursiveGaussLegendreIntegrator(IntegrandType* integrand, double error) :
37 m_integrand(integrand), m_error(error), m_depth(1) {
38 assert(m_integrand);
39 }
41
42 virtual double integrate(const eflowRange& range) {
43 /* Try a 5th and 6th order Gauss-Legendre quadrature */
44 double I5 = DoGaussLegendreIntegration(range, 5);
45 double I6 = DoGaussLegendreIntegration(range, 6);
46
47 /* If results agree within m_error threshold, return the mean... */
48 double Imean = (I5 + I6) / 2.0;
49 if (fabs(I5 - I6) / Imean <= m_error) {
50 return Imean;
51 } else {
52 /* ...else recursively part up the integration into n subRanges */
53 return RecurseIntegration(range, 2);
54 }
55 }
56
57 double getError() const {return m_error;}
58
59private:
60
61 double DoGaussLegendreIntegration(const eflowRange& range, int nOrder){
62 /* Perform nth order Gauss-Legendre quadrature, see
63 * http://en.wikipedia.org/wiki/Gaussian_quadrature#Gauss.E2.80.93Legendre_quadrature */
64
65 /* Array offset for legendre weights/roots */
66 int j = nOrder * (nOrder - 1) / 2;
67 const double* roots = legendreRoots + j;
68 const double* weights = legendreWeights + j;
69
70 double rangeCenter = range.getCenter();
71 double rangeHalfWidth = range.getWidth()/2;
72
73 double I = 0.0;
74 double x;
75 if(m_integrand){
76 for (int i = 0; i < nOrder; i++) {
77 x = rangeCenter + roots[i] * rangeHalfWidth;
78 I += weights[i] * m_integrand->evaluate(x);
79 }
80 }
81 else{
82 std::cerr << "eflowCellIntergrator::DoGaussLegendreIntegration ERROR : Invalid pointer to IntegrandType and cannot perform integration" << std::endl;
83 }
84 if (I < 0. || rangeHalfWidth < 0.) std::cerr << "eflowCellIntergrator::DoGaussLegendreIntegration WARNING: I = " << I << "\trange = " << range.print() << std::endl;
85
86 return I * rangeHalfWidth;
87 }
88
89 double RecurseIntegration(const eflowRange& range, int nSubRanges){
90 m_depth++;
91 /* Part up integration range in n subRanges */
92 double subRangeWidth = range.getWidth() / nSubRanges;
93 eflowRange subRange(range.getMin(), range.getMin() + subRangeWidth);
94
95 /* Loop over subranges */
96 double Integral(0.0);
97 for (int i = 0; i < nSubRanges; i++) {
98 Integral += integrate(subRange);
99 subRange.shift(subRangeWidth);
100 }
101 m_depth--;
102
103 return Integral;
104 }
105
106 IntegrandType* m_integrand;
107 double m_error;
109};
110
120template <Exp_t expType> class eflowCellIntegrand {
121public:
122 eflowCellIntegrand(double sigma): m_lookupExp(eflowLookupExp::getInstance()),
123 m_oneOverTwoSigmaSq(1 / (2 * sigma * sigma)), m_norm(m_oneOverTwoSigmaSq / M_PI), m_etaSq(NAN) { }
125
126 inline void setEtaSq(double xSq) { m_etaSq = xSq; }
127
128 inline double evaluateStdExp(double rSq) const { return m_norm * exp(-rSq * m_oneOverTwoSigmaSq); }
129 inline double evaluateLookupExp(double rSq) const { return m_lookupExp->evaluate(rSq * m_oneOverTwoSigmaSq)*m_norm; }
130
132 inline double evaluate(double y);
133
134private:
137 double m_norm;
138 double m_etaSq;
139};
140template<> inline double eflowCellIntegrand<stdExp>::evaluate(double phi) { return evaluateStdExp(m_etaSq+phi*phi); }
141template<> inline double eflowCellIntegrand<lookupExp>::evaluate(double phi) { return evaluateLookupExp(m_etaSq+phi*phi); }
142
154template <int expType> class eflowCellIntegrator {
155public:
156 eflowCellIntegrator(double stdDev, double error) : m_integrand2D(std::make_unique<eflowCellIntegrand<(Exp_t)expType> >(stdDev)), m_outerIntegrator(this, error), m_innerIntegrator(m_integrand2D.get(), error) { }
157 eflowCellIntegrator(const eflowCellIntegrator& original) : m_integrand2D(std::make_unique<eflowCellIntegrand<(Exp_t)expType> >(*(original.m_integrand2D.get()))), m_outerIntegrator(this, original.m_outerIntegrator.getError()), m_innerIntegrator( m_integrand2D.get(), original.m_innerIntegrator.getError()),
158 m_rangePhi (original.m_rangePhi)
159 {
160 }
162 m_integrand2D(std::make_unique<eflowCellIntegrand<(Exp_t)expType> >(*(original.m_integrand2D.get())));
165 m_rangePhi = original.m_rangePhi;
166 }
168
170 inline double integrate(const eflowRange& etaRange, const eflowRange& phiRange) {
171 /* Store the phi range for the inner integration */
172 m_rangePhi = phiRange;
173 /* Invoke the outer integration */
174 return m_outerIntegrator.integrate(etaRange);
175 }
176
178 inline double evaluate(double eta) {
179 /* Set the eta (square) value for the inner integration */
180 m_integrand2D->setEtaSq(eta*eta);
181 /* Perform the inner (i.e. phi) integration (will invoke the evaluate() method on m_innerIntegrator) */
182 return m_innerIntegrator.integrate(m_rangePhi);
183 }
184
185private:
186 std::unique_ptr<eflowCellIntegrand<(Exp_t)expType> > m_integrand2D;
190};
191
192
193#endif /* EFLOWCELLINTEGRATION_H_ */
#define M_PI
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
constexpr double legendreRoots[21]
Stores roots of 5th and 6th order Legendre polynomials.
constexpr double legendreWeights[21]
Stores weights of 5th and 6th order Legendre polynomials.
#define I(x, y, z)
Definition MD5.cxx:116
#define y
#define x
Class to represent the 2D Gaussian integrand.
void setEtaSq(double xSq)
double evaluateLookupExp(double rSq) const
eflowCellIntegrand(double sigma)
const eflowLookupExp * m_lookupExp
double evaluateStdExp(double rSq) const
double evaluate(double y)
The evaluate method for the integration.
double integrate(const eflowRange &etaRange, const eflowRange &phiRange)
Main method, which starts the integration.
eflowRecursiveGaussLegendreIntegrator< eflowCellIntegrand<(Exp_t) expType > > m_innerIntegrator
double evaluate(double eta)
Evaluate method for the outer (i.e.
eflowCellIntegrator & operator=(const eflowCellIntegrator &original)
eflowCellIntegrator(double stdDev, double error)
eflowCellIntegrator(const eflowCellIntegrator &original)
std::unique_ptr< eflowCellIntegrand<(Exp_t) expType > > m_integrand2D
eflowRecursiveGaussLegendreIntegrator< eflowCellIntegrator< expType > > m_outerIntegrator
Lookup-table based exponential function to save CPU time, which is used by eflowCellIntegrator.
void shift(double shift)
Definition eflowUtil.h:116
Class to perform a generic recursive Gauss-Legendre Integration, see http://en.wikipedia....
double DoGaussLegendreIntegration(const eflowRange &range, int nOrder)
virtual double integrate(const eflowRange &range)
double RecurseIntegration(const eflowRange &range, int nSubRanges)
eflowRecursiveGaussLegendreIntegrator(IntegrandType *integrand, double error)
Exp_t
Enum used as template argument to switch between std::exp and the lookup-table based version.
@ lookupExp
eflowRangeBase< double > eflowRange
Definition eflowUtil.h:139
T * get(TKey *tobj)
get a TObject* from a TKey* (why can't a TObject be a TKey?)
Definition hcg.cxx:130
STL namespace.