Manticore  Version 1.0
Physics of Molecular Clouds
Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | Static Protected Attributes | List of all members
manticore::Detector Class Referenceabstract

Detector base class. More...

#include <Detector.h>

Inheritance diagram for manticore::Detector:
Inheritance graph
[legend]

Public Member Functions

 Detector (bool extend=true) noexcept
 Default constructor. More...
 
 ~Detector ()
 Destructor. More...
 
virtual double ccf (const Graybody &gray, double T, double Sigma=1.0) const
 Color correction factor. More...
 
virtual double freqBand () const noexcept
 Reference band center frequency (Hz). More...
 
virtual std::pair< double, double > freqRange () const noexcept
 Frequency range containing significant response. More...
 
virtual double freqWidth () const noexcept
 Effective flux-conversion bandwidth \(\Delta\nu\) (Hz). More...
 
virtual void init ()
 Initializes detector characteristics. More...
 
bool isExtended () const noexcept
 Whether extended source detection is enabled. More...
 
virtual std::string name () const =0
 Official detector name. More...
 
virtual double response (double nu) const
 Absolute spectral response function. More...
 

Protected Types

using tableType = std::vector< std::pair< double, double > >
 Detector response table type. More...
 

Protected Member Functions

void calcWidth ()
 Calculates effective bandwidth using standard prescription. More...
 

Protected Attributes

double center_ = 0.0
 Effective band center (Hz). More...
 
bool extend_
 Default to extended (else point) source observation. More...
 
std::pair< double, double > freqRange_
 Response edge frequencies (Hz). More...
 
gsl_integration_workspace * integ_ = nullptr
 Integration workspace. More...
 
tableType resp_
 Spectral response table. More...
 
gsl_spline * spline_ = nullptr
 Response interpolator. More...
 
double width_ = 0.0
 Effective (flux) bandwidth (Hz). More...
 

Static Protected Attributes

static constexpr int intKey = GSL_INTEG_GAUSS41
 Integration rule key. More...
 
static constexpr size_t subLimit = 512
 Integral subdivision limit. More...
 

Detailed Description

Detector base class.

Concrete detector instances are characterized by the following:

By default this class will compute the reference bandwidth and color correction factor automatically from the response curve using the equations

\[ {\Delta\lambda\over \lambda_0} = {\Delta\nu\over \nu_0} = \int R(\nu){d\nu\over\nu} \]

and

\[ CCF = {\int F_\nu R(\nu)d\nu\over \Delta\nu F_{\nu_0}}, \]

where \(\lambda_0\) and \(\nu_0\) are the reference wavelength and frequency, respectively, and \(F_\nu\) is the model flux density under consideration. If the color correction prescription is inappropriate, concrete detector classes can override the calculation.

Note
This class assumes all integrals are performed in frequency space; but to match common reference, the nominal band center and width are given as wavelengths. Once initialized, an instance may be safely shared among multiple threads (but ccf() is thread-safe only if the Graybody instance is thread-local).
Warning
Clients must call init() on instances to populate detector parameters prior to any further usage.

Definition at line 50 of file Detector.h.

Member Typedef Documentation

◆ tableType

using manticore::Detector::tableType = std::vector<std::pair<double, double> >
protected

Detector response table type.

Definition at line 101 of file Detector.h.

Constructor & Destructor Documentation

◆ Detector()

manticore::Detector::Detector ( bool  extend = true)
inlinenoexcept

Default constructor.

Definition at line 53 of file Detector.h.

◆ ~Detector()

manticore::Detector::~Detector ( )
inline

Destructor.

Definition at line 56 of file Detector.h.

References integ_, and spline_.

Member Function Documentation

◆ calcWidth()

void manticore::Detector::calcWidth ( )
protected

Calculates effective bandwidth using standard prescription.

Uses the response() function to compute the reference bandwidth (see detailed description).

Postcondition
The width_ member is set to the result.

Definition at line 76 of file Detector.cc.

References manticore::c_light, center_, freqRange_, manticore::detector::fWidth(), integ_, intKey, MU_DEBUG, name(), subLimit, and width_.

Referenced by init().

◆ ccf()

double manticore::Detector::ccf ( const Graybody gray,
double  T,
double  Sigma = 1.0 
) const
virtual

Color correction factor.

Parameters
[in]grayGraybody model.
[in]TTemperature (K).
[in]SigmaMass surface density (g/cm^2).
Note
You must divide the reported detector flux density by the correction factor to obtain a fair estimate of the model flux density.

Definition at line 96 of file Detector.cc.

References center_, manticore::Graybody::F(), manticore::Graybody::Inu(), and width_.

Referenced by manticore::ccfCheck().

◆ freqBand()

virtual double manticore::Detector::freqBand ( ) const
inlinevirtualnoexcept

Reference band center frequency (Hz).

Definition at line 68 of file Detector.h.

References center_.

◆ freqRange()

virtual std::pair<double, double> manticore::Detector::freqRange ( ) const
inlinevirtualnoexcept

Frequency range containing significant response.

Definition at line 77 of file Detector.h.

References freqRange_.

Referenced by manticore::Graybody::Itot().

◆ freqWidth()

virtual double manticore::Detector::freqWidth ( ) const
inlinevirtualnoexcept

Effective flux-conversion bandwidth \(\Delta\nu\) (Hz).

The band-weighted total flux \(F\) is related to the reported flux density \(f_\nu\) by \(F = \Delta\nu\cdotf_\nu\).

Definition at line 74 of file Detector.h.

References width_.

◆ init()

void manticore::Detector::init ( )
virtual

Initializes detector characteristics.

Reimplemented in manticore::PACS, and manticore::SPIRE.

Definition at line 31 of file Detector.cc.

References manticore::c_light, calcWidth(), freqRange_, manticore::initSpline(), integ_, MU_DEBUG, name(), resp_, spline_, and subLimit.

Referenced by manticore::SPIRE::init(), and manticore::PACS::init().

◆ isExtended()

bool manticore::Detector::isExtended ( ) const
inlinenoexcept

Whether extended source detection is enabled.

Definition at line 81 of file Detector.h.

References extend_.

◆ name()

virtual std::string manticore::Detector::name ( ) const
pure virtual

Official detector name.

Implemented in manticore::PACS, and manticore::SPIRE.

Referenced by calcWidth(), and init().

◆ response()

double manticore::Detector::response ( double  nu) const
virtual

Absolute spectral response function.

Parameters
[in]nuFrequency (Hz).
Note
For extended source detection, this response includes a correction for the frequency-dependent beam size (scaled to unity at freqBand()). The default implementation is based on linear interpolation of a pre-calculated response table.

Definition at line 58 of file Detector.cc.

References freqRange_, and spline_.

Member Data Documentation

◆ center_

double manticore::Detector::center_ = 0.0
protected

Effective band center (Hz).

Definition at line 115 of file Detector.h.

Referenced by calcWidth(), ccf(), freqBand(), manticore::SPIRE::init(), and manticore::PACS::init().

◆ extend_

bool manticore::Detector::extend_
protected

Default to extended (else point) source observation.

Definition at line 113 of file Detector.h.

Referenced by manticore::SPIRE::init(), and isExtended().

◆ freqRange_

std::pair<double,double> manticore::Detector::freqRange_
protected

Response edge frequencies (Hz).

Definition at line 107 of file Detector.h.

Referenced by calcWidth(), freqRange(), init(), and response().

◆ integ_

gsl_integration_workspace* manticore::Detector::integ_ = nullptr
protected

Integration workspace.

Definition at line 119 of file Detector.h.

Referenced by calcWidth(), init(), and ~Detector().

◆ intKey

constexpr int manticore::Detector::intKey = GSL_INTEG_GAUSS41
staticprotected

Integration rule key.

Definition at line 98 of file Detector.h.

Referenced by calcWidth().

◆ resp_

tableType manticore::Detector::resp_
protected

Spectral response table.

Definition at line 110 of file Detector.h.

Referenced by manticore::SPIRE::init(), manticore::PACS::init(), and init().

◆ spline_

gsl_spline* manticore::Detector::spline_ = nullptr
protected

Response interpolator.

Definition at line 122 of file Detector.h.

Referenced by init(), response(), and ~Detector().

◆ subLimit

constexpr size_t manticore::Detector::subLimit = 512
staticprotected

Integral subdivision limit.

Definition at line 95 of file Detector.h.

Referenced by calcWidth(), and init().

◆ width_

double manticore::Detector::width_ = 0.0
protected

Effective (flux) bandwidth (Hz).

Definition at line 116 of file Detector.h.

Referenced by calcWidth(), ccf(), and freqWidth().


The documentation for this class was generated from the following files: