ode_solver*	ode_solver_alloc(ode_model* model){
  
  ode_solver* solver = (ode_solver*) malloc( sizeof(ode_solver) );				/* alloc */
  if (solver == 0){
    /* TODO: write a proper error handler */
    fprintf(stderr,"malloc failed to allocate memory for ode_solver\n");
    return 0;
  }
  
  solver->cvode_mem = CVodeCreate(CV_BDF,CV_NEWTON);								/* alloc */
  if( solver->cvode_mem == 0){
    /* TODO: write a proper error handler */
    fprintf(stderr,"CVodeCreate failed to allocate memory in ode_solver for cvode_mem.\n");
    
    free(solver);
    return 0;
  }
  
  int P = ode_model_getP(model);
  solver->params = (double*) malloc( sizeof(double) * P );						/* alloc */
  if( solver->params == 0 ){
    /* TODO: write a proper error handler */
    fprintf(stderr,"malloc failed to allocate memory in ode_solver for params.\n");
    CVodeFree(solver->cvode_mem);
    free(solver);
    return 0;
  }
  ode_model_get_default_params(model, solver->params, P);	
  
  int N = ode_model_getN(model);
  solver->odeModel = model;
  solver->y = N_VNewEmpty_Serial(N);												/* alloc */
  NV_DATA_S(solver->y) = solver->odeModel->v;
  solver->yS = 0;
  
  return solver;
}
Beispiel #2
0
/* Initializes MPI,
 * loads defaults, 
 *       command line arguments,
 *       hdf5 data,
 *       ode model from shared library @code dlopen@
 * allocates kernel, 
 *           ode model parameters
 *           MPI communivcation buffers
 * calls MCMC routines
 * finalizes and frees (most) structs
 */
int/*always returns success*/
main(int argc,/*count*/ char* argv[])/*array of strings*/ {
  int i=0;
  int warm_up=0; // sets the number of burn in points at command line
  char lib_name[BUFSZ];
  ode_model_parameters omp[1];
  omp->size=(problem_size*) malloc(sizeof(problem_size));

  char global_sample_filename_stem[BUFSZ]="Sample.h5"; // filename basis
  char rank_sample_file[BUFSZ]; // filename for sample output
  char resume_filename[BUFSZ]="resume.h5";
  double seed = 1;
  double gamma= 2;
  double t0=-1;
  int sampling_action=SMPL_FRESH;
  
  int start_from_prior=no;
  int sensitivity_approximation=no;

  main_options cnf_options=get_default_options(global_sample_filename_stem, lib_name);
  
  MPI_Init(&argc,&argv);
  int rank,R;
  MPI_Comm_size(MPI_COMM_WORLD,&R);
  MPI_Comm_rank(MPI_COMM_WORLD,&rank);
  char *h5file=NULL;

  gsl_set_error_handler_off();
  
  /* process command line arguments
   */  
  for (i=0;i<argc;i++){
    if (strcmp(argv[i],"-p")==0 || strcmp(argv[i],"--prior-start")==0) {
      start_from_prior=1;
    } else if (strcmp(argv[i],"-d")==0 || strcmp(argv[i],"--hdf5")==0) {
      h5file=argv[i+1];
    } else if (strcmp(argv[i],"-t")==0 || strcmp(argv[i],"--init-at-t")==0) {
      t0=strtod(argv[i+1],NULL);
      //printf("[main] t0=%f\n",t0);
    } else if (strcmp(argv[i],"-w")==0 || strcmp(argv[i],"--warm-up")==0) warm_up=strtol(argv[i+1],NULL,10);
    else if (strcmp(argv[i],"--resume")==0 || strcmp(argv[i],"-r")==0) sampling_action=SMPL_RESUME;
    else if (strcmp(argv[i],"--sens-approx")==0) sensitivity_approximation=1;
    else if (strcmp(argv[i],"-l")==0) strcpy(cnf_options.library_file,argv[i+1]);
    //    else if (strcmp(argv[i],"-n")==0) Tuning=0;
    else if (strcmp(argv[i],"-s")==0) cnf_options.sample_size=strtol(argv[i+1],NULL,0);
    else if (strcmp(argv[i],"-o")==0) strncpy(cnf_options.output_file,argv[i+1],BUFSZ);
    else if (strcmp(argv[i],"-a")==0) cnf_options.target_acceptance=strtod(argv[i+1],NULL);
    else if (strcmp(argv[i],"-i")==0 || strcmp(argv[i],"--initial-step-size")==0) cnf_options.initial_stepsize=strtod(argv[i+1],NULL);
    else if (strcmp(argv[i],"-m")==0 || strcmp(argv[i],"--initial-step-size-rank-multiplier")==0) cnf_options.initial_stepsize_rank_factor=strtod(argv[i+1],NULL);

    else if (strcmp(argv[i],"-g")==0) gamma=strtod(argv[i+1],NULL);
    else if (strcmp(argv[i],"--abs-tol")==0) cnf_options.abs_tol=strtod(argv[i+1],NULL);
    else if (strcmp(argv[i],"--rel-tol")==0) cnf_options.rel_tol=strtod(argv[i+1],NULL);
    else if (strcmp(argv[i],"--seed")==0) seed=strtod(argv[i+1],NULL);
    else if (strcmp(argv[i],"-h")==0 || strcmp(argv[i],"--help")==0) {
      print_help();
      MPI_Abort(MPI_COMM_WORLD,0);
    }
  }
  
  seed=seed*137+13*rank;

  /* load Data from hdf5 file
   */
  if (h5file){
    printf("# [main] (rank %i) reading hdf5 file, loading data.\n",rank);
    fflush(stdout);
    read_data(h5file,omp);
    fflush(stdout);
  } else {
    fprintf(stderr,"# [main] (rank %i) no data provided (-d option), exiting.\n",rank);
    MPI_Abort(MPI_COMM_WORLD,-1);
  }
    
  /* load model from shared library
   */
  ode_model *odeModel = ode_model_loadFromFile(lib_name);  /* alloc */
  if (!odeModel) {
    fprintf(stderr, "# [main] (rank %i) Library %s could not be loaded.\n",rank,lib_name);
    exit(1);
  } else printf( "# [main] (rank %i) Library %s loaded.\n",rank, lib_name);
  
  /* construct an output file from rank, library name, and user
   * supplied string.
   */
  char *dot;
  char *lib_base;
  lib_base=basename(lib_name);
  dot=strchr(lib_base,'.');
  dot[0]='\0';
  sprintf(resume_filename,"%s_resume_%02i.h5",lib_base,rank);
  sprintf(rank_sample_file,"mcmc_rank_%02i_of_%i_%s_%s",rank,R,lib_base,basename(cnf_options.output_file));
  cnf_options.output_file=rank_sample_file;
  cnf_options.resume_file=resume_filename;
  
  /* allocate a solver for each experiment for possible parallelization
   */
  ode_solver **solver;
  int c,C=omp->size->C;
  int c_success=0;
  solver=malloc(sizeof(ode_solver*)*C);
  for (c=0;c<C;c++){
    solver[c]=ode_solver_alloc(odeModel);
    if (solver[c]) c_success++;
  }
  if (c_success==C) {
    printf("# [main] Solver[0:%i] for «%s» created.\n",C,lib_base);
  } else {
    fprintf(stderr, "# [main] Solvers for «%s» could not be created.\n",lib_base);
    ode_model_free(odeModel);
    MPI_Abort(MPI_COMM_WORLD,-1);
  }

  /* sensitivity analysis is not feasible for large models. So, it can
   *  be turned off.
   */
  if (sensitivity_approximation){
    //printf("# [main] experimental: Sensitivity approximation activated.\n");
    for (c=0;c<C;c++) ode_solver_disable_sens(solver[c]);
    /* also: make sensitivity function unavailable; that way
     * ode_model_has_sens(model) will return «FALSE»;
     */
    odeModel->vf_sens=NULL;
  }
  
  /* init solver 
   */
  realtype solver_param[3] = {cnf_options.abs_tol, cnf_options.rel_tol, 0};

  const char **x_name=ode_model_get_var_names(odeModel);
  const char **p_name=ode_model_get_param_names(odeModel);
  const char **f_name=ode_model_get_func_names(odeModel);
  
  /* local variables for parameters and inital conditions as presented
     in ode model lib: */
  int N = ode_model_getN(odeModel);
  int P = ode_model_getP(odeModel);
  int F = ode_model_getF(odeModel);

  /* save in ode model parameter struct: */
  set_number_of_state_variables(omp,N);
  set_number_of_model_parameters(omp,P);
  set_number_of_model_outputs(omp,F);

  omp->t0=t0;
  /* ode model parameter struct has pointers for sim results that need
     memory allocation: */
  ode_model_parameters_alloc(omp);
  ode_model_parameters_link(omp);
  fflush(stdout);

  /* get default parameters from the model file
   */
  double p[P];
  gsl_vector_view p_view=gsl_vector_view_array(p,P);
  ode_model_get_default_params(odeModel, p, P);
  if (rank==0)  gsl_printf("default parameters",&(p_view.vector),GSL_IS_DOUBLE | GSL_IS_VECTOR);
  omp->solver=solver;

  /* All MCMC meta-parameters (like stepsize) here are positive (to
   * make sense). Some command line arguments can override parameters
   * read from files; but, input files are processed after the command
   * line parameters. So, to check whether default parameters were
   * altered by the command line, the variable declaration defaults
   * are negative at first. Alterations to some meta-parameter p can
   * be checked by: if (cnf_options.p<0)
   * cnf_options.p=read_from_file(SOME FILE);
   */
  cnf_options.initial_stepsize=fabs(cnf_options.initial_stepsize);
  cnf_options.target_acceptance=fabs(cnf_options.target_acceptance);
  cnf_options.sample_size=fabs(cnf_options.sample_size);

  /* load default initial conditions
   */
  double y[N];
  gsl_vector_view y_view=gsl_vector_view_array(y,N);
  ode_model_get_initial_conditions(odeModel, y, N);
  
  print_experiment_information(rank,R,omp,&(y_view.vector));

  /* initialize the ODE solver with initial time t, default ODE
   * parameters p and default initial conditions of the state y; In
   * addition error tolerances are set and sensitivity initialized.
   */
  //printf("# [main] (rank %i) init ivp: t0=%g\n",rank,omp->t0);
  for (c=0;c<C;c++){
    ode_solver_init(solver[c], omp->t0, omp->E[c]->init_y->data, N, p, P);
    //printf("# [main] solver initialised.\n");    
    ode_solver_setErrTol(solver[c], solver_param[1], &solver_param[0], 1);
    if (ode_model_has_sens(odeModel)) {
      ode_solver_init_sens(solver[c], omp->E[0]->yS0->data, P, N);
    }
  }
  /* An smmala_model is a struct that contains the posterior
   * probablity density function and a pointer to its parameters and
   * pre-allocated work-memory.
   */
  smmala_model* model = smmala_model_alloc(LogPosterior, NULL, omp);
  if (model){
    printf("[main] (rank %i) smmala_model allocated.\n",rank);
  }else{
    fprintf(stderr,"[main] (rank %i) smmala_model could not be allocated.\n",rank);
    MPI_Abort(MPI_COMM_WORLD,-1);
  }
  
  /* initial parameter values; after allocating an mcmc_kernel of the
   * right dimensions we set the initial Markov chain state from
   * either the model's default parametrization p, the prior's μ, or the
   * state of a previously completed mcmc run (resume).
   */
  int D=omp->size->D;
  double init_x[D];
  double beta=assign_beta(rank,R,round(gamma));
  double tgac=cnf_options.target_acceptance;
  double m=cnf_options.initial_stepsize_rank_factor;
  double step=cnf_options.initial_stepsize;
  if (m>1.0 && rank>0) step*=gsl_pow_int(m,rank);
  pdf_normalisation_constant(omp);
  printf("[main] (rank %i) likelihood log(normalisation constant): %g\n",rank,omp->pdf_lognorm);
  mcmc_kernel* kernel = smmala_kernel_alloc(beta,D,step,model,seed,tgac);
  
  int resume_load_status;
  if (sampling_action==SMPL_RESUME){
    resume_load_status=load_resume_state(resume_filename, rank, R, kernel);
    assert(resume_load_status==EXIT_SUCCESS);
    for (i=0;i<D;i++) init_x[i]=kernel->x[i];
  } else if (start_from_prior){     
    if (rank==0) printf("# [main] setting initial mcmc vector to prior mean.\n");
    for (i=0;i<D;i++) init_x[i]=gsl_vector_get(omp->prior->mu,i);
  } else {
    if (rank==0) printf("# [main] setting mcmc initial value to log(default parameters)\n");
    for (i=0;i<D;i++) init_x[i]=gsl_sf_log(p[i]);
  }
  fflush(stdout);
  //display_prior_information(omp->prior);
  
  /* here we initialize the mcmc_kernel; this makes one test
   * evaluation of the log-posterior density function. 
   */
  
  /*
  if (rank==0){
    printf("# [main] initializing MCMC.\n");
    printf("# [main] init_x:");
    for (i=0;i<D;i++) printf(" %g ",init_x[i]);
    printf("\n");
  }
  */
  
  mcmc_init(kernel, init_x);
  /* display the results of that test evaluation
   *
   */
  if (rank==0){
    printf("# [main] rank %i init complete .\n",rank);
    display_test_evaluation_results(kernel);
    ode_solver_print_stats(solver[0], stdout);
    fflush(stdout);
    fflush(stderr);  
  }
  
  size_t SampleSize = cnf_options.sample_size;  
  
  /* in parallel tempering th echains can swap their positions;
   * this buffers the communication between chains.
   */
  void *buffer=(void *) smmala_comm_buffer_alloc(D);
  
  /* Initialization of burin in length
   */
  size_t BurnInSampleSize;
  if (warm_up==0){
    BurnInSampleSize = 7 * (int) sqrt(cnf_options.sample_size);
  } else {
    BurnInSampleSize=warm_up;
  }
  if (rank==0){
    printf("# Performing Burn-In with step-size (%g) tuning: %lu iterations\n",get_step_size(kernel),BurnInSampleSize);
    fflush(stdout);
  }
  /* Burn In: these iterations are not recorded, but are used to find
   * an acceptable step size for each temperature regime.
   */
  int mcmc_error;
  mcmc_error=burn_in_foreach(rank,R, BurnInSampleSize, omp, kernel, buffer);
  assert(mcmc_error==EXIT_SUCCESS);
  if (rank==0){
    fprintf(stdout, "\n# Burn-in complete, sampling from the posterior.\n");
  }
  /* this struct contains all necessary id's and size arrays
   * for writing sample data to an hdf5 file in chunks
   */
  hdf5block_t *h5block = h5block_init(cnf_options.output_file,
				      omp,SampleSize,
				      x_name,p_name,f_name);
  
  /* The main loop of MCMC sampling
   * these iterations are recorded and saved to an hdf5 file
   * the file is set up and identified via the h5block variable.
   */  
  mcmc_error=mcmc_foreach(rank, R, SampleSize, omp, kernel, h5block, buffer, &cnf_options);
  assert(mcmc_error==EXIT_SUCCESS);
  append_meta_properties(h5block,&seed,&BurnInSampleSize, h5file, lib_base);
  h5block_close(h5block);

  /* clear memory */
  smmala_model_free(model);
  mcmc_free(kernel);
  ode_model_parameters_free(omp);
  MPI_Finalize();
  return EXIT_SUCCESS;
}