Manticore  Version 1.0
Physics of Molecular Clouds
Detector.cc
Go to the documentation of this file.
1 
8 #include "Detector.h"
9 #include "Graybody.h"
10 
11 namespace mu = mutils;
12 
13 namespace manticore {
14 
16 namespace detector {
17 
18 extern "C" {
19 
23 double fWidth(double nu, void *d)
24 { return static_cast<Detector*>(d)->response(nu)/nu; }
25 
26 }
27 
28 } // namescape detector
29 
30 
32 {
33  const char *const fn = "Detector::init";
34 
35  // Effective band limits.
36  freqRange_ = {resp_.front().first, resp_.back().first};
37  MU_DEBUG(fn, "%s response range is [%.2f, %.2f] microns", name(),
38  1e4*c_light/freqRange_.second, 1e4*c_light/freqRange_.first);
39 
40  // Initialize integrator workspace.
41  integ_ = gsl_integration_workspace_alloc(subLimit);
42 
43  // Initialize response interpolator.
45 
46  // Compute and set reference bandwidth.
47  calcWidth();
48 }
49 
50 
58 double Detector::response(double nu) const
59 {
60  if (nu >= freqRange_.first && nu <= freqRange_.second) {
61  return gsl_spline_eval(spline_, nu, nullptr);
62  } else {
63  // Truncate respone at specified band edges.
64  return 0.0;
65  }
66 }
67 
68 
77 {
78  const char *const fn = "Detector::calcWidth";
79  double result, abserr;
80  gsl_function f = {detector::fWidth, this};
81  gsl_integration_qag(&f, freqRange_.first, freqRange_.second,
82  0.0, 1e-4, subLimit, intKey,
83  integ_, &result, &abserr);
84 
85  width_ = result*center_;
86  MU_DEBUG(fn, "%s reference bandwidth: %.2f microns, %.3e Hz",
88 }
89 
90 
96 double Detector::ccf(const Graybody &gray, double T, double Sigma) const
97 { return gray.F(T, Sigma, this)/(width_*gray.Inu(center_, T, Sigma)); }
98 
99 } // namespace manticore
std::pair< double, double > freqRange_
Response edge frequencies (Hz).
Definition: Detector.h:107
double center_
Effective band center (Hz).
Definition: Detector.h:115
virtual double response(double nu) const
Absolute spectral response function.
Definition: Detector.cc:58
virtual void init()
Initializes detector characteristics.
Definition: Detector.cc:31
gsl_integration_workspace * integ_
Integration workspace.
Definition: Detector.h:119
gsl_spline * initSpline(const tableType &table, bool ln)
Initializes cubic spline interpolator.
Definition: manticore.cc:74
#define MU_DEBUG(src, msg,...)
Log a message at level Log::DEBUG.
Definition: Log.h:68
Package namespace.
Definition: Detector.cc:13
virtual std::string name() const =0
Official detector name.
virtual double ccf(const Graybody &gray, double T, double Sigma=1.0) const
Color correction factor.
Definition: Detector.cc:96
static constexpr size_t subLimit
Integral subdivision limit.
Definition: Detector.h:95
constexpr double c_light
Speed of light (cm/s).
Definition: manticore.h:20
MathUtils package.
Definition: CommandLine.cc:10
double fWidth(double nu, void *d)
Reference bandwidth integrand.
Definition: Detector.cc:23
double F(double T, double Sigma, const Detector *detect=nullptr, double sr=1.0) const
Graybody integrated flux (erg/cm^2/s).
Definition: Graybody.h:97
Detector base class.
Definition: Detector.h:50
Graybody emission model.
Definition: Graybody.h:33
double width_
Effective (flux) bandwidth (Hz).
Definition: Detector.h:116
gsl_spline * spline_
Response interpolator.
Definition: Detector.h:122
double Inu(double nu, double T, double Sigma) const
Graybody specific intensity (erg/cm^2/s/Hz/sr).
Definition: Graybody.h:73
tableType resp_
Spectral response table.
Definition: Detector.h:110
static constexpr int intKey
Integration rule key.
Definition: Detector.h:98
void calcWidth()
Calculates effective bandwidth using standard prescription.
Definition: Detector.cc:76