Manticore  Version 1.5.3
Physics of Molecular Clouds
Detector.cc
Go to the documentation of this file.
1 
8 #include "mutils/nmath.h"
9 
10 #include "Detector.h"
11 #include "Graybody.h"
12 
13 namespace mu = mutils;
14 
15 namespace manticore {
16 
18 namespace detector {
19 
20 extern "C" {
21 
25 double fWidth(double nu, void *d)
26 { return static_cast<Detector*>(d)->response(nu)/nu; }
27 
28 }
29 
30 } // namescape detector
31 
32 
34 {
35  const char *const fn = "Detector::init";
36 
37  // Effective band limits.
38  freqRange_ = {resp_.front().first, resp_.back().first};
39  MU_DEBUG(fn, "%s response range is [%.2f, %.2f] microns", name(),
41 
42  // Initialize response interpolator.
44 
45  // Compute and set reference bandwidth.
46  calcWidth();
47 }
48 
49 
57 double Detector::response(double nu, gsl_interp_accel *acc) const
58 {
59  if (nu >= freqRange_.first && nu <= freqRange_.second) {
60  return gsl_spline_eval(spline_, nu, acc);
61  } else {
62  // Truncate respone at specified band edges.
63  return 0.0;
64  }
65 }
66 
67 
76 {
77  const char *const fn = "Detector::calcWidth";
80  MU_DEBUG(fn, "%s reference bandwidth: %.2f microns, %.3e Hz",
82 }
83 
84 
90 double Detector::ccf(const Graybody &gray, double T, double Sigma) const
91 { return gray.F(T, Sigma, this)/(width_*gray.Inu(center_, T, Sigma)); }
92 
93 } // namespace manticore
std::pair< double, double > freqRange_
Response edge frequencies (Hz).
Definition: Detector.h:107
double center_
Effective band center (Hz).
Definition: Detector.h:118
MathUtils numerical mathematics library.
virtual void init()
Initializes detector characteristics.
Definition: Detector.cc:33
T front(T... args)
#define MU_DEBUG(src, msg,...)
Log a message at level Log::DEBUG.
Definition: Log.h:68
Package namespace.
Definition: Detector.cc:15
virtual std::string name() const =0
Official detector name (INSTRUMENT-BAND).
virtual double ccf(const Graybody &gray, double T, double Sigma=1.0) const
Color correction factor.
Definition: Detector.cc:90
double F(double T, double Sigma, const Detector *detect, double sr=1.0) const
Single-temperature graybody integrated flux (erg/cm^2/s).
Definition: Graybody.h:198
constexpr double c_light
Speed of light (cm/s).
Definition: manticore.h:24
MathUtils package.
Definition: CommandLine.cc:10
virtual double response(double nu, gsl_interp_accel *acc=nullptr) const
Absolute spectral response function.
Definition: Detector.cc:57
double fWidth(double nu, void *d)
Reference bandwidth integrand.
Definition: Detector.cc:25
Graybody emission model.
Definition: Graybody.h:40
double width_
Effective (flux) bandwidth (Hz).
Definition: Detector.h:119
gsl_spline * spline_
Response interpolator.
Definition: Detector.h:122
T back(T... args)
gsl_spline * initSpline(const tableType &table, bool ln, unsigned stride)
Initializes cubic spline interpolator.
Definition: manticore.cc:81
double Inu(double nu, double T, double Sigma) const
Single-temperature graybody specific intensity (erg/cm^2/s/Hz/sr).
Definition: Graybody.h:92
double kgauss_r(double(*f_r)(double x, void *), void *fdata, double a, double b, double err)
Definition: kgauss.cc:410
tableType resp_
Spectral response table.
Definition: Detector.h:113
void calcWidth()
Calculates effective bandwidth using standard prescription.
Definition: Detector.cc:75