int main(){ //bind together a Poisson and a Normal; //make a draw producing a 2-element vector apop_model *m1 = apop_model_set_parameters(apop_poisson, 3); apop_model *m2 = apop_model_set_parameters(apop_normal, -5, 1); apop_model *mm = apop_model_stack(m1, m2); int len = 1e5; gsl_rng *r = apop_rng_alloc(1); apop_data *draws = apop_data_alloc(len, 2); for (int i=0; i< len; i++){ Apop_row (draws, i, onev); apop_draw(onev->data, r, mm); assert((int)onev->data[0] == onev->data[0]); assert(onev->data[1]<0); } //The rest of the test script recovers the parameters. //First, set up a two-page data set: poisson data on p1, Normal on p2: apop_data *comeback = apop_data_alloc(); Apop_col(draws, 0,fishdraws) comeback->vector = apop_vector_copy(fishdraws); apop_data_add_page(comeback, apop_data_alloc(), "p2"); Apop_col(draws, 1, meandraws) comeback->more->vector = apop_vector_copy(meandraws); //set up the un-parameterized stacked model, including //the name at which to split the data set apop_model *estme = apop_model_stack(apop_model_copy(apop_poisson), apop_model_copy(apop_normal)); Apop_settings_add(estme, apop_stack, splitpage, "p2"); apop_model *ested = apop_estimate(comeback, *estme); //test that the parameters are as promised. apop_model *m1back = apop_settings_get(ested, apop_stack, model1); apop_model *m2back = apop_settings_get(ested, apop_stack, model2); assert(fabs(apop_data_get(m1back->parameters, .col=-1) - 3) < 1e-2); assert(fabs(apop_data_get(m2back->parameters, .col=-1) - -5) < 1e-2); assert(fabs(apop_data_get(m2back->parameters, .col=-1, .row=1) - 1) < 1e-2); }
long double apop_linear_constraint(gsl_vector *beta, apop_data * constraint, double margin){ #else apop_varad_head(long double, apop_linear_constraint){ static threadlocal apop_data *default_constraint; gsl_vector * apop_varad_var(beta, NULL); double apop_varad_var(margin, 0); apop_data * apop_varad_var(constraint, NULL); Apop_assert(beta, "The vector to be checked is NULL."); if (!constraint){ if (default_constraint && beta->size != default_constraint->vector->size){ apop_data_free(default_constraint); default_constraint = NULL; } if (!default_constraint){ default_constraint = apop_data_alloc(0,beta->size, beta->size); default_constraint->vector = gsl_vector_calloc(beta->size); gsl_matrix_set_identity(default_constraint->matrix); } constraint = default_constraint; } return apop_linear_constraint_base(beta, constraint, margin); } long double apop_linear_constraint_base(gsl_vector *beta, apop_data * constraint, double margin){ #endif static threadlocal gsl_vector *closest_pt = NULL; static threadlocal gsl_vector *candidate = NULL; static threadlocal gsl_vector *fix = NULL; int constraint_ct = constraint->matrix->size1; int bindlist[constraint_ct]; int i, bound = 0; /* For added efficiency, keep a scratch vector or two on hand. */ if (closest_pt==NULL || closest_pt->size != constraint->matrix->size2){ closest_pt = gsl_vector_calloc(beta->size); candidate = gsl_vector_alloc(beta->size); fix = gsl_vector_alloc(beta->size); closest_pt->data[0] = GSL_NEGINF; } /* Do any constraints bind?*/ memset(bindlist, 0, sizeof(int)*constraint_ct); for (i=0; i< constraint_ct; i++){ Apop_row_v(constraint, i, c); bound += bindlist[i] = binds(beta, apop_data_get(constraint, i, -1), c, margin); } if (!bound) return 0; //All constraints met. gsl_vector *base_beta = apop_vector_copy(beta); /* With only one constraint, it's easy. */ if (constraint->vector->size==1){ Apop_row_v(constraint, 0, c); find_nearest_point(base_beta, constraint->vector->data[0], c, beta); goto add_margin; } /* Finally, multiple constraints, at least one binding. For each surface, pick a candidate point. Check whether the point meets the other constraints. if not, translate to a new point that works. [Do this by maintaining a pseudopoint that translates by the necessary amount.] Once you have a candidate point, compare its distance to the current favorite; keep the best. */ for (i=0; i< constraint_ct; i++) if (bindlist[i]){ get_candiate(base_beta, constraint, i, candidate, margin); if(apop_vector_distance(base_beta, candidate) < apop_vector_distance(base_beta, closest_pt)) gsl_vector_memcpy(closest_pt, candidate); } gsl_vector_memcpy(beta, closest_pt); add_margin: for (i=0; i< constraint_ct; i++){ if(bindlist[i]){ Apop_row_v(constraint, i, c); gsl_vector_memcpy(fix, c); gsl_vector_scale(fix, magnitude(fix)); gsl_vector_scale(fix, margin); gsl_vector_add(beta, fix); } } long double out = apop_vector_distance(base_beta, beta); gsl_vector_free(base_beta); return out; }