Manticore  Version 1.0
Physics of Molecular Clouds
Namespaces | Classes | Typedefs | Functions | Variables
manticore Namespace Reference

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...
 

Detailed Description

Package namespace.

Typedef Documentation

◆ mapDataType

using manticore::mapDataType = typedef std::valarray<float>

Output FITS maps data type (cf. mapDataCode).

Definition at line 11 of file manticore.h.

◆ tableType

using manticore::tableType = typedef std::vector<std::pair<double,double> >

Basic table data.

Definition at line 17 of file manticore.h.

Function Documentation

◆ asctime()

std::string manticore::asctime ( )

Current time a la std::asctime().

Definition at line 114 of file solve.cc.

Referenced by solve(), and solveStage1().

◆ ccfCheck()

void manticore::ccfCheck ( )

◆ df_gray()

int manticore::df_gray ( const gsl_vector *  x,
void *  p,
gsl_matrix *  J 
)

Observation fitting function Jacobian.

Parameters
[in]xModel parameters {T, Sigma}.
[in]pFitting data structure.
[out]JOutput Jacobian.
Returns
GSL_SUCCESS.

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().

◆ f_gray()

int manticore::f_gray ( const gsl_vector *  x,
void *  p,
gsl_vector *  f 
)

Observation fitting function.

Parameters
[in]xModel parameters {T, Sigma}.
[in]pFitting data structure.
[out]fOutput functions.
Returns
GSL_SUCCESS.

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().

◆ fdf_gray()

int manticore::fdf_gray ( const gsl_vector *  x,
void *  p,
gsl_vector *  f,
gsl_matrix *  J 
)

Observation fitting function plus Jacobian.

Parameters
[in]xModel parameters {T, Sigma}.
[in]pFitting data structure.
[out]fOutput functions.
[out]JOutput Jacobian.
Returns
GSL_SUCCESS.

Definition at line 109 of file solve.cc.

References df_gray(), and f_gray().

Referenced by solveStage1().

◆ getError()

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.

Parameters
[in]dataInput image map data.
[in]hduInput image map FITS header.
[in]cliCommand line options and arguments.
Returns
Image 1-sigma error map.

Definition at line 30 of file process.cc.

References mutils::CommandLine::getDouble(), mutils::CommandLine::getOption(), MU_EXCEPTION_IF, and MU_INFO.

Referenced by process().

◆ initSpline()

gsl_spline * manticore::initSpline ( const tableType table,
bool  ln 
)

Initializes cubic spline interpolator.

Note
The input table is not required for interpolation following initialization.
Parameters
[in]tableInterpolation table.
[in]lnWhether to take natural logarithm of the table ordinates.
Returns
Initialized GSL spline instance.

Definition at line 74 of file manticore.cc.

Referenced by manticore::Detector::init(), and manticore::Dust::setModel().

◆ printVersion()

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().

◆ process()

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:

  • reading the input FITS images and headers
  • reading or generating input map uncertainties
  • performing input compatibility checks
  • populating the output FITS header and data structures
  • writing the output FITS file to disk
Parameters
[in]cliCommand 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().

◆ setKeys()

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.

Parameters
[in,out]hduOutput FITS header.
[in]hdu0Fiducial input image header.
[in]cliCommand line options and arguments.

Definition at line 195 of file process.cc.

References MU_EXCEPTION, and version.

Referenced by process().

◆ solve()

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:

  • T (gas temperature)
  • dT (1-sigma uncertainty in T)
  • N(H2) (N(H2) gas column density)
  • dN(H2) (1-sigma uncertainty in N)
  • Chi2 (fit chi-squared per degree of freedom)
Parameters
[out]outMapOutput result data.
[in,out]outHDUOutput FITS headers.
[in,out]outFITSOutput FITS file.
[in]bandDataInput band data.
[in]bandErrorInput band (1-sigma) uncertainties.
[in]bandHDUInput band FITS headers.
[in]cliCommand 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().

◆ solveStage1()

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.

Parameters
[out]outMapOutput result data.
[in]bandDataInput band data.
[in]bandErrorInput band (1-sigma) uncertainties.
[in]cliCommand line options and arguments.
[in]gdata0Model parameters.
[in]idThread ID number (0-based).
[in]nThreadsNumber 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().

Variable Documentation

◆ c_light

constexpr double manticore::c_light = 2.99792458e10

◆ h_Planck

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().

◆ k_Boltzman

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().

◆ m_H

constexpr double manticore::m_H = 1.6737237e-24

Mass of hygdrogen atom (g).

Definition at line 29 of file manticore.h.

Referenced by solveStage1().

◆ mapDataCode

constexpr int manticore::mapDataCode = FLOAT_IMG

Output FITS maps data code (cf. mapDataType).

Definition at line 14 of file manticore.h.

Referenced by process(), and solve().

◆ reldate

constexpr const char* const manticore::reldate = "Feb 09 2018"

Release date (if any).

Definition at line 7 of file version.h.

Referenced by printVersion().

◆ version

constexpr const char* const manticore::version = "1.0"

Version string.

Definition at line 6 of file version.h.

Referenced by printVersion(), and setKeys().