Exemple #1
0
void KinsolSolver::initialize(ComputeSystemFunction pComputeSystem,
                              double *pParameters, int pSize, void *pUserData)
{
    if (mSolver)
        // The solver has already been initialised, so reset things...

        reset();

    // Initialise the ODE solver itself

    OpenCOR::CoreSolver::CoreNlaSolver::initialize(pComputeSystem, pParameters, pSize);

    // Create some vectors

    mParametersVector = N_VMake_Serial(pSize, pParameters);
    mOnesVector = N_VNew_Serial(pSize);

    N_VConst(1.0, mOnesVector);

    // Create the KINSOL solver

    mSolver = KINCreate();

    // Use our own error handler

    KINSetErrHandlerFn(mSolver, errorHandler, this);

    // Initialise the KINSOL solver

    KINInit(mSolver, systemFunction, mParametersVector);

    // Set some user data

    mUserData = new KinsolSolverUserData(pUserData, pComputeSystem);

    KINSetUserData(mSolver, mUserData);

    // Set the linear solver

    KINDense(mSolver, pSize);
}
int nlsKinsolAllocate(int size, NONLINEAR_SYSTEM_DATA *nlsData, int linearSolverMethod)
{
  int i, flag, printLevel;

  NLS_KINSOL_DATA *kinsolData = (NLS_KINSOL_DATA*) malloc(sizeof(NLS_KINSOL_DATA));

  /* allocate system data */
  nlsData->solverData = (void*)kinsolData;

  kinsolData->size = size;
  kinsolData->linearSolverMethod = linearSolverMethod;
  kinsolData->solved = 0;

  kinsolData->fnormtol  = sqrt(newtonFTol);     /* function tolerance */
  kinsolData->scsteptol = sqrt(newtonXTol);     /* step tolerance */

  kinsolData->initialGuess = N_VNew_Serial(size);
  kinsolData->xScale = N_VNew_Serial(size);
  kinsolData->fScale = N_VNew_Serial(size);
  kinsolData->fRes = N_VNew_Serial(size);

  kinsolData->kinsolMemory = KINCreate();

  /* setup user defined functions */
  KINSetErrHandlerFn(kinsolData->kinsolMemory, nlsKinsolErrorPrint, kinsolData);
  KINSetInfoHandlerFn(kinsolData->kinsolMemory, nlsKinsolInfoPrint, kinsolData);
  KINSetUserData(kinsolData->kinsolMemory, (void*)&(kinsolData->userData));
  flag = KINInit(kinsolData->kinsolMemory, nlsKinsolResiduals, kinsolData->initialGuess);
  if (checkReturnFlag(flag)){
    errorStreamPrint(LOG_STDOUT, 0, "##KINSOL## Something goes wrong while initialize KINSOL solver!");
  }

  /* Specify linear solver and/or corresponding jacobian function*/
  if (kinsolData->linearSolverMethod == 3)
  {
    if(nlsData->isPatternAvailable)
    {
      kinsolData->nnz = nlsData->sparsePattern.numberOfNoneZeros;
      flag = KINKLU(kinsolData->kinsolMemory, size, kinsolData->nnz);
      if (checkReturnFlag(flag)){
        errorStreamPrint(LOG_STDOUT, 0, "##KINSOL## Something goes wrong while initialize KINSOL solver!");
      }
      flag = KINSlsSetSparseJacFn(kinsolData->kinsolMemory, nlsSparseJac);
      if (checkReturnFlag(flag)){
        errorStreamPrint(LOG_STDOUT, 0, "##KINSOL## Something goes wrong while initialize KINSOL Sparse Solver!");
      }
    }
    else
    {
      flag = KINDense(kinsolData->kinsolMemory, size);
      if (checkReturnFlag(flag)){
        errorStreamPrint(LOG_STDOUT, 0, "##KINSOL## Something goes wrong while initialize KINSOL solver!");
      }
    }
  }
  else if (kinsolData->linearSolverMethod == 1)
  {
    flag = KINDense(kinsolData->kinsolMemory, size);
    if (checkReturnFlag(flag)){
      errorStreamPrint(LOG_STDOUT, 0, "##KINSOL## Something goes wrong while initialize KINSOL solver!");
    }
  }
  else if (kinsolData->linearSolverMethod == 2)
  {
    flag = KINDense(kinsolData->kinsolMemory, size);
    if (checkReturnFlag(flag)){
      errorStreamPrint(LOG_STDOUT, 0, "##KINSOL## Something goes wrong while initialize KINSOL solver!");
    }
    flag = KINDlsSetDenseJacFn(kinsolData->kinsolMemory, nlsDenseJac);
    if (checkReturnFlag(flag)){
      errorStreamPrint(LOG_STDOUT, 0, "##KINSOL## Something goes wrong while initialize KINSOL Sparse Solver!");
    }
  }

  /* configuration */
  nlsKinsolConfigSetup(kinsolData);

  /* debug print level of kinsol */
  if (ACTIVE_STREAM(LOG_NLS))
    printLevel = 1;
  else if (ACTIVE_STREAM(LOG_NLS_V))
    printLevel = 3;
  else
    printLevel = 0;
  KINSetPrintLevel(kinsolData->kinsolMemory, printLevel);

  return 0;
}
Exemple #3
0
void FKIN_DENSE(int *neq, int *ier)
{
  *ier = KINDense(KIN_kinmem, *neq);
  KIN_ls = KIN_LS_DENSE;
}
int main()
{
  UserData data;
  realtype fnormtol, scsteptol;
  N_Vector u1, u2, u, s, c;
  int glstr, mset, flag;
  void *kmem;

  u1 = u2 = u = NULL;
  s = c = NULL;
  kmem = NULL;
  data = NULL;

  /* User data */

  data = (UserData)malloc(sizeof *data);
  data->lb[0] = PT25;       data->ub[0] = ONE;
  data->lb[1] = ONEPT5;     data->ub[1] = TWO*PI;

  /* Create serial vectors of length NEQ */
  u1 = N_VNew_Serial(NEQ);
  if (check_flag((void *)u1, "N_VNew_Serial", 0)) return(1);

  u2 = N_VNew_Serial(NEQ);
  if (check_flag((void *)u2, "N_VNew_Serial", 0)) return(1);

  u = N_VNew_Serial(NEQ);
  if (check_flag((void *)u, "N_VNew_Serial", 0)) return(1);

  s = N_VNew_Serial(NEQ);
  if (check_flag((void *)s, "N_VNew_Serial", 0)) return(1);

  c = N_VNew_Serial(NEQ);
  if (check_flag((void *)c, "N_VNew_Serial", 0)) return(1);

  SetInitialGuess1(u1,data);
  SetInitialGuess2(u2,data);

  N_VConst_Serial(ONE,s); /* no scaling */

  Ith(c,1) =  ZERO;   /* no constraint on x1 */
  Ith(c,2) =  ZERO;   /* no constraint on x2 */
  Ith(c,3) =  ONE;    /* l1 = x1 - x1_min >= 0 */
  Ith(c,4) = -ONE;    /* L1 = x1 - x1_max <= 0 */
  Ith(c,5) =  ONE;    /* l2 = x2 - x2_min >= 0 */
  Ith(c,6) = -ONE;    /* L2 = x2 - x22_min <= 0 */
  
  fnormtol=FTOL; scsteptol=STOL;


  kmem = KINCreate();
  if (check_flag((void *)kmem, "KINCreate", 0)) return(1);

  flag = KINSetUserData(kmem, data);
  if (check_flag(&flag, "KINSetUserData", 1)) return(1);
  flag = KINSetConstraints(kmem, c);
  if (check_flag(&flag, "KINSetConstraints", 1)) return(1);
  flag = KINSetFuncNormTol(kmem, fnormtol);
  if (check_flag(&flag, "KINSetFuncNormTol", 1)) return(1);
  flag = KINSetScaledStepTol(kmem, scsteptol);
  if (check_flag(&flag, "KINSetScaledStepTol", 1)) return(1);

  flag = KINInit(kmem, func, u);
  if (check_flag(&flag, "KINInit", 1)) return(1);

  /* Call KINDense to specify the linear solver */

  flag = KINDense(kmem, NEQ);
  if (check_flag(&flag, "KINDense", 1)) return(1);

  /* Print out the problem size, solution parameters, initial guess. */
  PrintHeader(fnormtol, scsteptol);

  /* --------------------------- */

  printf("\n------------------------------------------\n");
  printf("\nInitial guess on lower bounds\n");
  printf("  [x1,x2] = ");
  PrintOutput(u1);

  N_VScale_Serial(ONE,u1,u);
  glstr = KIN_NONE;
  mset = 1;
  SolveIt(kmem, u, s, glstr, mset);

  /* --------------------------- */

  N_VScale_Serial(ONE,u1,u);
  glstr = KIN_LINESEARCH;
  mset = 1;
  SolveIt(kmem, u, s, glstr, mset);

  /* --------------------------- */

  N_VScale_Serial(ONE,u1,u);
  glstr = KIN_NONE;
  mset = 0;
  SolveIt(kmem, u, s, glstr, mset);

  /* --------------------------- */

  N_VScale_Serial(ONE,u1,u);
  glstr = KIN_LINESEARCH;
  mset = 0;
  SolveIt(kmem, u, s, glstr, mset);



  /* --------------------------- */

  printf("\n------------------------------------------\n");
  printf("\nInitial guess in middle of feasible region\n");
  printf("  [x1,x2] = ");
  PrintOutput(u2);

  N_VScale_Serial(ONE,u2,u);
  glstr = KIN_NONE;
  mset = 1;
  SolveIt(kmem, u, s, glstr, mset);

  /* --------------------------- */

  N_VScale_Serial(ONE,u2,u);
  glstr = KIN_LINESEARCH;
  mset = 1;
  SolveIt(kmem, u, s, glstr, mset);

  /* --------------------------- */

  N_VScale_Serial(ONE,u2,u);
  glstr = KIN_NONE;
  mset = 0;
  SolveIt(kmem, u, s, glstr, mset);

  /* --------------------------- */

  N_VScale_Serial(ONE,u2,u);
  glstr = KIN_LINESEARCH;
  mset = 0;
  SolveIt(kmem, u, s, glstr, mset);




  /* Free memory */

  N_VDestroy_Serial(u);
  N_VDestroy_Serial(s);
  N_VDestroy_Serial(c);
  KINFree(&kmem);
  free(data);

  return(0);
}
Exemple #5
0
int main()
{
  realtype fnormtol, scsteptol;
  N_Vector y, scale, constraints;
  int mset, flag, i;
  void *kmem;

  y = scale = constraints = NULL;
  kmem = NULL;

  printf("\nRobot Kinematics Example\n");
  printf("8 variables; -1 <= x_i <= 1\n");
  printf("KINSOL problem size: 8 + 2*8 = 24 \n\n");

  /* Create vectors for solution, scales, and constraints */

  y = N_VNew_Serial(NEQ);
  if (check_flag((void *)y, "N_VNew_Serial", 0)) return(1);

  scale = N_VNew_Serial(NEQ);
  if (check_flag((void *)scale, "N_VNew_Serial", 0)) return(1);

  constraints = N_VNew_Serial(NEQ);
  if (check_flag((void *)constraints, "N_VNew_Serial", 0)) return(1);

  /* Initialize and allocate memory for KINSOL */

  kmem = KINCreate();
  if (check_flag((void *)kmem, "KINCreate", 0)) return(1);

  flag = KINInit(kmem, func, y); /* y passed as a template */
  if (check_flag(&flag, "KINInit", 1)) return(1);

  /* Set optional inputs */

  N_VConst_Serial(ZERO,constraints);
  for (i = NVAR+1; i <= NEQ; i++) Ith(constraints, i) = ONE;
  
  flag = KINSetConstraints(kmem, constraints);
  if (check_flag(&flag, "KINSetConstraints", 1)) return(1);

  fnormtol  = FTOL; 
  flag = KINSetFuncNormTol(kmem, fnormtol);
  if (check_flag(&flag, "KINSetFuncNormTol", 1)) return(1);

  scsteptol = STOL;
  flag = KINSetScaledStepTol(kmem, scsteptol);
  if (check_flag(&flag, "KINSetScaledStepTol", 1)) return(1);

  /* Attach dense linear solver */

  flag = KINDense(kmem, NEQ);
  if (check_flag(&flag, "KINDense", 1)) return(1);

  flag = KINDlsSetDenseJacFn(kmem, jac);
  if (check_flag(&flag, "KINDlsSetDenseJacFn", 1)) return(1);

  /* Indicate exact Newton */

  mset = 1;
  flag = KINSetMaxSetupCalls(kmem, mset);
  if (check_flag(&flag, "KINSetMaxSetupCalls", 1)) return(1);

  /* Initial guess */

  N_VConst_Serial(ONE, y);
  for(i = 1; i <= NVAR; i++) Ith(y,i) = SQRT(TWO)/TWO;

  printf("Initial guess:\n");
  PrintOutput(y);

  /* Call KINSol to solve problem */

  N_VConst_Serial(ONE,scale);
  flag = KINSol(kmem,           /* KINSol memory block */
                y,              /* initial guess on input; solution vector */
                KIN_LINESEARCH, /* global stragegy choice */
                scale,          /* scaling vector, for the variable cc */
                scale);         /* scaling vector for function values fval */
  if (check_flag(&flag, "KINSol", 1)) return(1);

  printf("\nComputed solution:\n");
  PrintOutput(y);

  /* Print final statistics and free memory */  

  PrintFinalStats(kmem);

  N_VDestroy_Serial(y);
  N_VDestroy_Serial(scale);
  N_VDestroy_Serial(constraints);
  KINFree(&kmem);

  return(0);
}
  /*! \fn kinsol_initialization
   *
   *  \param [ref] [data]
   *  \param [in]  [initData]
   *  \param [ref] [initialResiduals]
   *
   *  \author lochel
   */
  int kinsol_initialization(DATA* data, INIT_DATA* initData, int useScaling)
  {
    long i, indz;
    KINSOL_DATA* kdata = NULL;
    double fnormtol  = 1.e-9;     /* function tolerance */
    double scsteptol = 1.e-9;     /* step tolerance */

    long int nni, nfe, nje, nfeD;

    N_Vector z = NULL;
    N_Vector sVars = NULL;
    N_Vector sEqns = NULL;
    N_Vector c = NULL;

    int glstr = KIN_NONE;   /* globalization strategy applied to the Newton method. It must be one of KIN_NONE or KIN_LINESEARCH */
    long int mset = 1;      /* maximum number of nonlinear iterations without a call to the preconditioner setup function. Pass 0 to indicate the default [10]. */
    void *kmem = NULL;
    int error_code = -1;

    ASSERT(data->modelData.nInitResiduals == initData->nz, "The number of initial equations are not consistent with the number of unfixed variables. Select a different initialization.");

    do /* Try it first with KIN_NONE. If that fails, try it with KIN_LINESEARCH. */
    {
      if(mset == 1 && glstr == KIN_NONE)
        DEBUG_INFO(LOG_INIT, "using exact Newton");
      else if(mset == 1)
        DEBUG_INFO(LOG_INIT, "using exact Newton with line search");
      else if(glstr == KIN_NONE)
        DEBUG_INFO(LOG_INIT, "using modified Newton");
      else
        DEBUG_INFO(LOG_INIT, "using modified Newton with line search");

      DEBUG_INFO_AL1(LOG_INIT, "| mset               = %10ld", mset);
      DEBUG_INFO_AL1(LOG_INIT, "| function tolerance = %10.6g", fnormtol);
      DEBUG_INFO_AL1(LOG_INIT, "| step tolerance     = %10.6g", scsteptol);

      kdata = (KINSOL_DATA*)malloc(sizeof(KINSOL_DATA));
      ASSERT(kdata, "out of memory");

      kdata->initData = initData;
      kdata->data = data;

      z = N_VNew_Serial(3*initData->nz);
      ASSERT(z, "out of memory");

      sVars = N_VNew_Serial(3*initData->nz);
      ASSERT(sVars, "out of memory");

      sEqns = N_VNew_Serial(3*initData->nz);
      ASSERT(sEqns, "out of memory");

      c = N_VNew_Serial(3*initData->nz);
      ASSERT(c, "out of memory");

      /* initial guess */
      for(i=0; i<initData->nz; ++i)
      {
        NV_Ith_S(z, i) = initData->start[i];
        NV_Ith_S(z, initData->nInitResiduals+2*i+0) = NV_Ith_S(z, i) - initData->min[i];
        NV_Ith_S(z, initData->nInitResiduals+2*i+1) = NV_Ith_S(z, i) - initData->max[i];
      }

      kdata->useScaling=useScaling;
      if(useScaling)
      {
        computeInitialResidualScalingCoefficients(data, initData);
        for(i=0; i<initData->nz; ++i)
        {
          NV_Ith_S(sVars, i) = 1.0 / (initData->nominal[i] == 0.0 ? 1.0 : initData->nominal[i]);
          NV_Ith_S(sVars, initData->nInitResiduals+2*i+0) = NV_Ith_S(sVars, i);
          NV_Ith_S(sVars, initData->nInitResiduals+2*i+1) = NV_Ith_S(sVars, i);

          NV_Ith_S(sEqns, i) = 1.0 / (initData->residualScalingCoefficients[i] == 0.0 ? 1.0 : initData->residualScalingCoefficients[i]);
          NV_Ith_S(sEqns, initData->nInitResiduals+2*i+0) = NV_Ith_S(sEqns, i);
          NV_Ith_S(sEqns, initData->nInitResiduals+2*i+1) = NV_Ith_S(sEqns, i);
        }
      }
      else
      {
        N_VConst_Serial(1.0, sVars);        /* no scaling */
        N_VConst_Serial(1.0, sEqns);        /* no scaling */
      }

      for(i=0; i<initData->nz; ++i)
      {
        NV_Ith_S(c, i) =  0.0;        /* no constraint on z[i] */
        NV_Ith_S(c, initData->nInitResiduals+2*i+0) = 1.0;
        NV_Ith_S(c, initData->nInitResiduals+2*i+1) = -1.0;
      }

      kmem = KINCreate();
      ASSERT(kmem, "out of memory");

      KINSetErrHandlerFn(kmem, kinsol_errorHandler, NULL);
      KINSetUserData(kmem, kdata);
      KINSetConstraints(kmem, c);
      KINSetFuncNormTol(kmem, fnormtol);
      KINSetScaledStepTol(kmem, scsteptol);
      KINInit(kmem, kinsol_residuals, z);

      /* Call KINDense to specify the linear solver */
      KINDense(kmem, 3*initData->nz);

      KINSetMaxSetupCalls(kmem, mset);
      /*KINSetNumMaxIters(kmem, 2000);*/

      globalInitialResiduals = initData->initialResiduals;

      error_code = KINSol(kmem,           /* KINSol memory block */
             z,              /* initial guess on input; solution vector */
             glstr,          /* global stragegy choice */
             sVars,          /* scaling vector, for the variable cc */
             sEqns);         /* scaling vector for function values fval */

      globalInitialResiduals = NULL;

      KINGetNumNonlinSolvIters(kmem, &nni);
      KINGetNumFuncEvals(kmem, &nfe);
      KINDlsGetNumJacEvals(kmem, &nje);
      KINDlsGetNumFuncEvals(kmem, &nfeD);

      DEBUG_INFO(LOG_INIT, "final kinsol statistics");
      DEBUG_INFO_AL1(LOG_INIT, "| KINGetNumNonlinSolvIters = %5ld", nni);
      DEBUG_INFO_AL1(LOG_INIT, "| KINGetNumFuncEvals       = %5ld", nfe);
      DEBUG_INFO_AL1(LOG_INIT, "| KINDlsGetNumJacEvals     = %5ld", nje);
      DEBUG_INFO_AL1(LOG_INIT, "| KINDlsGetNumFuncEvals    = %5ld", nfeD);

      /* Free memory */
      N_VDestroy_Serial(z);
      N_VDestroy_Serial(sVars);
      N_VDestroy_Serial(sEqns);
      N_VDestroy_Serial(c);
      KINFree(&kmem);
      free(kdata);

      if(error_code < 0)
        glstr++;  /* try next globalization strategy */
    }while(error_code < 0 && glstr <= KIN_LINESEARCH);

    /* debug output */
    indz = 0;
    DEBUG_INFO(LOG_INIT, "solution");
    for(i=0; i<data->modelData.nStates; ++i)
      if(data->modelData.realVarsData[i].attribute.fixed==0)
        DEBUG_INFO_AL2(LOG_INIT, "| %s = %g", initData->name[indz++], data->localData[0]->realVars[i]);

    for(i=0; i<data->modelData.nParametersReal; ++i)
      if(data->modelData.realParameterData[i].attribute.fixed == 0)
        DEBUG_INFO_AL2(LOG_INIT, "| %s = %g", initData->name[indz++], data->simulationInfo.realParameter[i]);

    if(error_code < 0)
      THROW("kinsol failed. see last warning. use [-lv LOG_INIT] for more output.");

    return 0;
  }
static void KIM_Malloc(int nlhs, mxArray *plhs[], int nrhs, const mxArray *prhs[])
{
    int status;

    mxArray *mx_in[3], *mx_out[2];

    int mxiter, msbset, msbsetsub, etachoice, mxnbcf;
    double eta, egamma, ealpha, mxnewtstep, relfunc, fnormtol, scsteptol;
    booleantype verbose, noInitSetup, noMinEps;

    double *constraints;
    N_Vector NVconstraints;

    int ptype;
    int mudq, mldq, mupper, mlower;
    int maxl, maxrs;
    double dqrely;

    /*
     * -----------------------------
     * Find out the vector type and
     * then pass it to the vector
     * library.
     * -----------------------------
     */

    /* Send vec_type and mx_comm */

    InitVectors();

    /*
     * -----------------------------
     * Extract stuff from arguments:
     * - SYS function
     * - problem dimension
     * - solver options
     * - user data
     * -----------------------------
     */

    /* Matlab user-provided function */

    mxDestroyArray(mx_SYSfct);
    mx_SYSfct = mxDuplicateArray(prhs[0]);

    /* problem dimension */

    N = (int) mxGetScalar(prhs[1]);

    /* Solver Options -- optional argument */

    status = get_SolverOptions(prhs[2],
                               &verbose,
                               &mxiter, &msbset, &msbsetsub, &etachoice, &mxnbcf,
                               &eta, &egamma, &ealpha, &mxnewtstep,
                               &relfunc, &fnormtol, &scsteptol,
                               &constraints,
                               &noInitSetup, &noMinEps);


    /* User data -- optional argument */

    mxDestroyArray(mx_data);
    mx_data = mxDuplicateArray(prhs[3]);

    /*
     * -----------------------------------------------------
     * Set solution vector (used as a template to KINMAlloc)
     * -----------------------------------------------------
     */

    y = NewVector(N);

    /*
     * ----------------------------------------
     * Create kinsol object and allocate memory
     * ----------------------------------------
     */

    kin_mem = KINCreate();

    /* attach error handler function */
    status = KINSetErrHandlerFn(kin_mem, mtlb_KINErrHandler, NULL);

    if (verbose) {
        status = KINSetPrintLevel(kin_mem,3);
        /* attach info handler function */
        status = KINSetInfoHandlerFn(kin_mem, mtlb_KINInfoHandler, NULL);
        /* initialize the output window */
        mx_in[0] = mxCreateScalarDouble(0);
        mx_in[1] = mxCreateScalarDouble(0); /* ignored */
        mx_in[2] = mxCreateScalarDouble(0); /* ignored */
        mexCallMATLAB(1,mx_out,3,mx_in,"kim_info");
        fig_handle = (int)*mxGetPr(mx_out[0]);
    }

    /* Call KINMalloc */

    status = KINMalloc(kin_mem, mtlb_KINSys, y);

    /* Redirect output */
    status = KINSetErrFile(kin_mem, stdout);

    /* Optional inputs */

    status = KINSetNumMaxIters(kin_mem,mxiter);
    status = KINSetNoInitSetup(kin_mem,noInitSetup);
    status = KINSetNoMinEps(kin_mem,noMinEps);
    status = KINSetMaxSetupCalls(kin_mem,msbset);
    status = KINSetMaxSubSetupCalls(kin_mem,msbsetsub);
    status = KINSetMaxBetaFails(kin_mem,mxnbcf);
    status = KINSetEtaForm(kin_mem,etachoice);
    status = KINSetEtaConstValue(kin_mem,eta);
    status = KINSetEtaParams(kin_mem,egamma,ealpha);
    status = KINSetMaxNewtonStep(kin_mem,mxnewtstep);
    status = KINSetRelErrFunc(kin_mem,relfunc);
    status = KINSetFuncNormTol(kin_mem,fnormtol);
    status = KINSetScaledStepTol(kin_mem,scsteptol);
    if (constraints != NULL) {
        NVconstraints = N_VCloneEmpty(y);
        N_VSetArrayPointer(constraints, NVconstraints);
        status = KINSetConstraints(kin_mem,NVconstraints);
        N_VDestroy(NVconstraints);
    }

    status = get_LinSolvOptions(prhs[2],
                                &mupper, &mlower,
                                &mudq, &mldq, &dqrely,
                                &ptype, &maxrs, &maxl);

    switch (ls) {

    case LS_NONE:

        mexErrMsgTxt("KINMalloc:: no linear solver specified.");

        break;

    case LS_DENSE:

        status = KINDense(kin_mem, N);
        if (!mxIsEmpty(mx_JACfct))
            status = KINDenseSetJacFn(kin_mem, mtlb_KINDenseJac, NULL);

        break;

    case LS_BAND:

        status = KINBand(kin_mem, N, mupper, mlower);
        if (!mxIsEmpty(mx_JACfct))
            status = KINBandSetJacFn(kin_mem, mtlb_KINBandJac, NULL);

        break;

    case LS_SPGMR:

        switch(pm) {
        case PM_NONE:
            status = KINSpgmr(kin_mem, maxl);
            if (!mxIsEmpty(mx_PSOLfct)) {
                if (!mxIsEmpty(mx_PSETfct))
                    status = KINSpilsSetPreconditioner(kin_mem, mtlb_KINSpilsPset, mtlb_KINSpilsPsol, NULL);
                else
                    status = KINSpilsSetPreconditioner(kin_mem, NULL, mtlb_KINSpilsPsol, NULL);
            }
            break;
        case PM_BBDPRE:
            if (!mxIsEmpty(mx_GCOMfct))
                bbd_data = KINBBDPrecAlloc(kin_mem, N, mudq, mldq, mupper, mlower, dqrely, mtlb_KINGloc, mtlb_KINGcom);
            else
                bbd_data = KINBBDPrecAlloc(kin_mem, N, mudq, mldq, mupper, mlower, dqrely, mtlb_KINGloc, NULL);
            status = KINBBDSpgmr(kin_mem, maxl, bbd_data);
            break;
        }

        status = KINSpilsSetMaxRestarts(kin_mem, maxrs);

        if (!mxIsEmpty(mx_JACfct))
            status = KINSpilsSetJacTimesVecFn(kin_mem, mtlb_KINSpilsJac, NULL);

        break;

    case LS_SPBCG:

        switch(pm) {
        case PM_NONE:
            status = KINSpbcg(kin_mem, maxl);
            if (!mxIsEmpty(mx_PSOLfct)) {
                if (!mxIsEmpty(mx_PSETfct))
                    status = KINSpilsSetPreconditioner(kin_mem, mtlb_KINSpilsPset, mtlb_KINSpilsPsol, NULL);
                else
                    status = KINSpilsSetPreconditioner(kin_mem, NULL, mtlb_KINSpilsPsol, NULL);
            }
            break;
        case PM_BBDPRE:
            if (!mxIsEmpty(mx_GCOMfct))
                bbd_data = KINBBDPrecAlloc(kin_mem, N, mudq, mldq, mupper, mlower, dqrely, mtlb_KINGloc, mtlb_KINGcom);
            else
                bbd_data = KINBBDPrecAlloc(kin_mem, N, mudq, mldq, mupper, mlower, dqrely, mtlb_KINGloc, NULL);
            status = KINBBDSpbcg(kin_mem, maxl, bbd_data);
            break;
        }

        if (!mxIsEmpty(mx_JACfct))
            status = KINSpilsSetJacTimesVecFn(kin_mem, mtlb_KINSpilsJac, NULL);

        break;

    case LS_SPTFQMR:

        switch(pm) {
        case PM_NONE:
            status = KINSptfqmr(kin_mem, maxl);
            if (!mxIsEmpty(mx_PSOLfct)) {
                if (!mxIsEmpty(mx_PSETfct))
                    status = KINSpilsSetPreconditioner(kin_mem, mtlb_KINSpilsPset, mtlb_KINSpilsPsol, NULL);
                else
                    status = KINSpilsSetPreconditioner(kin_mem, NULL, mtlb_KINSpilsPsol, NULL);
            }
            break;
        case PM_BBDPRE:
            if (!mxIsEmpty(mx_GCOMfct))
                bbd_data = KINBBDPrecAlloc(kin_mem, N, mudq, mldq, mupper, mlower, dqrely, mtlb_KINGloc, mtlb_KINGcom);
            else
                bbd_data = KINBBDPrecAlloc(kin_mem, N, mudq, mldq, mupper, mlower, dqrely, mtlb_KINGloc, NULL);
            status = KINBBDSptfqmr(kin_mem, maxl, bbd_data);
            break;
        }

        if (!mxIsEmpty(mx_JACfct))
            status = KINSpilsSetJacTimesVecFn(kin_mem, mtlb_KINSpilsJac, NULL);

        break;

    }

    return;
}
Exemple #8
0
bool kinsol_solve(void) {
	double fnormtol;
	N_Vector y, scale, constraints;
	int mset, flag;
	void *kmem;

	y = scale = constraints = NULL;
	kmem = NULL;

	/* Create vectors for solution, scales, and constraints */

	y = N_VMake_Serial(nvariables, variables);
	if (y == NULL) goto cleanup;

	scale = N_VNew_Serial(nvariables);
	if (scale == NULL) goto cleanup;

	constraints = N_VNew_Serial(nvariables);
	if (constraints == NULL) goto cleanup;

	/* Initialize and allocate memory for KINSOL */

	kmem = KINCreate();
	if (kmem == NULL) goto cleanup;

	flag = KINInit(kmem, _f, y); 
	if (flag < 0) goto cleanup;

	// This constrains metabolites >= 0
	N_VConst_Serial(1.0, constraints);
	flag = KINSetConstraints(kmem, constraints);
	if (flag < 0) goto cleanup;

	fnormtol  = ARGS.ABSTOL; 
	flag = KINSetFuncNormTol(kmem, fnormtol);
	if (flag < 0) goto cleanup;

	//scsteptol = ARGS.RELTOL;
	//flag = KINSetScaledStepTol(kmem, scsteptol);
	//if (flag < 0) goto cleanup;

	/* Attach dense linear solver */

	flag = KINDense(kmem, nvariables);
	if (flag < 0) goto cleanup;

	flag = KINDlsSetDenseJacFn(kmem, _dfdy);
	if (flag < 0) goto cleanup;

	/* Indicate exact Newton */

	//This makes sure that the user-supplied Jacobian gets evaluated onevery step.
	mset = 1;
	flag = KINSetMaxSetupCalls(kmem, mset);
	if (flag < 0) goto cleanup;

	/* Initial guess is the intial variable concentrations */

	/* Call KINSol to solve problem */

	N_VConst_Serial(1.0,scale);
	flag = KINSol(kmem,           /* KINSol memory block */
	              y,              /* initial guess on input; solution vector */
	              KIN_NONE, /* basic Newton iteration */
	              scale,          /* scaling vector, for the variable cc */
	              scale);         /* scaling vector for function values fval */
	if((ARGS.MODE != ABC) and (ARGS.MODE != ABCSIM) and (ARGS.MODE != Prior))
		OUT_nums();
	if (flag < 0) goto cleanup;


	N_VDestroy_Serial(y);
	N_VDestroy_Serial(scale);
	N_VDestroy_Serial(constraints);
	KINFree(&kmem);

	return true;
	
	 cleanup:
	error("Could not integrate!\n");
	if (y!=NULL) N_VDestroy_Serial(y);
	if (scale!=NULL) N_VDestroy_Serial(scale);
	if (constraints!=NULL) N_VDestroy_Serial(constraints);
	if (kmem!=NULL) KINFree(&kmem);
	return false;
}