ATLAS Offline Software
Public Member Functions | Private Attributes | List of all members
eflowCellIntegrator< expType > Class Template Reference

Class that controls the 2D integration. More...

#include <eflowCellIntegrator.h>

Collaboration diagram for eflowCellIntegrator< expType >:

Public Member Functions

 eflowCellIntegrator (double stdDev, double error)
 
 eflowCellIntegrator (const eflowCellIntegrator &original)
 
eflowCellIntegratoroperator= (const eflowCellIntegrator &original)
 
 ~eflowCellIntegrator ()
 
double integrate (const eflowRange &etaRange, const eflowRange &phiRange)
 Main method, which starts the integration. More...
 
double evaluate (double eta)
 Evaluate method for the outer (i.e. More...
 

Private Attributes

std::unique_ptr< eflowCellIntegrand<(Exp_t) expType > > m_integrand2D
 
eflowRecursiveGaussLegendreIntegrator< eflowCellIntegrator< expType > > m_outerIntegrator
 
eflowRecursiveGaussLegendreIntegrator< eflowCellIntegrand<(Exp_t) expType > > m_innerIntegrator
 
eflowRange m_rangePhi
 

Detailed Description

template<int expType>
class eflowCellIntegrator< expType >

Class that controls the 2D integration.

Holds two eflowRecursiveGaussLegendreIntegrators, one to perform the outer (i.e. eta) and one to perform the inner (i.e. phi) integration. The class itself acts as the integrand of the outer integration. Its evaluate method performs a full integral in phi at the given eta value. (This is how the 2D integration is achieved.) The integrand of the inner (i.e. phi) integration is an eflowCellIntegrand.

To evaluate the actual Gaussian, either std::exp or eflowLookupExp is used, depending on the value of the template parameter (see eflowCellIntegrand). Templates are used instead of polymorphism here due to the better CPU performance.

Definition at line 154 of file eflowCellIntegrator.h.

Constructor & Destructor Documentation

◆ eflowCellIntegrator() [1/2]

template<int expType>
eflowCellIntegrator< expType >::eflowCellIntegrator ( double  stdDev,
double  error 
)
inline

Definition at line 156 of file eflowCellIntegrator.h.

◆ eflowCellIntegrator() [2/2]

template<int expType>
eflowCellIntegrator< expType >::eflowCellIntegrator ( const eflowCellIntegrator< expType > &  original)
inline

Definition at line 157 of file eflowCellIntegrator.h.

158  m_rangePhi (original.m_rangePhi)
159  {
160  }

◆ ~eflowCellIntegrator()

template<int expType>
eflowCellIntegrator< expType >::~eflowCellIntegrator ( )
inline

Definition at line 167 of file eflowCellIntegrator.h.

167 {}

Member Function Documentation

◆ evaluate()

template<int expType>
double eflowCellIntegrator< expType >::evaluate ( double  eta)
inline

Evaluate method for the outer (i.e.

eta) integration (invoked by m_outerIntegrator)

Definition at line 178 of file eflowCellIntegrator.h.

178  {
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) */
183  }

◆ integrate()

template<int expType>
double eflowCellIntegrator< expType >::integrate ( const eflowRange etaRange,
const eflowRange phiRange 
)
inline

Main method, which starts the integration.

Definition at line 170 of file eflowCellIntegrator.h.

170  {
171  /* Store the phi range for the inner integration */
173  /* Invoke the outer integration */
175  }

◆ operator=()

template<int expType>
eflowCellIntegrator& eflowCellIntegrator< expType >::operator= ( const eflowCellIntegrator< expType > &  original)
inline

Definition at line 161 of file eflowCellIntegrator.h.

161  {
162  m_integrand2D(std::make_unique<eflowCellIntegrand<(Exp_t)expType> >(*(original.m_integrand2D.get())));
163  m_outerIntegrator(this, original.m_outerIntegrator.getError());
165  m_rangePhi = original.m_rangePhi;
166  }

Member Data Documentation

◆ m_innerIntegrator

template<int expType>
eflowRecursiveGaussLegendreIntegrator<eflowCellIntegrand<(Exp_t)expType> > eflowCellIntegrator< expType >::m_innerIntegrator
private

Definition at line 188 of file eflowCellIntegrator.h.

◆ m_integrand2D

template<int expType>
std::unique_ptr<eflowCellIntegrand<(Exp_t)expType> > eflowCellIntegrator< expType >::m_integrand2D
private

Definition at line 186 of file eflowCellIntegrator.h.

◆ m_outerIntegrator

template<int expType>
eflowRecursiveGaussLegendreIntegrator<eflowCellIntegrator<expType> > eflowCellIntegrator< expType >::m_outerIntegrator
private

Definition at line 187 of file eflowCellIntegrator.h.

◆ m_rangePhi

template<int expType>
eflowRange eflowCellIntegrator< expType >::m_rangePhi
private

Definition at line 189 of file eflowCellIntegrator.h.


The documentation for this class was generated from the following file:
Exp_t
Exp_t
Enum used as template argument to switch between std::exp and the lookup-table based version.
Definition: eflowCellIntegrator.h:25
make_unique
std::unique_ptr< T > make_unique(Args &&... args)
Definition: SkimmingToolEXOT5.cxx:23
eta
Scalar eta() const
pseudorapidity method
Definition: AmgMatrixBasePlugin.h:83
LArCellBinning.etaRange
etaRange
Filling Eta range.
Definition: LArCellBinning.py:128
eflowRecursiveGaussLegendreIntegrator::integrate
virtual double integrate(const eflowRange &range)
Definition: eflowCellIntegrator.h:42
eflowCellIntegrator::m_innerIntegrator
eflowRecursiveGaussLegendreIntegrator< eflowCellIntegrand<(Exp_t) expType > > m_innerIntegrator
Definition: eflowCellIntegrator.h:188
eflowRecursiveGaussLegendreIntegrator::getError
double getError() const
Definition: eflowCellIntegrator.h:57
eflowCellIntegrator::m_rangePhi
eflowRange m_rangePhi
Definition: eflowCellIntegrator.h:189
eflowCellIntegrator::m_integrand2D
std::unique_ptr< eflowCellIntegrand<(Exp_t) expType > > m_integrand2D
Definition: eflowCellIntegrator.h:186
eflowCellIntegrator::m_outerIntegrator
eflowRecursiveGaussLegendreIntegrator< eflowCellIntegrator< expType > > m_outerIntegrator
Definition: eflowCellIntegrator.h:187
LArCellBinning.phiRange
phiRange
Filling Phi ranges.
Definition: LArCellBinning.py:107
error
Definition: IImpactPoint3dEstimator.h:70
eflowCellIntegrand
Class to represent the 2D Gaussian integrand.
Definition: eflowCellIntegrator.h:120