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) {
char* data;
if (argc > 1) {
data = argv[1];
} else {
data = "";
}
bs_model_rng* model = bs_construct(data, 123, 0);
if (!model) {
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));
return bs_destruct(model);
}
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_rng *bs_construct(char *data_file, unsigned int seed, unsigned int chain_id)#
Construct an instance of a model and pseudorandom number generator (PRNG) wrapper. Data must be encoded in JSON as indicated in the CmdStan Reference Manual.
- Parameters
data_file – [in] C-style string. This is either a path to JSON-encoded data file (must end with “.json”), or a JSON string literal.
seed – [in] seed for PRNG
chain_id – [in] identifier for concurrent sequence of PRNG draws
- Returns
pointer to constructed model or
nullptr
if construction fails
-
int bs_destruct(bs_model_rng *mr)#
Destroy the model and return 0 for success and -1 if there is an exception while freeing memory.
- Parameters
mr – [in] pointer to model and RNG structure
- Returns
0 for success and -1 if there is an exception freeing one of the model components.
-
const char *bs_name(bs_model_rng *mr)#
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 and RNG wrapper is destroyed.
- Parameters
mr – [in] pointer to model and RNG structure
- Returns
name of model
-
const char *bs_model_info(bs_model_rng *mr)#
Return information about the compiled model as a C-style string.
The returned string should not be modified; it is freed when the model and RNG wrapper is destroyed.
- Parameters
mr – [in] pointer to model and RNG structure
- Returns
Information about the model including Stan version, Stan defines, and compiler flags.
-
const char *bs_param_names(bs_model_rng *mr, 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 writtena.3
andb[2, 3]
asb.2.3
. The numbering follows Stan and is indexed from 1.The returned string should not be modified; it is freed when the model and RNG wrapper is destroyed.
- Parameters
mr – [in] pointer to model and RNG structure
include_tp – [in]
true
to include transformed parametersinclude_gq – [in]
true
to include generated quantities
- Returns
CSV-separated, indexed, parameter names
-
const char *bs_param_unc_names(bs_model_rng *mr)#
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 writtena.3
andb[2, 3]
asb.2.3
. The numbering follows Stan and is indexed from 1.The returned string should not be modified; it is freed when the model and RNG wrapper is destroyed.
- Parameters
mr – [in] pointer to model and RNG structure
- Returns
CSV-separated, indexed, unconstrained parameter names
-
int bs_param_num(bs_model_rng *mr, 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
mr – [in] pointer to model and RNG structure
include_tp – [in]
true
to include transformed parametersinclude_gq – [in]
true
to include generated quantities
- Returns
number of parameters
-
int bs_param_unc_num(bs_model_rng *mr)#
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
mr – [in] pointer to model and RNG structure
- Returns
number of unconstrained parameters
-
int bs_param_constrain(bs_model_rng *mr, bool include_tp, bool include_gq, const double *theta_unc, double *theta)#
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
mr – [in] pointer to model and RNG structure
include_tp – [in]
true
to include transformed parametersinclude_gq – [in]
true
to include generated quantitiestheta_unc – [in] sequence of unconstrained parameters
theta – [out] sequence of constrained parameters
- Returns
code 0 if successful and code -1 if there is an exception in the underlying Stan code
-
int bs_param_unconstrain(bs_model_rng *mr, const double *theta, double *theta_unc)#
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
mr – [in] pointer to model and RNG structure
theta – [in] sequence of constrained parameters
theta_unc – [out] sequence of unconstrained parameters
- Returns
code 0 if successful and code -1 if there is an exception in the underlying Stan code
-
int bs_param_unconstrain_json(bs_model_rng *mr, const char *json, double *theta_unc)#
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. Parameter order is as declared in the Stan program, with multivariate parameters given in last-index-major order. The JSON schema assumed is fully defined in the CmdStan Reference Manual.
- Parameters
mr – [in] pointer to model and RNG structure
json – [in] JSON-encoded constrained parameters
theta_unc – [out] sequence of unconstrained parameters
- Returns
code 0 if successful and code -1 if there is an exception in the underlying Stan code
-
int bs_log_density(bs_model_rng *mr, bool propto, bool jacobian, const double *theta, double *lp)#
Set the log density of the specified parameters, dropping constants if
propto
istrue
and including the Jacobian terms resulting from constraining parameters ifjacobian
istrue
, and return a return code of 0 for success and -1 if there is an exception executing the Stan program.- Parameters
mr – [in] pointer to model and RNG structure
propto – [in]
true
to discard constant termsjacobian – [in]
true
to include change-of-variables termstheta – [in] unconstrained parameters
lp – [out] log density to be set
- Returns
code 0 if successful and code -1 if there is an exception in the underlying Stan code
-
int bs_log_density_gradient(bs_model_rng *mr, bool propto, bool jacobian, const double *theta, double *val, double *grad)#
Set the log density and gradient of the specified parameters, dropping constants if
propto
istrue
and including the Jacobian terms resulting from constraining parameters ifjacobian
istrue
, 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
mr – [in] pointer to model and RNG structure
propto – [in]
true
to discard constant termsjacobian – [in]
true
to include change-of-variables termstheta – [in] unconstrained parameters
val – [out] log density to be set
grad – [out] gradient to set
- Returns
code 0 if successful and code -1 if there is an exception in the underlying Stan code
-
int bs_log_density_hessian(bs_model_rng *mr, bool propto, bool jacobian, const double *theta, double *val, double *grad, double *hessian)#
Set the log density, gradient, and Hessian of the specified parameters, dropping constants if
propto
istrue
and including the Jacobian terms resulting from constraining parameters ifjacobian
istrue
, and return a return code of 0 for success and -1 if there is an exception executing the Stan program. The pointergrad
must have enough space to hold the gradient. The pointerHessian
must have enough space to hold the Hessian.The gradients are computed using automatic differentiation. the Hessians are
- Parameters
mr – [in] pointer to model and RNG structure
propto – [in]
true
to discard constant termsjacobian – [in]
true
to include change-of-variables termstheta – [in] unconstrained parameters
val – [out] log density to be set
grad – [out] gradient to set
hessian – [out] hessian to set
- Returns
code 0 if successful and code -1 if there is an exception in the underlying Stan code
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 handeled via pointers.
Functions
-
void bs_construct_R(char **data, int *rng, int *chain, bs_model_rng **ptr_out)#
-
void bs_destruct_R(bs_model_rng **model, int *return_code)#
-
void bs_name_R(bs_model_rng **model, char const **name_out)#
-
void bs_model_info_R(bs_model_rng **model, char const **info_out)#
-
void bs_param_names_R(bs_model_rng **model, int *include_tp, int *include_gq, char const **name_out)#
-
void bs_param_unc_names_R(bs_model_rng **model, char const **name_out)#
-
void bs_param_num_R(bs_model_rng **model, int *include_tp, int *include_gq, int *num_out)#
-
void bs_param_unc_num_R(bs_model_rng **model, int *num_out)#
-
void bs_param_constrain_R(bs_model_rng **model, int *include_tp, int *include_gq, const double *theta_unc, double *theta, int *return_code)#
-
void bs_param_unconstrain_R(bs_model_rng **model, const double *theta, double *theta_unc, int *return_code)#
-
void bs_param_unconstrain_json_R(bs_model_rng **model, char const **json, double *theta_unc, int *return_code)#
-
void bs_log_density_R(bs_model_rng **model, int *propto, int *jacobian, const double *theta, double *val, int *return_code)#
-
void bs_log_density_gradient_R(bs_model_rng **model, int *propto, int *jacobian, const double *theta, double *val, double *grad, int *return_code)#
-
void bs_log_density_hessian_R(bs_model_rng **model, int *propto, int *jacobian, const double *theta, double *val, double *grad, double *hess, int *return_code)#