![]() |
|
Set some of the parameters of a model to fixed values.
| apop_model* apop_model_fix_params | ( | apop_model * | model_in | ) |
Produce a model based on another model, but with some of the parameters fixed at a given value.
As well as a pointer to the model whose parameters are to be fixed, I need two sets of data.
For the fixed parameters, I need to know the values at which they'll be fixed. Send me an apop_data set that has the same shape as the parameters of your model; at the positions of the fixed parameters, give the values to which they will be fixed. For the free parameters, I (mostly) don't care what value they have. This set of parameters can be either set as the model_in->parameters element, or as an argument to the model. [If you give me both, I will use the one explicitly sent in rather than the one attached to the input model.]
I also need to know which parameters to fix, which requires a mask that I can hold over the paramete set. Again, the mask is an apop_data set of the same size and shape as your data. Where there is a nonzero marker, I will fix the parameter.
You again have two options for giving me this information. You can use the parameter matrix as the mask: just set parameters to be left free to a nonzero value (including GSL_NAN). Or, you can explicitly send in a mask, with ones at params to be fixed and zero elsewhere. Again, I will try the explicit mask first.
You also need to input the base model, which I will copy (along with its settings groups) to form the new model.
The output is an apop_model that can be estimated, Bayesian updated, et cetera.
estimate method always uses an MLE, and it never calls the base model's estimate method.estimate method. Otherwise, I'll set my own.Here is a sample program. It produces a few thousand draws from a multivariate normal distribution, and then tries to recover the means given a var/covar matrix fixed at the correct variance.
#include <apop.h> int main(){ gsl_rng *r = apop_rng_alloc(10); size_t i, ct = 1e4; apop_data *d = apop_data_alloc(0,ct,2); double draws[2]; apop_data *params = apop_data_alloc(2,2,2); apop_data_fill(params, 8, 1, 0.5, 2, 0.5, 1); apop_model *pvm = apop_model_copy(apop_multivariate_normal); pvm->parameters = apop_data_copy(params); for(i=0; i< ct; i++){ apop_draw(draws, r, pvm); apop_data_set(d, i, 0, draws[0]); apop_data_set(d, i, 1, draws[1]); } gsl_vector_set_all(pvm->parameters->vector, GSL_NAN); apop_model *mep1 = apop_model_fix_params(pvm); apop_model *e1 = apop_estimate(d, *mep1); printf("original params: "); apop_vector_show(params->vector); printf("estimated params: "); apop_vector_show(e1->parameters->vector); return 0; }
| model_in | The base model |