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
31
void
Detector::init
()
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.
44
spline_
=
initSpline
(
resp_
);
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
76
void
Detector::calcWidth
()
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"
,
87
name
(),
width_
*1e4*
c_light
/(
center_
*
center_
),
width_
);
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
manticore::Detector::freqRange_
std::pair< double, double > freqRange_
Response edge frequencies (Hz).
Definition:
Detector.h:107
manticore::Detector::center_
double center_
Effective band center (Hz).
Definition:
Detector.h:115
manticore::Detector::response
virtual double response(double nu) const
Absolute spectral response function.
Definition:
Detector.cc:58
manticore::Detector::init
virtual void init()
Initializes detector characteristics.
Definition:
Detector.cc:31
Graybody.h
manticore::Detector::integ_
gsl_integration_workspace * integ_
Integration workspace.
Definition:
Detector.h:119
manticore::initSpline
gsl_spline * initSpline(const tableType &table, bool ln)
Initializes cubic spline interpolator.
Definition:
manticore.cc:74
MU_DEBUG
#define MU_DEBUG(src, msg,...)
Log a message at level Log::DEBUG.
Definition:
Log.h:68
manticore
Package namespace.
Definition:
Detector.cc:13
manticore::Detector::name
virtual std::string name() const =0
Official detector name.
manticore::Detector::ccf
virtual double ccf(const Graybody &gray, double T, double Sigma=1.0) const
Color correction factor.
Definition:
Detector.cc:96
manticore::Detector::subLimit
static constexpr size_t subLimit
Integral subdivision limit.
Definition:
Detector.h:95
manticore::c_light
constexpr double c_light
Speed of light (cm/s).
Definition:
manticore.h:20
mutils
MathUtils package.
Definition:
CommandLine.cc:10
manticore::detector::fWidth
double fWidth(double nu, void *d)
Reference bandwidth integrand.
Definition:
Detector.cc:23
manticore::Graybody::F
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
manticore::Detector
Detector base class.
Definition:
Detector.h:50
manticore::Graybody
Graybody emission model.
Definition:
Graybody.h:33
manticore::Detector::width_
double width_
Effective (flux) bandwidth (Hz).
Definition:
Detector.h:116
manticore::Detector::spline_
gsl_spline * spline_
Response interpolator.
Definition:
Detector.h:122
manticore::Graybody::Inu
double Inu(double nu, double T, double Sigma) const
Graybody specific intensity (erg/cm^2/s/Hz/sr).
Definition:
Graybody.h:73
manticore::Detector::resp_
tableType resp_
Spectral response table.
Definition:
Detector.h:110
manticore::Detector::intKey
static constexpr int intKey
Integration rule key.
Definition:
Detector.h:98
Detector.h
manticore::Detector::calcWidth
void calcWidth()
Calculates effective bandwidth using standard prescription.
Definition:
Detector.cc:76
Generated on Fri Feb 9 2018 17:14:46 for Manticore by
1.8.14