Manticore  Version 1.0
Physics of Molecular Clouds
Detector.h
Go to the documentation of this file.
1 #ifndef DETECTOR_H
2 #define DETECTOR_H
3 
4 #include <gsl/gsl_integration.h>
5 #include <gsl/gsl_spline.h>
6 
7 #include <mutils/util.h>
8 
9 namespace manticore {
10 
11 // Forward reference.
12 class Graybody;
13 
14 
50 class Detector {
51 public:
53  Detector(bool extend = true) noexcept : extend_(extend) { }
54 
57  if (integ_) { gsl_integration_workspace_free(integ_); }
58  if (spline_) { gsl_spline_free(spline_); }
59  }
60 
62  virtual std::string name() const = 0;
63 
65  virtual void init();
66 
68  virtual double freqBand() const noexcept { return center_; }
69 
74  virtual double freqWidth() const noexcept { return width_; }
75 
77  virtual std::pair<double, double> freqRange() const noexcept
78  { return freqRange_; }
79 
81  bool isExtended() const noexcept { return extend_; }
82 
85  virtual double response(double nu) const;
86 
91  virtual double ccf(const Graybody &gray, double T, double Sigma = 1.0) const;
92 
93 protected:
95  static constexpr size_t subLimit = 512;
96 
98  static constexpr int intKey = GSL_INTEG_GAUSS41;
99 
101  using tableType = std::vector<std::pair<double, double>>;
102 
104  void calcWidth();
105 
107  std::pair<double,double> freqRange_;
108 
111 
113  bool extend_;
114 
115  double center_ = 0.0,
116  width_ = 0.0;
117 
119  gsl_integration_workspace *integ_ = nullptr;
120 
122  gsl_spline *spline_ = nullptr;
123 };
124 
125 } // namespace manticore
126 
127 #endif // DETECTOR_H
std::pair< double, double > freqRange_
Response edge frequencies (Hz).
Definition: Detector.h:107
~Detector()
Destructor.
Definition: Detector.h:56
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
virtual double freqBand() const noexcept
Reference band center frequency (Hz).
Definition: Detector.h:68
virtual double freqWidth() const noexcept
Effective flux-conversion bandwidth (Hz).
Definition: Detector.h:74
Package namespace.
Definition: Detector.cc:13
std::vector< std::pair< double, double > > tableType
Detector response table type.
Definition: Detector.h:101
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
The MathUtils miscellaneous utilities library.
virtual std::pair< double, double > freqRange() const noexcept
Frequency range containing significant response.
Definition: Detector.h:77
Detector base class.
Definition: Detector.h:50
bool extend_
Default to extended (else point) source observation.
Definition: Detector.h:113
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
bool isExtended() const noexcept
Whether extended source detection is enabled.
Definition: Detector.h:81
Detector(bool extend=true) noexcept
Default constructor.
Definition: Detector.h:53
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