Manticore
Version 1.0
Physics of Molecular Clouds
|
Package namespace. More...
Namespaces | |
detector | |
Detector support. | |
dust | |
Dust models. | |
graybody | |
Detector support. | |
Classes | |
class | Detector |
Detector base class. More... | |
class | Dust |
Dust model. More... | |
class | Graybody |
Graybody emission model. More... | |
struct | GrayData |
Observation fitting data structure. More... | |
class | PACS |
PACS detector class. More... | |
class | SPIRE |
SPIRE detector class. More... | |
Typedefs | |
using | mapDataType = std::valarray< float > |
Output FITS maps data type (cf. mapDataCode). More... | |
using | tableType = std::vector< std::pair< double, double > > |
Basic table data. More... | |
Functions | |
std::string | asctime () |
Current time a la std::asctime(). More... | |
void | ccfCheck () |
int | df_gray (const gsl_vector *x, void *p, gsl_matrix *J) |
Observation fitting function Jacobian. More... | |
int | f_gray (const gsl_vector *x, void *p, gsl_vector *f) |
Observation fitting function. More... | |
int | fdf_gray (const gsl_vector *x, void *p, gsl_vector *f, gsl_matrix *J) |
Observation fitting function plus Jacobian. More... | |
mapDataType | getError (const mapDataType &data, CCfits::PHDU &hdu, const mu::CommandLine &cli) |
Generates image error map. More... | |
gsl_spline * | initSpline (const tableType &table, bool ln) |
Initializes cubic spline interpolator. More... | |
void | printVersion () |
Prints version information to the log. More... | |
void | process (const mu::CommandLine &cli) |
Processes multi-wavelength maps into temperature/column density maps. More... | |
void | setKeys (CCfits::ExtHDU &hdu, const CCfits::PHDU &hdu0) |
Sets header keywords for an output map. More... | |
void | solve (std::vector< mapDataType > &outMap, std::vector< CCfits::ExtHDU * > &outHDU, CCfits::FITS *outFITS, const std::vector< mapDataType > &bandData, const std::vector< mapDataType > &bandError, const std::vector< CCfits::PHDU * > &bandHDU, const mu::CommandLine &cli) |
Performs least-squares fitting of band SEDs. More... | |
void | solveStage1 (const std::vector< mapDataType > &bandData, const std::vector< mapDataType > &bandError, const mu::CommandLine &cli, const GrayData &gdata0, unsigned id, unsigned nThreads) |
Stage 1 solver. More... | |
Variables | |
constexpr double | c_light = 2.99792458e10 |
Speed of light (cm/s). More... | |
constexpr double | h_Planck = 6.626070040e-27 |
Planck's constant (erg*s). More... | |
constexpr double | k_Boltzman = 1.38064852e-16 |
Planck's constant (erg/K). More... | |
constexpr double | m_H = 1.6737237e-24 |
Mass of hygdrogen atom (g). More... | |
constexpr int | mapDataCode = FLOAT_IMG |
Output FITS maps data code (cf. mapDataType). More... | |
constexpr const char *const | reldate = "Feb 09 2018" |
Release date (if any). More... | |
constexpr const char *const | version = "1.0" |
Version string. More... | |
Package namespace.
using manticore::mapDataType = typedef std::valarray<float> |
Output FITS maps data type (cf. mapDataCode).
Definition at line 11 of file manticore.h.
using manticore::tableType = typedef std::vector<std::pair<double,double> > |
Basic table data.
Definition at line 17 of file manticore.h.
std::string manticore::asctime | ( | ) |
Current time a la std::asctime().
Definition at line 114 of file solve.cc.
Referenced by solve(), and solveStage1().
void manticore::ccfCheck | ( | ) |
Definition at line 17 of file manticore.cc.
References manticore::Detector::ccf(), manticore::Graybody::dF_dS(), manticore::Graybody::dF_dT(), manticore::Graybody::F(), manticore::PACS::init(), manticore::Graybody::Inu(), manticore::Dust::setModel(), and mutils::Log::stream.
Referenced by main().
int manticore::df_gray | ( | const gsl_vector * | x, |
void * | p, | ||
gsl_matrix * | J | ||
) |
Observation fitting function Jacobian.
[in] | x | Model parameters {T, Sigma}. |
[in] | p | Fitting data structure. |
[out] | J | Output Jacobian. |
Definition at line 83 of file solve.cc.
References manticore::GrayData::band, manticore::GrayData::Ferr, and manticore::GrayData::gray.
Referenced by fdf_gray(), and solveStage1().
int manticore::f_gray | ( | const gsl_vector * | x, |
void * | p, | ||
gsl_vector * | f | ||
) |
Observation fitting function.
[in] | x | Model parameters {T, Sigma}. |
[in] | p | Fitting data structure. |
[out] | f | Output functions. |
Definition at line 59 of file solve.cc.
References manticore::GrayData::band, manticore::GrayData::Ferr, manticore::GrayData::Fobs, and manticore::GrayData::gray.
Referenced by fdf_gray(), and solveStage1().
int manticore::fdf_gray | ( | const gsl_vector * | x, |
void * | p, | ||
gsl_vector * | f, | ||
gsl_matrix * | J | ||
) |
Observation fitting function plus Jacobian.
[in] | x | Model parameters {T, Sigma}. |
[in] | p | Fitting data structure. |
[out] | f | Output functions. |
[out] | J | Output Jacobian. |
Definition at line 109 of file solve.cc.
References df_gray(), and f_gray().
Referenced by solveStage1().
mapDataType manticore::getError | ( | const mapDataType & | data, |
CCfits::PHDU & | hdu, | ||
const mu::CommandLine & | cli | ||
) |
Generates image error map.
If any --error-map
options were supplied and an entry matching the BAND and TELESCOP keywords of the input is found, the error map is read from the corresponding image. Otherwise, the value of --error-pct
(default: 3.0) is multiplied by data to create the error map.
[in] | data | Input image map data. |
[in] | hdu | Input image map FITS header. |
[in] | cli | Command line options and arguments. |
Definition at line 30 of file process.cc.
References mutils::CommandLine::getDouble(), mutils::CommandLine::getOption(), MU_EXCEPTION_IF, and MU_INFO.
Referenced by process().
gsl_spline * manticore::initSpline | ( | const tableType & | table, |
bool | ln | ||
) |
Initializes cubic spline interpolator.
table
is not required for interpolation following initialization.[in] | table | Interpolation table. |
[in] | ln | Whether to take natural logarithm of the table ordinates. |
Definition at line 74 of file manticore.cc.
Referenced by manticore::Detector::init(), and manticore::Dust::setModel().
void manticore::printVersion | ( | ) |
Prints version information to the log.
Definition at line 86 of file manticore.cc.
References reldate, mutils::Log::stream, and version.
Referenced by main().
void manticore::process | ( | const mu::CommandLine & | cli | ) |
Processes multi-wavelength maps into temperature/column density maps.
This routine handles the higher-level management required to generate mass and temperature images from the input continuum data, including:
[in] | cli | Command line options and arguments. |
Definition at line 84 of file process.cc.
References mutils::CommandLine::arguments(), getError(), mutils::CommandLine::getOption(), mutils::CommandLine::getString(), mapDataCode, MU_ERROR_IF, MU_EXCEPTION_IF, MU_INFO, setKeys(), and solve().
Referenced by main().
void manticore::setKeys | ( | CCfits::ExtHDU & | hdu, |
const CCfits::PHDU & | hdu0 | ||
) |
Sets header keywords for an output map.
Combines information from the input FITS headers and the command line arguments to populate keywords in the output FITS file.
[in,out] | hdu | Output FITS header. |
[in] | hdu0 | Fiducial input image header. |
[in] | cli | Command line options and arguments. |
Definition at line 195 of file process.cc.
References MU_EXCEPTION, and version.
Referenced by process().
void manticore::solve | ( | std::vector< mapDataType > & | outMap, |
std::vector< CCfits::ExtHDU * > & | outHDU, | ||
CCfits::FITS * | outFITS, | ||
const std::vector< mapDataType > & | bandData, | ||
const std::vector< mapDataType > & | bandError, | ||
const std::vector< CCfits::PHDU * > & | bandHDU, | ||
const mu::CommandLine & | cli | ||
) |
Performs least-squares fitting of band SEDs.
Given input band continuum fluxes and uncertainties, calculates the best-fit temperature and column density maps (including uncertainties and reduced chi-square information). This routine does not hard-wire a particular number of temperature/density components; an appropriate fitting engine is selected based on the --components NCOMP
option (default: 1). Currently only NCOMP
= 1 is supported.
For NCOMP
= 1 the output FITS file will contain the following named image extensions:
[out] | outMap | Output result data. |
[in,out] | outHDU | Output FITS headers. |
[in,out] | outFITS | Output FITS file. |
[in] | bandData | Input band data. |
[in] | bandError | Input band (1-sigma) uncertainties. |
[in] | bandHDU | Input band FITS headers. |
[in] | cli | Command line options and arguments. |
Definition at line 286 of file solve.cc.
References asctime(), manticore::GrayData::band, mutils::CommandLine::getDouble(), mutils::CommandLine::getInteger(), mutils::CommandLine::getOption(), mutils::CommandLine::getString(), mapDataCode, MU_INFO, and solveStage1().
Referenced by process().
void manticore::solveStage1 | ( | const std::vector< mapDataType > & | bandData, |
const std::vector< mapDataType > & | bandError, | ||
const mu::CommandLine & | cli, | ||
const GrayData & | gdata0, | ||
unsigned | id, | ||
unsigned | nThreads | ||
) |
Stage 1 solver.
Single-plane mass/temperature solution.
[out] | outMap | Output result data. |
[in] | bandData | Input band data. |
[in] | bandError | Input band (1-sigma) uncertainties. |
[in] | cli | Command line options and arguments. |
[in] | gdata0 | Model parameters. |
[in] | id | Thread ID number (0-based). |
[in] | nThreads | Number of threads. |
Definition at line 138 of file solve.cc.
References asctime(), df_gray(), f_gray(), fdf_gray(), mutils::CommandLine::getDouble(), mutils::CommandLine::getInteger(), manticore::GrayData::gray, m_H, MU_DEBUG, MU_INFO_IF, MU_WARN_IF, and manticore::GrayData::outMap.
Referenced by solve().
constexpr double manticore::c_light = 2.99792458e10 |
Speed of light (cm/s).
Definition at line 20 of file manticore.h.
Referenced by manticore::Detector::calcWidth(), manticore::SPIRE::init(), manticore::PACS::init(), manticore::Detector::init(), main(), and manticore::Dust::setModel().
constexpr double manticore::h_Planck = 6.626070040e-27 |
Planck's constant (erg*s).
Definition at line 23 of file manticore.h.
Referenced by manticore::Graybody::Bnu(), manticore::Graybody::dBnu_dT(), and manticore::Graybody::Itot().
constexpr double manticore::k_Boltzman = 1.38064852e-16 |
Planck's constant (erg/K).
Definition at line 26 of file manticore.h.
Referenced by manticore::Graybody::Bnu(), manticore::Graybody::dBnu_dT(), and manticore::Graybody::Itot().
constexpr double manticore::m_H = 1.6737237e-24 |
constexpr int manticore::mapDataCode = FLOAT_IMG |
Output FITS maps data code (cf. mapDataType).
Definition at line 14 of file manticore.h.
constexpr const char* const manticore::reldate = "Feb 09 2018" |
constexpr const char* const manticore::version = "1.0" |
Version string.
Definition at line 6 of file version.h.
Referenced by printVersion(), and setKeys().