C API#


Installation#

Please follow the Getting Started guide to install BridgeStan’s pre-requisites and downloaded a copy of the BridgeStan source code.

Example Program#

An example program is provided alongside the BridgeStan source in c-example/. Details for building the example can be found in c-example/Makefile.

Show example.c
#include "bridgestan.h"
#include <stdio.h>

int main(int argc, char** argv) {
  printf("Using BridgeStan version %d.%d.%d\n", bs_major_version,
         bs_minor_version, bs_patch_version);

  char* data;
  if (argc > 1) {
    data = argv[1];
  } else {
    data = NULL;
  }

  // this could potentially error, and we may get information back about why.
  char* err;
  bs_model* model = bs_model_construct(data, 123, &err);
  if (!model) {
    if (err) {
      printf("Error: %s", err);
      bs_free_error_msg(err);
    }
    return 1;
  }

  printf("This model's name is %s.\n", bs_name(model));
  printf("It has %d parameters.\n", bs_param_num(model, 0, 0));

  bs_model_destruct(model);
  return 0;
}

API Reference#

The following are the C functions exposed by the BridgeStan library in bridgestan.h. These are wrapped in the various high-level interfaces.

These functions are implemented in C++, see How It Works for more details.

Functions

bs_model *bs_model_construct(const char *data, unsigned int seed, char **error_msg)#

Construct an instance of a model wrapper. Data must be encoded in JSON as indicated in the CmdStan Reference Manual.

Parameters:
  • data[in] C-style string. This is either a path to JSON-encoded data file (must end with “.json”), a JSON string literal, or nullptr. An empty string or null pointer are both interpreted as no data.

  • seed[in] seed for PRNG used during model construction. This PRNG is used for RNG functions in the transformed data block of the model, and then discarded.

  • error_msg[out] a pointer to a string that will be allocated if there is an error. This must later be freed by calling bs_free_error_msg().

Returns:

pointer to constructed model or nullptr if construction fails

void bs_model_destruct(bs_model *m)#

Destroy the model.

Parameters:

m[in] pointer to model structure

void bs_free_error_msg(char *error_msg)#

Free the error messages created by other methods.

Parameters:

error_msg[in] pointer to error message

const char *bs_name(const bs_model *m)#

Return the name of the specified model as a C-style string.

The returned string should not be modified; it is freed when the model wrapper is destroyed.

Parameters:

m[in] pointer to model and RNG structure

Returns:

name of model

const char *bs_model_info(const bs_model *m)#

Return information about the compiled model as a C-style string.

The returned string should not be modified; it is freed when the model wrapper is destroyed.

Parameters:

m[in] pointer to model structure

Returns:

Information about the model including Stan version, Stan defines, and compiler flags.

const char *bs_param_names(const bs_model *m, bool include_tp, bool include_gq)#

Return a comma-separated sequence of indexed parameter names, including the transformed parameters and/or generated quantities as specified.

The parameters are returned in the order they are declared. Multivariate parameters are return in column-major (more generally last-index major) order. Parameters are separated with periods (.). For example, a[3] is written a.3 and b[2, 3] as b.2.3. The numbering follows Stan and is indexed from 1.

The returned string should not be modified; it is freed when the model wrapper is destroyed.

Parameters:
  • m[in] pointer to model structure

  • include_tp[in] true to include transformed parameters

  • include_gq[in] true to include generated quantities

Returns:

CSV-separated, indexed, parameter names

const char *bs_param_unc_names(const bs_model *m)#

Return a comma-separated sequence of unconstrained parameters. Only parameters are unconstrained, so there are no unconstrained transformed parameters or generated quantities.

The parameters are returned in the order they are declared. Multivariate parameters are return in column-major (more generally last-index major) order. Parameters are separated with periods (.). For example, a[3] is written a.3 and b[2, 3] as b.2.3. The numbering follows Stan and is indexed from 1.

The returned string should not be modified; it is freed when the model wrapper is destroyed.

Parameters:

m[in] pointer to model structure

Returns:

CSV-separated, indexed, unconstrained parameter names

int bs_param_num(const bs_model *m, bool include_tp, bool include_gq)#

Return the number of scalar parameters, optionally including the number of transformed parameters and/or generated quantities. For example, a 2 x 3 matrix counts as 6 scalar parameters.

Parameters:
  • m[in] pointer to model structure

  • include_tp[in] true to include transformed parameters

  • include_gq[in] true to include generated quantities

Returns:

number of parameters

int bs_param_unc_num(const bs_model *m)#

Return the number of unconstrained parameters. The number of unconstrained parameters might be smaller than the number of parameters if the unconstrained space has fewer dimensions than the constrained (e.g., for simplexes or correlation matrices).

Parameters:

m[in] pointer to model structure

Returns:

number of unconstrained parameters

int bs_param_constrain(const bs_model *m, bool include_tp, bool include_gq, const double *theta_unc, double *theta, bs_rng *rng, char **error_msg)#

Set the sequence of constrained parameters based on the specified unconstrained parameters, including transformed parameters and/or generated quantities as specified, and return a return code of 0 for success and -1 for failure. Parameter order is as declared in the Stan program, with multivariate parameters given in last-index-major order.

Parameters:
  • m[in] pointer to model structure

  • include_tp[in] true to include transformed parameters

  • include_gq[in] true to include generated quantities

  • theta_unc[in] sequence of unconstrained parameters

  • theta[out] sequence of constrained parameters

  • rng[in] pointer to pseudorandom number generator, should be created by bs_rng_construct(). This is only required when include_gq is true, otherwise it can be null.

  • error_msg[out] a pointer to a string that will be allocated if there is an error. This must later be freed by calling bs_free_error_msg().

Returns:

code 0 if successful and code -1 if there is an exception in the underlying Stan code

int bs_param_unconstrain(const bs_model *m, const double *theta, double *theta_unc, char **error_msg)#

Set the sequence of unconstrained parameters based on the specified constrained parameters, and return a return code of 0 for success and -1 for failure. Parameter order is as declared in the Stan program, with multivariate parameters given in last-index-major order.

Parameters:
  • m[in] pointer to model structure

  • theta[in] sequence of constrained parameters

  • theta_unc[out] sequence of unconstrained parameters

  • error_msg[out] a pointer to a string that will be allocated if there is an error. This must later be freed by calling bs_free_error_msg().

Returns:

code 0 if successful and code -1 if there is an exception in the underlying Stan code

int bs_param_unconstrain_json(const bs_model *m, const char *json, double *theta_unc, char **error_msg)#

Set the sequence of unconstrained parameters based on the JSON specification of the constrained parameters, and return a return code of 0 for success and -1 for failure. The JSON schema assumed is fully defined in the CmdStan Reference Manual.

Parameters:
  • m[in] pointer to model structure

  • json[in] JSON-encoded constrained parameters

  • theta_unc[out] sequence of unconstrained parameters

  • error_msg[out] a pointer to a string that will be allocated if there is an error. This must later be freed by calling bs_free_error_msg().

Returns:

code 0 if successful and code -1 if there is an exception in the underlying Stan code

int bs_log_density(const bs_model *m, bool propto, bool jacobian, const double *theta_unc, double *lp, char **error_msg)#

Set the log density of the specified parameters, dropping constants if propto is true and including the Jacobian terms resulting from constraining parameters if jacobian is true, and return a return code of 0 for success and -1 if there is an exception executing the Stan program.

Parameters:
  • m[in] pointer to model structure

  • propto[in] true to discard constant terms

  • jacobian[in] true to include change-of-variables terms

  • theta_unc[in] unconstrained parameters

  • lp[out] log density to be set

  • error_msg[out] a pointer to a string that will be allocated if there is an error. This must later be freed by calling bs_free_error_msg().

Returns:

code 0 if successful and code -1 if there is an exception in the underlying Stan code

int bs_log_density_gradient(const bs_model *m, bool propto, bool jacobian, const double *theta_unc, double *val, double *grad, char **error_msg)#

Set the log density and gradient of the specified parameters, dropping constants if propto is true and including the Jacobian terms resulting from constraining parameters if jacobian is true, and return a return code of 0 for success and -1 if there is an exception executing the Stan program. The gradient must have enough space to hold the gradient.

The gradients are computed using automatic differentiation.

Parameters:
  • m[in] pointer to model structure

  • propto[in] true to discard constant terms

  • jacobian[in] true to include change-of-variables terms

  • theta_unc[in] unconstrained parameters

  • val[out] log density to be set

  • grad[out] gradient to set

  • error_msg[out] a pointer to a string that will be allocated if there is an error. This must later be freed by calling bs_free_error_msg().

Returns:

code 0 if successful and code -1 if there is an exception in the underlying Stan code

int bs_log_density_hessian(const bs_model *m, bool propto, bool jacobian, const double *theta_unc, double *val, double *grad, double *hessian, char **error_msg)#

Set the log density, gradient, and Hessian of the specified parameters, dropping constants if propto is true and including the Jacobian terms resulting from constraining parameters if jacobian is true, and return a return code of 0 for success and -1 if there is an exception executing the Stan program. The pointer grad must have enough space to hold the gradient. The pointer hessian must have enough space to hold the Hessian.

The gradients are computed using automatic differentiation. Hessians are computed using nested automatic differentiation if BRIDGESTAN_AD_HESSIAN is defined, otherwise they are computed using central finite differences of size(theta_unc) calculations of gradient.

Parameters:
  • m[in] pointer to model structure

  • propto[in] true to discard constant terms

  • jacobian[in] true to include change-of-variables terms

  • theta_unc[in] unconstrained parameters

  • val[out] log density to be set

  • grad[out] gradient to set

  • hessian[out] hessian to set

  • error_msg[out] a pointer to a string that will be allocated if there is an error. This must later be freed by calling bs_free_error_msg().

Returns:

code 0 if successful and code -1 if there is an exception in the underlying Stan code

int bs_log_density_hessian_vector_product(const bs_model *m, bool propto, bool jacobian, const double *theta_unc, const double *vector, double *val, double *hvp, char **error_msg)#

Calculate the log density and the product of the Hessian with the specified vector for the specified unconstrained parameters and write it into the specified value pointer and Hessian-vector product pointer, dropping constants it propto is true and including the Jacobian adjustment if jacobian is true. Returns a return code of 0 for success and -1 if there is an exception executing the Stan program. The pointer hvp must have enough space to hold the product.

Hessian-vector-products are computed using nested automatic differentiation if BRIDGESTAN_AD_HESSIAN is defined, otherwise they are computed using central finite differences of the gradient of theta_unc perturbed in the direction of vector. This approximates the Hessian-vector product using two gradient evaluations, but at a lower accuracy than the nested automatic differentiation.

Parameters:
  • m[in] pointer to model structure

  • propto[in] true to drop constant terms

  • jacobian[in] true to include Jacobian adjustment for constrained parameter transforms

  • theta_unc[in] unconstrained parameters

  • vector[in] vector to multiply Hessian by

  • val[out] log density to set

  • hvp[out] Hessian-vector to set

  • error_msg[out] a pointer to a string that will be allocated if there is an error. This must later be freed by calling bs_free_error_msg().

Returns:

code 0 if successful and code -1 if there is an exception in the underlying Stan code

bs_rng *bs_rng_construct(unsigned int seed, char **error_msg)#

Construct an PRNG object to be used in bs_param_constrain(). This object is not thread safe and should be constructed and destructed for each thread.

Parameters:
  • seed[in] seed for the RNG

  • error_msg[out] a pointer to a string that will be allocated if there is an error. This must later be freed by calling bs_free_error_msg().

void bs_rng_destruct(bs_rng *rng)#

Destruct an RNG object.

Parameters:

rng[in] pointer to RNG object

int bs_set_print_callback(STREAM_CALLBACK callback, char **error_msg)#

Provide a function for printing. This will be called when the Stan model prints output. The default is to print to stdout.

Parameters:
  • callback[in] function to call when the Stan model prints. This function will be guarded by a mutex, so it need not be thread safe. It must never propagate an exception. Passing NULL will redirect printing back to stdout.

  • error_msg[out] a pointer to a string that will be allocated if there is an error. This must later be freed by calling bs_free_error_msg().

Returns:

code 0 if successful and code -1 if there is an exception

Typedefs

typedef struct bs_model bs_model#

Opaque type for model.

typedef struct bs_rng bs_rng#

Opaque type for RNG.

typedef void (*STREAM_CALLBACK)(const char *data, size_t size)#

Type signature for optional print callback

Variables

int bs_major_version#

Version information for the BridgeStan library.

Note

These are not the version of the wrapped Stan library.

Note

These were not available pre-2.0.0, so their absence implies the library is in the 1.0.x series

int bs_minor_version#
int bs_patch_version#

R-compatible functions#

To support calling these functions from R without including R-specific headers into the project, the following functions are exposed in bridgestanR.h.

These are small shims which call the above functions. All arguments and return values must be handled via pointers.

Functions

void bs_model_construct_R(char **data, int *rng, bs_model **ptr_out, char **err_msg, void **err_ptr)#

See bs_model_construct() for more details.

void bs_version_R(int *major, int *minor, int *patch)#

See bs_major_version for more details.

void bs_model_destruct_R(bs_model **model)#

See bs_model_destruct() for more details.

void bs_free_error_msg_R(void **err_msg)#

Free error message allocated in C++. Because R performs copies at the boundary on char**s, this uses void** pointing to the same memory.

See bs_free_error_msg() for more details.

void bs_name_R(bs_model **model, char const **name_out)#

See bs_name() for more details.

void bs_model_info_R(bs_model **model, char const **info_out)#

See bs_model_info() for more details.

void bs_param_names_R(bs_model **model, int *include_tp, int *include_gq, char const **name_out)#

See bs_param_names() for more details.

void bs_param_unc_names_R(bs_model **model, char const **name_out)#

See bs_param_unc_names() for more details.

void bs_param_num_R(bs_model **model, int *include_tp, int *include_gq, int *num_out)#

See bs_param_num() for more details.

void bs_param_unc_num_R(bs_model **model, int *num_out)#

See bs_param_unc_num() for more details.

void bs_param_constrain_R(bs_model **model, int *include_tp, int *include_gq, const double *theta_unc, double *theta, bs_rng **rng, int *return_code, char **err_msg, void **err_ptr)#

See bs_param_constrain() for more details.

void bs_param_unconstrain_R(bs_model **model, const double *theta, double *theta_unc, int *return_code, char **err_msg, void **err_ptr)#

See bs_param_unconstrain() for more details.

void bs_param_unconstrain_json_R(bs_model **model, char const **json, double *theta_unc, int *return_code, char **err_msg, void **err_ptr)#

See bs_param_unconstrain_json() for more details.

void bs_log_density_R(bs_model **model, int *propto, int *jacobian, const double *theta_unc, double *val, int *return_code, char **err_msg, void **err_ptr)#

See bs_log_density() for more details.

void bs_log_density_gradient_R(bs_model **model, int *propto, int *jacobian, const double *theta_unc, double *val, double *grad, int *return_code, char **err_msg, void **err_ptr)#

See bs_log_density_gradient() for more details.

void bs_log_density_hessian_R(bs_model **model, int *propto, int *jacobian, const double *theta_unc, double *val, double *grad, double *hess, int *return_code, char **err_msg, void **err_ptr)#

See bs_log_density_hessian() for more details.

void bs_log_density_hessian_vector_product_R(bs_model **model, int *propto, int *jacobian, const double *theta_unc, const double *vector, double *val, double *Hvp, int *return_code, char **err_msg, void **err_ptr)#

See bs_log_density_hessian_vector_product() for more details.

void bs_rng_construct_R(int *seed, bs_rng **ptr_out, char **err_msg, void **err_ptr)#

See bs_rng_construct() for more details.

void bs_rng_destruct_R(bs_rng **rng)#

See bs_rng_destruct() for more details.