mrcal C API

Table of Contents

The C API consists of 3 headers:

Most usages would simply #include <mrcal.h>, and this would include all the headers. This is a C (not C++) library, so X macros are used in several places for templating.

The best documentation for the C interfaces is the comments in the headers themselves. The available functions are broken down into categories, and described in a bit more detail here.

Geometry structures

We have 3 structures in basic_geometry.h:

  • mrcal_point2_t: a vector containing 2 double-precision floating-point values. The elements can be accessed individually as .x and .y or as an array .xy[]
  • mrcal_point3_t: exactly like mrcal_point2_t, but 3-dimensional. A vector containing 3 double-precision floating-point values. The elements can be accessed individually as .x and .y and .z or as an array .xyz[]
  • mrcal_pose_t: an unconstrained 6-DOF pose. Contains two sub-structures:

Geometry functions

A number of utility functions are defined in poseutils.h. Each routine has two forms:

  • A mrcal_..._full() function that supports a non-contiguous memory layout for each input and output
  • A convenience mrcal_...() macro that wraps mrcal_..._full(), and expects contiguous data. This has many fewer arguments, and is easier to call

Each data argument (input or output) has several items in the argument list:

  • double* xxx: a pointer to the first element in the array
  • int xxx_stride0, int xxx_stride1, …: the strides, one per dimension

The strides are given in bytes, and work as expected. For a (for instance) 3-dimensional xxx, the element at xxx[i,j,k] would be accessible as

*(double*) &((char*)xxx)[ i*xxx_stride0 +
                          j*xxx_stride1 +
                          k*xxx_stride2 ]

These all have direct Python bindings. For instance mrcal.rt_from_Rt().

The poseutils.h header serves as the listing of available functions.

Lens models

The lens model structures are defined here:

  • mrcal_lensmodel_type_t: an enum decribing the lens model type. No configuration is stored here.
  • mrcal_lensmodel_t: a lens model type and the configuration parameters. The configuration lives in a union supporting all the known lens models
  • mrcal_lensmodel_metadata_t: some metadata that decribes a model type. These are inherent properties of a particular model type; answers questions like: Can this model project behind the camera? Does it have an intrinsics core? Does it have gradients implemented?

The Python API describes a lens model with a string that contains the model type and the configuration, while the C API stores the same information in a mrcal_lensmodel_t. So much of the functionality here is used to convert between the two. The listing of available functions is best given with the commented header (with the extraneous bits removed, and the X-macros expanded):

// parametric models have no extra configuration
typedef struct {} mrcal_LENSMODEL_PINHOLE__config_t;
typedef struct {} mrcal_LENSMODEL_OPENCV4__config_t;
// ... and the same for all the other configuration-less models

// Configuration for the splined stereographic models. Generated by an X-macro
typedef struct
{
    /* Maximum degree of each 1D polynomial. This is almost certainly 2 */
    /* (quadratic splines, C1 continuous) or 3 (cubic splines, C2 continuous) */
    uint16_t order;
    /* We have a Nx by Ny grid of control points */
    uint16_t Nx;
    uint16_t Ny;
    /* The horizontal field of view. Not including fov_y. It's proportional with */
    /* Ny and Nx */
    uint16_t fov_x_deg;
} mrcal_LENSMODEL_SPLINED_STEREOGRAPHIC__config_t;


// This lensmodel type selects the lens model, but does NOT provide the
// configuration. mrcal_lensmodel_t does that.
typedef enum
{ MRCAL_LENSMODEL_INVALID           = -2,
  MRCAL_LENSMODEL_INVALID_BADCONFIG = -1,
  // The rest, starting with 0

  // Generated by an X-macro
  // ...,
  MRCAL_LENSMODEL_PINHOLE,
  // ...,
  MRCAL_LENSMODEL_OPENCV4,
  // ...,
  MRCAL_LENSMODEL_SPLINED_STEREOGRAPHIC,
  // ... and so on for the other models
} mrcal_lensmodel_type_t;


// Defines a lens model: the type AND the configuration values
typedef struct
{
    // The type of lensmodel. This is an enum, selecting elements of
    // MRCAL_LENSMODEL_LIST (with "MRCAL_" prepended)
    mrcal_lensmodel_type_t type;

    // A union of all the possible configuration structures. We pick the
    // structure type based on the value of "type
    union
    {
        // Generated by an X-macro
        mrcal_LENSMODEL_PINHOLE__config_t               LENSMODEL_PINHOLE__config;
        mrcal_LENSMODEL_OPENCV4__config_t               LENSMODEL_OPENCV4__config;
        mrcal_LENSMODEL_SPLINED_STEREOGRAPHIC__config_t LENSMODEL_SPLINED_STEREOGRAPHIC__config;
        // ... and so on for the other models
    };
} mrcal_lensmodel_t;


// Return an array of strings listing all the available lens models
//
// These are all "unconfigured" strings that use "..." placeholders for any
// configuration values. Each returned string is a \0-terminated const char*. The
// end of the list is signified by a NULL pointer
const char* const* mrcal_supported_lensmodel_names( void ); // NULL-terminated array of char* strings


// Return true if the given mrcal_lensmodel_type_t specifies a valid lens model
bool mrcal_lensmodel_type_is_valid(mrcal_lensmodel_type_t t);


// Return a string describing a lens model.
//
// This function returns a static string. For models with no configuration, this
// is the FULL string for that model. For models with a configuration, the
// configuration values have "..." placeholders. These placeholders mean that
// the resulting strings do not define a lens model fully, and cannot be
// converted to a mrcal_lensmodel_t with mrcal_lensmodel_from_name()
//
// This is the inverse of mrcal_lensmodel_type_from_name()
const char* mrcal_lensmodel_name_unconfigured( mrcal_lensmodel_t model );


// Return a CONFIGURED string describing a lens model.
//
// This function generates a fully-configured string describing the given lens
// model. For models with no configuration, this is just the static string
// returned by mrcal_lensmodel_name_unconfigured(). For models that have a
// configuration, however, the configuration values are filled-in. The resulting
// string may be converted back into a mrcal_lensmodel_t by calling
// mrcal_lensmodel_from_name().
//
// This function writes the string into the given buffer "out". The size of the
// buffer is passed in the "size" argument. The meaning of "size" is as with
// snprintf(), which is used internally. Returns true on success
//
// This is the inverse of mrcal_lensmodel_from_name()
bool mrcal_lensmodel_name( char* out, int size, mrcal_lensmodel_t model );


// Parse the lens model type from a lens model name string
//
// The configuration is ignored. Thus this function works even if the
// configuration is missing or unparseable. Unknown model names return
// MRCAL_LENSMODEL_INVALID
//
// This is the inverse of mrcal_lensmodel_name_unconfigured()
mrcal_lensmodel_type_t mrcal_lensmodel_type_from_name( const char* name );


// Parse the full configured lens model from a lens model name string
//
// The lens mode type AND the configuration are read into a mrcal_lensmodel_t
// structure, which this function returns. Strings with valid model names but
// missing or unparseable configuration return
//
//   {.type = MRCAL_LENSMODEL_INVALID_BADCONFIG}.
//
// Any other errors result in some other invalid lensmodel.type values, which
// can be checked with mrcal_lensmodel_type_is_valid(lensmodel->type)
//
// This is the inverse of mrcal_lensmodel_name()
mrcal_lensmodel_t mrcal_lensmodel_from_name( const char* name );


// Each lens model type has some metadata that describes its inherent
// properties. These properties can be queried by calling
// mrcal_lensmodel_metadata().
typedef struct
{
    // generated by an X-macro

    /* If true, this model contains an "intrinsics core". This is described */
    /* in mrcal_intrinsics_core_t. If present, the 4 core parameters ALWAYS */
    /* appear at the start of a model's parameter vector                    */
    bool has_core :1;


    /* Whether a model is able to project points behind the camera          */
    /* (z<0 in the camera coordinate system). Models based on a pinhole     */
    /* projection (pinhole, OpenCV, CAHVOR(E)) cannot do this. models based */
    /* on a stereographic projection (stereographic, splined stereographic) */
    /* can                                                                  */
    bool can_project_behind_camera :1;
} mrcal_lensmodel_metadata_t;


// Return a structure containing a model's metadata
//
// The available metadata is described in the definition of the
// MRCAL_LENSMODEL_META_LIST() macro
mrcal_lensmodel_metadata_t mrcal_lensmodel_metadata( const mrcal_lensmodel_t m );


// Return the number of parameters required to specify a given lens model
//
// For models that have a configuration, the parameter count value generally
// depends on the configuration. For instance, splined models use the model
// parameters as the spline control points, so the spline density (specified in
// the configuration) directly affects how many parameters such a model requires
int mrcal_lensmodel_num_params( const mrcal_lensmodel_t m );


// Return the number of parameters needed in optimizing the given lens model
//
// This is identical to mrcal_lensmodel_num_params(), but takes into account the
// problem selections. Any intrinsics parameters locked down in the
// mrcal_problem_selections_t do NOT count towards the optimization parameters
int mrcal_num_intrinsics_optimization_params( mrcal_problem_selections_t problem_selections,
                                              mrcal_lensmodel_t m );


// Return the locations of x and y spline knots

// Splined models are defined by the locations of their control points. These
// are arranged in a grid, the size and density of which is set by the model
// configuration. We fill-in the x knot locations into ux[] and the y locations
// into uy[]. ux[] and uy[] must be large-enough to hold configuration->Nx and
// configuration->Ny values respectively.
//
// This function applies to splined models only. Returns true on success
bool mrcal_knots_for_splined_models( double* ux, double* uy,
                                     mrcal_lensmodel_t lensmodel);

Projections

The fundamental functions for projection and unprojection are defined here. mrcal_project() is the main routine that implements the "forward" direction, and is available for every camera model. This function can return gradients in respect to the coordinates of the point being projected and/or in respect to the intrinsics vector.

mrcal_unproject() is the reverse direction, and is implemented as a numerical optimization to reverse the projection operation. Naturally, this is much slower than mrcal_project(), and has no gradient reporting. Models that have no gradients implemented (CAHVORE only, as of this writing) do not support mrcal_unproject(). They may have a Python mrcal.unproject() implementation available that uses a slower optimization routine that uses numerical differences instead of analytical gradients.

mrcal_project_stereographic() and mrcal_unproject_stereographic() are available as special-case routines. These are used in analysis and not to represent any actual lenses.

The listing of available functions is best given with the commented header:

// Project the given camera-coordinate-system points
//
// Compute a "projection", a mapping of points defined in the camera coordinate
// system to their observed pixel coordinates. If requested, gradients are
// computed as well.
//
// We project N 3D points p to N 2D pixel coordinates q using the given
// lensmodel and intrinsics parameter values.
//
// if (dq_dp != NULL) we report the gradient dq/dp in a dense (N,2,3) array
// ((N,2) mrcal_point3_t objects).
//
// if (dq_dintrinsics != NULL) we report the gradient dq/dintrinsics in a dense
// (N,2,Nintrinsics) array. Note that splined models have very high Nintrinsics
// and very sparse gradients. THIS function reports the gradients densely,
// however, so it is inefficient for splined models.
//
// This function supports CAHVORE distortions only if we don't ask for any
// gradients
//
// Projecting out-of-bounds points (beyond the field of view) returns undefined
// values. Generally things remain continuous even as we move off the imager
// domain. Pinhole-like projections will work normally if projecting a point
// behind the camera. Splined projections clamp to the nearest spline segment:
// the projection will fly off to infinity quickly since we're extrapolating a
// polynomial, but the function will remain continuous.
bool mrcal_project( // out
                   mrcal_point2_t* q,
                   mrcal_point3_t* dq_dp,
                   double*         dq_dintrinsics,

                   // in
                   const mrcal_point3_t* p,
                   int N,
                   mrcal_lensmodel_t lensmodel,
                   // core, distortions concatenated
                   const double* intrinsics);


// Unproject the given pixel coordinates
//
// Compute an "unprojection", a mapping of pixel coordinates to the camera
// coordinate system.
//
// We unproject N 2D pixel coordinates q to N 3D direction vectors v using the
// given lensmodel and intrinsics parameter values. The returned vectors v are
// not normalized, and may have any length.

// This is the "reverse" direction, so an iterative nonlinear optimization is
// performed internally to compute this result. This is much slower than
// mrcal_project(). For OpenCV models specifically, OpenCV has
// cvUndistortPoints() (and cv2.undistortPoints()), but these are unreliable:
// https://github.com/opencv/opencv/issues/8811
//
// This function does NOT support CAHVORE
bool mrcal_unproject( // out
                     mrcal_point3_t* v,

                     // in
                     const mrcal_point2_t* q,
                     int N,
                     mrcal_lensmodel_t lensmodel,
                     // core, distortions concatenated
                     const double* intrinsics);


// Project the given camera-coordinate-system points using a stereographic model
//
// Compute a "projection", a mapping of points defined in the camera coordinate
// system to their observed pixel coordinates. If requested, gradients are
// computed as well.
//
// We project N 3D points p to N 2D pixel coordinates q using the stereographic
// model with the given intrinsics core.
//
// if (dq_dp != NULL) we report the gradient dq/dp in a dense (N,2,3) array
// ((N,2) mrcal_point3_t objects).
//
// This is a special case of mrcal_project(). Useful as part of data analysis,
// not to represent any real-world lens
void mrcal_project_stereographic( // output
                                 mrcal_point2_t* q,
                                 mrcal_point3_t* dq_dp,

                                  // input
                                 const mrcal_point3_t* p,
                                 int N,
                                 double fx, double fy,
                                 double cx, double cy);


// Unproject the given pixel coordinates using a stereographic model
//
// Compute an "unprojection", a mapping pixel coordinates to the camera
// coordinate system.
//
// We project N 2D pixel coordinates q to N 3D direction vectors v using the
// stereographic model with the given intrinsics core. The returned vectors v
// are not normalized, and may have any length.
//
// if (dv_dq != NULL) we report the gradient dv/dq in a dense (N,3,2) array
// ((N,3) mrcal_point2_t objects).
//
// This is a special case of mrcal_unproject(). Useful as part of data analysis,
// not to represent any real-world lens
void mrcal_unproject_stereographic( // output
                                   mrcal_point3_t* v,
                                   mrcal_point2_t* dv_dq,

                                   // input
                                   const mrcal_point2_t* q,
                                   int N,
                                   double fx, double fy,
                                   double cx, double cy);

Layout of the measurement and state vectors

The optimization routine tries to minimize the 2-norm of the measurement vector \(\vec x\) by moving around the state vector \(\vec p\).

We select which parts of the optimization problem we're solving by setting bits in the mrcal_problem_selections_t structure. This defines

  • Which elements of the optimization vector are locked-down, and which are given to the optimizer to adjust
  • Whether we apply regularization to stabilize the solution
  • Whether the chessboard should be assumed flat, or if we should optimize deformation factors

This structure is defined like this:

// Bits indicating which parts of the optimization problem being solved. We can
// ask mrcal to solve for ALL the lens parameters and ALL the geometry and
// everything else. OR we can ask mrcal to lock down some part of the
// optimization problem, and to solve for the rest. If any variables are locked
// down, we use their initial values passed-in to mrcal_optimize()
typedef struct
{
    // If true, we solve for the intrinsics core. Applies only to those models
    // that HAVE a core (fx,fy,cx,cy)
    bool do_optimize_intrinsics_core        : 1;

    // If true, solve for the non-core lens parameters
    bool do_optimize_intrinsics_distortions : 1;

    // If true, solve for the geometry of the cameras
    bool do_optimize_extrinsics             : 1;

    // If true, solve for the poses of the calibration object
    bool do_optimize_frames                 : 1;

    // If true, apply the regularization terms in the solver
    bool do_apply_regularization            : 1;

    // If true, optimize the shape of the calibration object
    bool do_optimize_calobject_warp         : 1;

    // Whether to try to find NEW outliers. The outliers given on
    // input are respected regardless
    bool do_apply_outlier_rejection         : 1;
} mrcal_problem_selections_t;

Thus the state vector may contain any of

  • The lens parameters
  • The geometry of the cameras
  • The geometry of the observed chessboards and discrete points
  • The chessboard shape

The measurement vector may contain

  • The errors in observations of the chessboards
  • The errors in observations of discrete points
  • The penalties in the solved point positions
  • The regularization terms

Given mrcal_problem_selections_t and a vector \(\vec p\) or \(\vec x\), it is useful to know where specific quantities lie inside those vectors. Here we have 4 sets of functions to answer such questions:

  • int mrcal_state_index_THING(): Returns the index in the state vector \(\vec p\) where the contiguous block of values describing the THING begins. THING is any of
    • intrinsics
    • extrinsics
    • frames
    • points
    • calobject_warp
  • int mrcal_num_states_THING(): Returns the number of values in the contiguous block in the state vector \(\vec p\) that describe the given THING. THING is any of
    • intrinsics
    • extrinsics
    • frames
    • points
    • calobject_warp
  • int mrcal_measurement_index_THING(): Returns the index in the measurement vector \(\vec x\) where the contiguous block of values describing the THING begins. THING is any of
    • boards
    • points
    • regularization
  • int mrcal_num_measurements_THING(): Returns the number of values in the contiguous block in the measurement vector \(\vec x\) that describe the given THING. THING is any of
    • boards
    • points
    • regularization

The function listing:

int mrcal_measurement_index_boards(int i_observation_board,
                                   int Nobservations_board,
                                   int Nobservations_point,
                                   int calibration_object_width_n,
                                   int calibration_object_height_n);
int mrcal_num_measurements_boards(int Nobservations_board,
                                  int calibration_object_width_n,
                                  int calibration_object_height_n);
int mrcal_measurement_index_points(int i_observation_point,
                                   int Nobservations_board,
                                   int Nobservations_point,
                                   int calibration_object_width_n,
                                   int calibration_object_height_n);
int mrcal_num_measurements_points(int Nobservations_point);
int mrcal_measurement_index_regularization(int Nobservations_board,
                                           int Nobservations_point,
                                           int calibration_object_width_n,
                                           int calibration_object_height_n);
int mrcal_num_measurements_regularization(int Ncameras_intrinsics, int Ncameras_extrinsics,
                                          int Nframes,
                                          int Npoints, int Npoints_fixed, int Nobservations_board,
                                          mrcal_problem_selections_t problem_selections,
                                          mrcal_lensmodel_t lensmodel);

int mrcal_num_measurements(int Nobservations_board,
                           int Nobservations_point,
                           int calibration_object_width_n,
                           int calibration_object_height_n,
                           int Ncameras_intrinsics, int Ncameras_extrinsics,
                           int Nframes,
                           int Npoints, int Npoints_fixed,
                           mrcal_problem_selections_t problem_selections,
                           mrcal_lensmodel_t lensmodel);

int mrcal_num_states(int Ncameras_intrinsics, int Ncameras_extrinsics,
                     int Nframes,
                     int Npoints, int Npoints_fixed, int Nobservations_board,
                     mrcal_problem_selections_t problem_selections,
                     mrcal_lensmodel_t lensmodel);
int mrcal_state_index_intrinsics(int icam_intrinsics,
                                 int Ncameras_intrinsics, int Ncameras_extrinsics,
                                 int Nframes,
                                 int Npoints, int Npoints_fixed, int Nobservations_board,
                                 mrcal_problem_selections_t problem_selections,
                                 mrcal_lensmodel_t lensmodel);
int mrcal_num_states_intrinsics(int Ncameras_intrinsics,
                                mrcal_problem_selections_t problem_selections,
                                mrcal_lensmodel_t lensmodel);
int mrcal_state_index_extrinsics(int icam_extrinsics,
                                 int Ncameras_intrinsics, int Ncameras_extrinsics,
                                 int Nframes,
                                 int Npoints, int Npoints_fixed, int Nobservations_board,
                                 mrcal_problem_selections_t problem_selections,
                                 mrcal_lensmodel_t lensmodel);
int mrcal_num_states_extrinsics(int Ncameras_extrinsics,
                                mrcal_problem_selections_t problem_selections);
int mrcal_state_index_frames(int iframe,
                             int Ncameras_intrinsics, int Ncameras_extrinsics,
                             int Nframes,
                             int Npoints, int Npoints_fixed, int Nobservations_board,
                             mrcal_problem_selections_t problem_selections,
                             mrcal_lensmodel_t lensmodel);
int mrcal_num_states_frames(int Nframes,
                            mrcal_problem_selections_t problem_selections);
int mrcal_state_index_points(int i_point,
                             int Ncameras_intrinsics, int Ncameras_extrinsics,
                             int Nframes,
                             int Npoints, int Npoints_fixed, int Nobservations_board,
                             mrcal_problem_selections_t problem_selections,
                             mrcal_lensmodel_t lensmodel);
int mrcal_num_states_points(int Npoints, int Npoints_fixed,
                            mrcal_problem_selections_t problem_selections);
int mrcal_state_index_calobject_warp(int Ncameras_intrinsics, int Ncameras_extrinsics,
                                     int Nframes,
                                     int Npoints, int Npoints_fixed, int Nobservations_board,
                                     mrcal_problem_selections_t problem_selections,
                                     mrcal_lensmodel_t lensmodel);
int mrcal_num_states_calobject_warp(mrcal_problem_selections_t problem_selections,
                                    int Nobservations_board);

State packing

The optimization routine works in the space of scaled parameters, and several functions are available to pack/unpack the state vector \(\vec p\):

// Scales a state vector to the packed, unitless form used by the optimizer
//
// In order to make the optimization well-behaved, we scale all the variables in
// the state and the gradients before passing them to the optimizer. The internal
// optimization library thus works only with unitless (or "packed") data.
//
// This function takes an (Nstate,) array of full-units values p[], and scales
// it to produce packed data. This function applies the scaling directly to the
// input array; the input is modified, and nothing is returned.
//
// This is the inverse of mrcal_unpack_solver_state_vector()
void mrcal_pack_solver_state_vector( // out, in
                                     double* p,

                                     // in
                                     int Ncameras_intrinsics, int Ncameras_extrinsics,
                                     int Nframes,
                                     int Npoints, int Npoints_fixed,
                                     mrcal_problem_selections_t problem_selections,
                                     const mrcal_lensmodel_t lensmodel);


// Scales a state vector from the packed, unitless form used by the optimizer
//
// In order to make the optimization well-behaved, we scale all the variables in
// the state and the gradients before passing them to the optimizer. The internal
// optimization library thus works only with unitless (or "packed") data.
//
// This function takes an (Nstate,) array of unitless values p[], and scales it
// to produce full-units data. This function applies the scaling directly to the
// input array; the input is modified, and nothing is returned.
//
// This is the inverse of mrcal_pack_solver_state_vector()
void mrcal_unpack_solver_state_vector( // out, in
                                       double* p, // unitless state on input,
                                                  // scaled, meaningful state on
                                                  // output

                                       // in
                                       int Ncameras_intrinsics, int Ncameras_extrinsics,
                                       int Nframes,
                                       int Npoints, int Npoints_fixed,
                                       mrcal_problem_selections_t problem_selections,
                                       const mrcal_lensmodel_t lensmodel);

Optimization

The mrcal optimization routines are defined in mrcal.h. There are two primary functions, each accessing a lot of functionality, and taking many arguments:

  • mrcal_optimize() is the entry point to the optimization routine. This function ingests the state, runs the optimization, and returns the optimal state in the same variables. The optimization routine tries out different values of the state vector by calling an optimization callback function to evaluate each one.
  • mrcal_optimizer_callback() provides access to the optimization callback function standalone, without being wrapped into the optimization loop

Helper structures

We define some structures to organize the input to these functions. Each observation has a mrcal_camera_index_t to identify the observing camera:

// Used to specify which camera is making an observation. The "intrinsics" index
// is used to identify a specific camera, while the "extrinsics" index is used
// to locate a camera in space. If I have a camera that is moving over time, the
// intrinsics index will remain the same, while the extrinsics index will change
typedef struct
{
    // indexes the intrinsics array
    int  intrinsics;
    // indexes the extrinsics array. -1 means "at coordinate system reference"
    int  extrinsics;
} mrcal_camera_index_t;

When solving a vanilla calibration problem, we have a set of stationary cameras observing a moving scene. By convention, in such a problem we set the reference coordinate system to camera 0, so that camera has no extrinsics. So in a vanilla calibration problem mrcal_camera_index_t.intrinsics will be in \([0, N_\mathrm{cameras})\) and mrcal_camera_index_t.extrinsics will always be mrcal_camera_index_t.intrinsics - 1.

When solving a vanilla structure-from-motion problem, we have a set of moving cameras observing a stationary scene. Here mrcal_camera_index_t.intrinsics would be in \([0, N_\mathrm{cameras})\) and mrcal_camera_index_t.extrinsics would be specify the camera pose, unrelated to mrcal_camera_index_t.intrinsics.

These are the limiting cases; anything in-between is allowed.

A board observation is defined by a mrcal_observation_board_t:

// An observation of a calibration board. Each "observation" is ONE camera
// observing a board
typedef struct
{
    // which camera is making this observation
    mrcal_camera_index_t icam;

    // indexes the "frames" array to select the pose of the calibration object
    // being observed
    int                  iframe;
} mrcal_observation_board_t;

And an observation of a discrete point is defined by a mrcal_observation_point_t:

// An observation of a discrete point. Each "observation" is ONE camera
// observing a single point in space
typedef struct
{
    // which camera is making this observation
    mrcal_camera_index_t icam;

    // indexes the "points" array to select the position of the point being
    // observed
    int                  i_point;

    // Observed pixel coordinates. This works just like elements of
    // observations_board_pool:
    //
    // .x, .y are the pixel observations
    // .z is the weight of the observation. Most of the weights are expected to
    // be 1.0. Less precise observations have lower weights.
    // .z<0 indicates that this is an outlier. This is respected on
    // input
    //
    // Unlike observations_board_pool, outlier rejection is NOT YET IMPLEMENTED
    // for points, so outlier points will NOT be found and reported on output in
    // .z<0
    mrcal_point3_t px;
} mrcal_observation_point_t;

Note that the details of the handling of discrete points may change in the future.

We have mrcal_problem_constants_t to define some details of the optimization problem. These are similar to mrcal_problem_selections_t, but consist of numerical values, rather than just bits. Currently this structure contains valid ranges for interpretation of discrete points. These may change in the future.

// Constants used in a mrcal optimization. This is similar to
// mrcal_problem_selections_t, but contains numerical values rather than just
// bits
typedef struct
{
    // The minimum distance of an observed discrete point from its observing
    // camera. Any observation of a point below this range will be penalized to
    // encourage the optimizer to move the point further away from the camera
    double  point_min_range;


    // The maximum distance of an observed discrete point from its observing
    // camera. Any observation of a point abive this range will be penalized to
    // encourage the optimizer to move the point closer to the camera
    double  point_max_range;
} mrcal_problem_constants_t;

The optimization function returns most of its output in the same memory as its input variables. A few metrics that don't belong there are returned in a separate mrcal_stats_t structure:

// This structure is returned by the optimizer, and contains some statistics
// about the optimization
typedef struct
{
    // generated by an X-macro

    /* The RMS error of the optimized fit at the optimum. Generally the residual */
    /* vector x contains error values for each element of q, so N observed pixels */
    /* produce 2N measurements: len(x) = 2*N. And the RMS error is */
    /*   sqrt( norm2(x) / N ) */
    double rms_reproj_error__pixels;

    /* How many pixel observations were thrown out as outliers. Each pixel */
    /* observation produces two measurements. Note that this INCLUDES any */
    /* outliers that were passed-in at the start */
    int Noutliers;
} mrcal_stats_t;

This contains some statistics describing the discovered optimal solution.

Arguments

The full prototypes of the two optimization functions are:

mrcal_stats_t
mrcal_optimize( // out
                // Each one of these output pointers may be NULL
                // Shape (Nstate,)
                double* p_packed,
                // used only to confirm that the user passed-in the buffer they
                // should have passed-in. The size must match exactly
                int buffer_size_p_packed,

                // Shape (Nmeasurements,)
                double* x,
                // used only to confirm that the user passed-in the buffer they
                // should have passed-in. The size must match exactly
                int buffer_size_x,

                // out, in

                // These are a seed on input, solution on output

                // intrinsics is a concatenation of the intrinsics core and the
                // distortion params. The specific distortion parameters may
                // vary, depending on lensmodel, so this is a variable-length
                // structure
                double*             intrinsics,         // Ncameras_intrinsics * NlensParams
                mrcal_pose_t*       extrinsics_fromref, // Ncameras_extrinsics of these. Transform FROM the reference frame
                mrcal_pose_t*       frames_toref,       // Nframes of these.    Transform TO the reference frame
                mrcal_point3_t*     points,             // Npoints of these.    In the reference frame
                mrcal_point2_t*     calobject_warp,     // 1 of these. May be NULL if !problem_selections.do_optimize_calobject_warp

                // in
                int Ncameras_intrinsics, int Ncameras_extrinsics, int Nframes,
                int Npoints, int Npoints_fixed, // at the end of points[]

                const mrcal_observation_board_t* observations_board,
                const mrcal_observation_point_t* observations_point,
                int Nobservations_board,
                int Nobservations_point,

                // All the board pixel observations, in an array of shape
                //
                // ( Nobservations_board,
                //   calibration_object_height_n,
                //   calibration_object_width_n )
                //
                // .x, .y are the
                // pixel observations .z is the weight of the observation. Most
                // of the weights are expected to be 1.0. Less precise
                // observations have lower weights.
                //
                // .z<0 indicates that this is an outlier. This is respected on
                // input (even if !do_apply_outlier_rejection). New outliers are
                // marked with .z<0 on output, so this isn't const
                mrcal_point3_t* observations_board_pool,

                mrcal_lensmodel_t lensmodel,
                double observed_pixel_uncertainty,
                const int* imagersizes, // Ncameras_intrinsics*2 of these
                mrcal_problem_selections_t       problem_selections,
                const mrcal_problem_constants_t* problem_constants,
                double calibration_object_spacing,
                int calibration_object_width_n,
                int calibration_object_height_n,
                bool verbose,

                bool check_gradient);


// This is cholmod_sparse. I don't want to include the full header that defines
// it in mrcal.h, and I don't need to: mrcal.h just needs to know that it's a
// structure
struct cholmod_sparse_struct;

// Evaluate the value of the callback function at the given operating point
//
// The main optimization routine in mrcal_optimize() searches for optimal
// parameters by repeatedly calling a function to evaluate each hypothethical
// parameter set. This evaluation function is available by itself here,
// separated from the optimization loop. The arguments are largely the same as
// those to mrcal_optimize(), but the inputs are all read-only It is expected
// that this will be called from Python only.
bool mrcal_optimizer_callback(// out

                             // These output pointers may NOT be NULL, unlike
                             // their analogues in mrcal_optimize()

                             // Shape (Nstate,)
                             double* p_packed,
                             // used only to confirm that the user passed-in the buffer they
                             // should have passed-in. The size must match exactly
                             int buffer_size_p_packed,

                             // Shape (Nmeasurements,)
                             double* x,
                             // used only to confirm that the user passed-in the buffer they
                             // should have passed-in. The size must match exactly
                             int buffer_size_x,

                             // output Jacobian. May be NULL if we don't need
                             // it. This is the unitless Jacobian, used by the
                             // internal optimization routines
                             struct cholmod_sparse_struct* Jt,


                             // in

                             // intrinsics is a concatenation of the intrinsics core
                             // and the distortion params. The specific distortion
                             // parameters may vary, depending on lensmodel, so
                             // this is a variable-length structure
                             const double*             intrinsics,         // Ncameras_intrinsics * NlensParams
                             const mrcal_pose_t*       extrinsics_fromref, // Ncameras_extrinsics of these. Transform FROM the reference frame
                             const mrcal_pose_t*       frames_toref,       // Nframes of these.    Transform TO the reference frame
                             const mrcal_point3_t*     points,             // Npoints of these.    In the reference frame
                             const mrcal_point2_t*     calobject_warp,     // 1 of these. May be NULL if !problem_selections.do_optimize_calobject_warp

                             int Ncameras_intrinsics, int Ncameras_extrinsics, int Nframes,
                             int Npoints, int Npoints_fixed, // at the end of points[]

                             const mrcal_observation_board_t* observations_board,
                             const mrcal_observation_point_t* observations_point,
                             int Nobservations_board,
                             int Nobservations_point,

                             // All the board pixel observations, in an array of shape
                             //
                             // ( Nobservations_board,
                             //   calibration_object_height_n,
                             //   calibration_object_width_n )
                             //
                             // .x, .y are the pixel observations .z is the
                             // weight of the observation. Most of the weights
                             // are expected to be 1.0. Less precise
                             // observations have lower weights.
                             //
                             // .z<0 indicates that this is an outlier
                             const mrcal_point3_t* observations_board_pool,

                             mrcal_lensmodel_t lensmodel,
                             double observed_pixel_uncertainty,
                             const int* imagersizes, // Ncameras_intrinsics*2 of these
                             mrcal_problem_selections_t       problem_selections,
                             const mrcal_problem_constants_t* problem_constants,
                             double calibration_object_spacing,
                             int calibration_object_width_n,
                             int calibration_object_height_n,
                             bool verbose);

Most of the arguments to mrcal_optimize() and mrcal_optimizer_callback() represent an optimization state, so these two functions accept a very similar set of arguments.

The output buffers are given as two arguments: the buffer pointer itself and a int buffer_size_.... to describe the size of the given buffer. This exists purely for error-checking: mrcal knows how big these buffers should be, and it makes sure that the given buffer is of the correct size. If it doesn't match, an error is reported.

The resulting state and measurement vectors are returned in p_packed and x respectively.

mrcal_optimizer_callback() also returns the jacobian at the operating point in the Jt array. This is large and sparse, so it is stored in a cholmod_sparse structure from CHOLMOD in the suitesparse project. CHOLMOD uses a FORTRAN-style column-major representation, so from CHOLMOD's point of view we're returning \(J^T\) and not \(J\).

The optimization state is given in the intrinsics, extrinsics_fromref, frames_toref, points, calobject_warp, arguments. These are const inputs to mrcal_optimizer_callback(). In calls to mrcal_optimize() these are the optimization seed, and the optimized results are reported in the same arrays.

The integers Ncameras_intrinsics, Ncameras_extrinsics, Nframes, and Npoints, denote the respective lengths of the arrays intrinsics, extrinsics_fromref, frames and observations_point.

Npoints_fixed denotes how many points at the end of the points array are fixed, and not optimized. The usual value of 0 indicates that all points should be optimized. This logic is likely to change in the future.

The observations_board and observations_point arrays describe the observations. Each structure element indicate which camera (intrinsics and extrinsics) made the corresponding observation. The observations_point array contains the actual observed pixels and weights, which the observed chessboard pixels and weights live in a separate array: observations_board_pool.

The integers Nobservations_board and Nobservations_point set the sizes of the arrays observations_board and observations_point respectively. These are different from Nframes and Npoints because one frame (or point) may be observed by multiple cameras, producing larger Nobservations_board (or Nobservations_point).

The actual chessboard observations are given in observations_board_pool, an array of shape (Nobservations_board, calibration_object_height_n, calibration_object_width_n) containing mrcal_point3_t elements. The observed pixel coordinates are in .x and .y, and the observation weight is in .z. .z < 0 means this point should be ignored; this can be used to specify incomplete board observations. If we're calling mrcal_optimize() and mrcal_problem_selections_t.do_apply_outlier_rejection, then in addition to the ignored points on input the outlier rejection algorithm will be active, and upon return the outlier points will be marked with .z < 0.

The lens model of all the cameras is specified in the lensmodel argument.

The expected uncertainty of the pixel observations is given in observed_pixel_uncertainty. This is used primarily by the caller to estimate uncertainties. These functions currently use this value only in the outlier rejection routine.

The imagersizes array contains the (width,height) dimensions of the imager of all the cameras. Each camera is specified separately, so Ncameras_intrinsics*2 integers must be passed in. Currently these C routines use this value only in the regularization computation.

The problem_selections and problem_constants define the problem details, as described above.

The calibration_object_spacing and calibration_object_width_n and calibration_object_height_n arguments define the dimensions of the chessboard. A regular grid of points is expected.

If we want verbose reporting about what the optimizer is doing, pass verbose = true to mrcal_optimize().

Camera model reading/writing

Currently there's no support for reading/writing .cameramodel files in the C API. This is already partially implemented, and I will finish it when I need it or when somebody bugs me about it, whichever comes first.

Miscellaneous

When calibrating cameras, each observations is associated with some intrinsics and some extrinsics. Those two chunks of data live in different parts of the optimization vector, and are indexed independently. If we have stationary cameras, then each set of camera intrinsics is associated with exactly one set of camera extrinsics, and we can use this function to query this correspondence.

The arguments are the same as the ones to mrcal_optimize(). The output index is returned in icam_extrinsics. If this camera was used to define the reference coordinate system, this camera has no explicit extrinsics, and we set icam_extrinsics = -1.

We return true on success. If we have moving cameras, then a single physical camera would have one set of intrinsics but many different extrinsics, and this function will fail, returning false.

The prototype:

// Reports the icam_extrinsics corresponding to a given icam_intrinsics.
//
// If we're solving a calibration problem (stationary cameras observing a moving
// calibration object), each camera has a unique intrinsics vector and a unique
// extrinsics vector. And this function reports the latter, given the former. On
// success, the result is written to *icam_extrinsics, and we return true. If
// the given camera is at the reference coordinate system, it has no extrinsics,
// and we report -1.
//
// If we have moving cameras, there won't be a single icam_extrinsics for a
// given icam_intrinsics, and we report an error by returning false
bool mrcal_corresponding_icam_extrinsics(// out
                                         int* icam_extrinsics,

                                         // in
                                         int icam_intrinsics,
                                         int Ncameras_intrinsics,
                                         int Ncameras_extrinsics,
                                         int Nobservations_board,
                                         const mrcal_observation_board_t* observations_board,
                                         int Nobservations_point,
                                         const mrcal_observation_point_t* observations_point);