ModelObject

Overview

The heart of Imfit is the ModelObject class. This holds information about an image model and can compute images using the model and the current set of model parameter values; it can also hold a data image and its associated information, and can compute chi^2 and and other statistics by comparing the model and data images.

This class is instantiated once in each program, and holds (among other things) the following objects and data (”imfit only” denotes data used in imfit or imfit-mcmc, but not needed when used in makeimage)

  • Model image pixel values

  • Vector of FunctionObject instances which define the model

  • Array of current parameter values for the model

  • Image pixel data for PSF(s), if any

  • Convolver objects(s) for PSF convolution

  • imfit only: Data-image pixel values

  • imfit only: other characteristics of the data image (size, gain, read noise, etc.)

  • imfit only: Error image (either supplied by user or internally generated)

  • imfit only: Mask image, if supplied by user

  • imfit only: Internal weight image (from combination of error image and mask image).

Not all of these data members are initialized or used – for example, none of the data-related (”imfit only”) members are used by makeimage, and the error image is not generated or used if a fit uses Poisson-based statistics instead of chi^2.

The main functionality of the class includes:

  • Generating a model image using the vector of FunctionObjects and the current array of parameter values

  • Generating individual-function model images (i.e., the --output-functions option of makeimage)

  • Computing individual-function and total fluxes for current model

  • Computing deviances vector between current model image and data image (for use by Levenberg-Marquardt solver)

  • Computing the fit statistic (chi^2, etc.) from comparison of the current model image and the data image

  • Printing current parameter values in different formats

  • Generating a bootstrap resampling pixel-index vector when bootstrap resampling is being done

API

[Warning: This is currently somewhat incomplete!]

Files: core/model_object.h, core/model_object.cpp

class ModelObject

Main class holding data, model information, and code for generating model images, computing chi^2, etc.

Public Functions

ModelObject()

Constructor.

virtual ~ModelObject()

Destructor.

void SetVerboseLevel(int verbosity)

Set the verbosity level.

void SetDebugLevel(int debuggingLevel)

Set the debugging level (must be 0 [default] or larger).

void SetMaxThreads(int maxThreadNumber)

Specify the maximum number of OpenMP threads to use in computations; also sets maximum number FFTW threads for convolutions.

void SetOMPChunkSize(int chunkSize)

Sets the chunk size for OpenMP.

virtual int AddFunction(FunctionObject *newFunctionObj_ptr)

Adds a FunctionObject subclass to the model.

int SetupPsfInterpolation(int interpolationType = kInterpolator_bicubic)

Specify that PSF interpolation (by PointSource functions) will be used; causes an internal PsfInterpolator object of the appropriate subclass to be allocated and set up with previously supplied PSF data. This should only be called after AddPSFVector has been called to supply the PSF data.

virtual void DefineFunctionSets(vector<int> &functionStartIndices)

Given a vector of indices specifying the start of individual function sets, set up the internal fsetStartFlags array; also compute nParamsTot

virtual void SetZeroPoint(double zeroPointValue)

Set the internal zero-point (for printing flux values)

virtual double CashStatistic(double params[])

Computes and returns the Cash statistic for the current set of model parameters (computes model image and compares it with the data image)

virtual void PrintDescription()

Prints the number of data values (pixels) in the data image.

void GetFunctionNames(vector<string> &functionNames)

Adds the names of image functions in the model (calling each function’s GetShortName method) to the input vector of strings

void GetFunctionLabels(vector<string> &functionLabels)

Adds the names of image functions in the model (calling each function’s GetShortName method) to the input vector of strings

void GetImageOffsets(double params[])

Given an input array (all zeros), this method returns the array with locations corresponding to each function set’s imageOffset_X0 and imageOffset_Y0 values properly set (if there are no image offsets, then the corresponding values will remain zero). Intended for use in bootstrap_errors.cpp.

virtual string GetParamHeader()

Prints all function and parameter names in order all on one line; e.g., for use as header in bootstrap-parameters output file.

virtual int PrintModelParamsToStrings(vector<string> &stringVector, double params[], double errs[], const char *prefix = "", bool printLimits = false)

Prints description of model (e.g., best-fit result) as strings to the input vector of string. Optionally, the lower and upper limits defined in parameterInfo are also printed, OR associated lower and upper error bounds in errs can be printed.

If errs != NULL, then +/- errors are printed as well (only if printLimits is false)

If prefix != NULL, then the specified character (e.g., ‘#’) is prepended to each output line.

If printLimits == true, then lower and upper parameter limits will be printed for each parameter (or else “fixed” for fixed parameters)

virtual string PrintModelParamsHorizontalString(const double params[], const string &separator = "\)

Like PrintModelParamsToString, but prints parameter values all in one line to a string (without parameter names or limits or errors), which is returned. Meant to be used in printing results of bootstrap resampling (imfit) or MCMC chains (imfit-mcmc)

void PrintImage(double *pixelVector, int nColumns, int nRows)

Basic function which prints an image to stdout. Mainly meant to be called by PrintInputImage, PrintModelImage, and PrintWeights.

virtual void PrintInputImage()

Prints the input data image to stdout (for debugging purposes).

virtual void PrintModelImage()

Prints the current computed model image to stdout (for debugging purposes).

virtual void PrintWeights()

Prints the current weight image to stdout (for debugging purposes).

virtual void PrintMask()

Prints the input mask image to stdout (for debugging purposes).

virtual void PopulateParameterNames()

Note that this is usually called by AddFunctions() in add_functions.cpp.

int GetNFunctions()

Prints the total number of image functions (instances of FunctionObject subclasses) making up the model.

int GetNParams()

Prints the total number of parameters making up the model.

long GetNDataValues()

Prints the number of data values (pixels, masked or unmasked) in the data image.

virtual long GetNValidPixels()

Prints the number of valid (i.e., unmasked) data values (pixels) in the data image.

bool HasPSF()

Returns true if the model has a PSF image.

bool HasOversampledPSF()

Returns true if the model has one or more oversampled regions (with oversampled PSFs).

bool HasMask()

Returns true if a mask image exists.

virtual double *GetModelImageVector()

Returns a pointer to the model image (matching the data image in size if convolution is being done). If the model image has not yet been computed, returns NULL.

double *GetExpandedModelImageVector()

This differs from GetModelImageVector() in that it always returns the full model image, even in the case of PSF convolution (where the full model image will be larger than the data image!)

double *GetResidualImageVector()

Computes residual image (data - model) using current model image (usually best-fit image), and returns pointer to residual-image vector. Returns NULL if memory allocation failed.

double *GetWeightImageVector()

Returns the weightVector converted to 1/sigma^2 (i.e., 1/variance) form. Returns NULL if memory allocation failed.

double *GetDataVector()

Returns a pointer to the data image.

double FindTotalFluxes(double params[], int xSize, int ySize, double individualFluxes[])

Estimate total fluxes for individual components (and entire model) by integrating over a very large image, with each component/function centered in the image. Total flux is returned by the function; fluxes for individual components are returned in individualFluxes.

virtual int UseBootstrap()

Tells ModelObject1d object that from now on we’ll operate in bootstrap resampling mode, so that bootstrapIndices vector is used to access the data and model values (and weight values, if any). Returns the status from MakeBootstrapSample(), which will be -1 if memory allocation for the bootstrap-indices vector failed.

virtual int MakeBootstrapSample()

Generate a new bootstrap resampling of the data (more precisely, this generate a bootstrap resampling of the data indices) Returns -1 if memory allocation for the bootstrap indices vector failed, otherwise returns 0.

Protected Functions

bool CheckParamVector(int nParams, double paramVector[])

Returns true if all values in the parameter vector are finite.

bool CheckWeightVector()

Returns true if all pixels in the weight vector are finite and nonnegative.

bool VetDataVector()

Returns true if all non-masked pixels in the image data vector are finite; returns false if one or more are not, and prints an error message to stderr. ALSO sets any masked pixels which are non-finite to 0.

More simply: the purpose of this method is to check the data vector (profile or image) to ensure that all non-masked pixels are finite. Any non-finite pixels which are masked will be set = 0.