static void KINForcingTerm(KINMem kin_mem, real fnormp) { real eta_max = POINT9, eta_min = POINTOHOHOHONE, eta_safe = 0.5, linmodel_norm; if (etaflag == ETACHOICE1) /* Choice 1 forcing terms. */ { /* compute the norm of f + J p , scaled L2 norm */ linmodel_norm = RSqrt(fnorm * fnorm + TWO * sfdotJp + sJpnorm * sJpnorm); /* Form the safeguarded Choice 1 eta. */ eta_safe = RPowerR(eta, ealpha); eta = ABS(fnormp - linmodel_norm) / fnorm; } if (etaflag == ETACHOICE2) /* Choice 2 forcing terms. */ { eta_safe = egamma * RPowerR(eta, ealpha); eta = egamma * RPowerR(fnormp / fnorm, ealpha); } /* Apply the safeguards. */ if (eta_safe < POINT1) eta_safe = ZERO; eta = MAX(eta, eta_safe); eta = MAX(eta, eta_min); eta = MIN(eta, eta_max); return; }
int KINSetScaledStepTol(void *kinmem, realtype scsteptol) { KINMem kin_mem; realtype uround; if (kinmem == NULL) { KINProcessError(NULL, KIN_MEM_NULL, "KINSOL", "KINSetScaledStepTol", MSG_NO_MEM); return(KIN_MEM_NULL); } kin_mem = (KINMem) kinmem; if (scsteptol < ZERO) { KINProcessError(NULL, KIN_ILL_INPUT, "KINSOL", "KINSetScaledStepTol", MSG_BAD_SCSTEPTOL); return(KIN_ILL_INPUT); } if (scsteptol == ZERO) { uround = kin_mem->kin_uround; kin_mem->kin_scsteptol = RPowerR(uround,TWOTHIRDS); } else { kin_mem->kin_scsteptol = scsteptol; } return(KIN_SUCCESS); }
int KINSetFuncNormTol(void *kinmem, realtype fnormtol) { KINMem kin_mem; realtype uround; if (kinmem == NULL) { KINProcessError(NULL, KIN_MEM_NULL, "KINSOL", "KINSetFuncNormTol", MSG_NO_MEM); return(KIN_MEM_NULL); } kin_mem = (KINMem) kinmem; if (fnormtol < ZERO) { KINProcessError(NULL, KIN_ILL_INPUT, "KINSOL", "KINSetFuncNormTol", MSG_BAD_FNORMTOL); return(KIN_ILL_INPUT); } if (fnormtol == ZERO) { uround = kin_mem->kin_uround; kin_mem->kin_fnormtol = RPowerR(uround,ONETHIRD); } else { kin_mem->kin_fnormtol = fnormtol; } return(KIN_SUCCESS); }
static int KINSolInit(void *kinmem, integer Neq, N_Vector uu, SysFn func, int globalstrategy, N_Vector uscale, N_Vector fscale, real fnormtol, real scsteptol, N_Vector constraints, boole optIn, long int iopt[], real ropt[], void *f_data) { boole ioptExists, roptExists, ioptBad, roptNeg; KINMem kin_mem; FILE *fp; /* Check for legal input parameters */ kin_mem = (KINMem)kinmem; fp = kin_mem->kin_msgfp; if (Neq <= 0) { fprintf(fp, MSG_BAD_NEQ, Neq); return(-1); } if (uu == NULL) { fprintf(fp, MSG_UU_NULL); return(-1); } if (func == NULL) { fprintf(fp, MSG_FUNC_NULL); return(-1); } if ((globalstrategy != INEXACT_NEWTON) && (globalstrategy != LINESEARCH)) { fprintf(fp, MSG_BAD_GLSTRAT, globalstrategy, INEXACT_NEWTON, LINESEARCH); return(-1); } if (uscale == NULL) { fprintf(fp, MSG_BAD_USCALE); return(-1); } if (N_VMin(uscale) <= ZERO) { fprintf(fp, MSG_USCALE_NONPOSITIVE); return(-1); } if (fscale == NULL) { fprintf(fp, MSG_BAD_FSCALE); return(-1); } if (N_VMin(fscale) <= ZERO) { fprintf(fp, MSG_FSCALE_NONPOSITIVE); return(-1); } if (fnormtol < ZERO) { fprintf(fp, MSG_BAD_FNORMTOL, fnormtol); return(-1); } if (scsteptol < ZERO) { fprintf(fp, MSG_BAD_SCSTEPTOL, scsteptol); return(-1); } if ((optIn != FALSE) && (optIn != TRUE)) { fprintf(fp, MSG_BAD_OPTIN, optIn, FALSE, TRUE); return(-1); } if ((optIn) && (iopt == NULL) && (ropt == NULL)) { fprintf(fp, MSG_BAD_OPT); return(-1); } kin_mem->kin_ioptExists = ioptExists = (iopt != NULL); kin_mem->kin_roptExists = roptExists = (ropt != NULL); /* Set the constraints flag */ if (constraints == NULL) constraintsSet = FALSE; else constraintsSet = TRUE; /* All error checking is complete at this point */ /* Copy the input parameters into KINSol state */ kin_mem->kin_uu = uu; kin_mem->kin_func = func; kin_mem->kin_f_data = f_data; kin_mem->kin_globalstrategy = globalstrategy; kin_mem->kin_fnormtol = fnormtol; kin_mem->kin_scsteptol = scsteptol; kin_mem->kin_iopt = iopt; kin_mem->kin_ropt = ropt; kin_mem->kin_constraints = constraints; kin_mem->kin_uscale = uscale; kin_mem->kin_fscale = fscale; kin_mem->kin_precondflag = FALSE; /* set to the correct state in KINSpgmr */ /* (additional readability constants are defined below/after KINMalloc) */ /* check the value of the two tolerances */ if (scsteptol <= ZERO) scsteptol = RPowerR(uround, TWOTHIRDS); if (fnormtol <= ZERO) fnormtol = RPowerR(uround, ONETHIRD); /* Handle the remaining optional inputs */ printfl = PRINTFL_DEFAULT; mxiter = MXITER_DEFAULT; if (ioptExists && optIn) { if (iopt[PRINTFL] > 0 && iopt[PRINTFL] < 4) printfl = iopt[PRINTFL]; if (iopt[MXITER] > 0) mxiter = iopt[MXITER]; ioptBad = (iopt[PRINTFL] < 0 || iopt[MXITER] < 0 || iopt[PRINTFL] > 3); if (ioptBad) { fprintf(fp, MSG_IOPT_OUT); free(kin_mem); return(-1); } } if (printfl > 0) fprintf(fp, "scsteptol used: %12.3g \n", scsteptol); if (printfl > 0) fprintf(fp, "fnormtol used: %12.3g \n", fnormtol); sqrt_relfunc = RSqrt(uround); /* calculate the default value for mxnewtstep (max Newton step) */ mxnewtstep = THOUSAND * N_VWL2Norm(uu, uscale); if (mxnewtstep < ONE) mxnewtstep = ONE; relu = RELU_DEFAULT; if (roptExists && optIn) { if (ropt[MXNEWTSTEP] > ZERO) mxnewtstep = ropt[MXNEWTSTEP]; if (ropt[RELFUNC] > ZERO) sqrt_relfunc = RSqrt(ropt[RELFUNC]); if (ropt[RELU] > ZERO) relu = ropt[RELU]; roptNeg = (ropt[MXNEWTSTEP] < ZERO || ropt[RELFUNC] < ZERO || ropt[RELU] < ZERO); if (roptNeg) { fprintf(fp, MSG_ROPT_OUT); free(kin_mem); return(-1); } } /* set up the coefficients for the eta calculation */ etaflag = ETACHOICE1; /* default choice */ if (ioptExists && optIn) { etaflag = iopt[ETACHOICE]; if (etaflag < ETACHOICE1 || etaflag > ETACONSTANT) { fprintf(fp, MSG_BAD_ETACHOICE); etaflag = ETACHOICE1; } } callForcingTerm = (etaflag != ETACONSTANT); if (etaflag == ETACHOICE1) /* this value ALWAYS used for Choice 1 */ ealpha = ONE + HALF * RSqrt(FIVE); if (etaflag == ETACHOICE2) { /* default values: */ ealpha = TWO; egamma = POINT9; } if (etaflag == ETACONSTANT) eta = POINT1; /* default value used --constant case */ else eta = HALF; /* Initial value for eta set to 0.5 for other than the ETACONSTANT option */ if (roptExists && optIn) { if (etaflag == ETACHOICE2) { ealpha = ropt[ETAALPHA]; egamma = ropt[ETAGAMMA]; if (ealpha <= ONE || ealpha > TWO) { /* alpha value out of range */ if (ealpha != ZERO) fprintf(fp, MSG_ALPHA_BAD); ealpha = 2.; } if (egamma <= ZERO || egamma > ONE) { /* gamma value out of range */ if (egamma != ZERO) fprintf(fp, MSG_GAMMA_BAD); egamma = POINT9; } } else if (etaflag == ETACONSTANT) eta = ropt[ETACONST]; if (eta <= ZERO || eta > ONE) { if (eta != ZERO) fprintf(fp, MSG_ETACONST_BAD); eta = POINT1; } } /* Initialize all the counters */ nfe = nnilpre = nni = nbcf = nbktrk = 0; /* Initialize optional output locations in iopt, ropt */ if (ioptExists) { iopt[NFE] = iopt[NNI] = 0; iopt[NBKTRK] = 0; } if (roptExists) { ropt[FNORM] = ZERO; ropt[STEPL] = ZERO; } /* check the initial guess uu against the constraints, if any, and * also see if the system func(uu) = 0 is satisfied by the initial guess uu */ if (constraintsSet) { if (!KINInitialConstraint(kin_mem)) { fprintf(fp, MSG_INITIAL_CNSTRNT); KINFreeVectors(kin_mem); free(kin_mem); return(-1); } } if (KINInitialStop(kin_mem)) return(+1); /* initialize the L2 norms of f for the linear iteration steps */ fnorm = N_VWL2Norm(fval, fscale); f1norm = HALF * fnorm * fnorm; if (printfl > 0) fprintf(fp, "KINSolInit nni= %4ld fnorm= %26.16g nfe=%6ld \n", nni, fnorm, nfe); /* Problem has been successfully initialized */ return(0); }
/* Main program */ int main() { UserData data; void *mem; N_Vector yy, yp, rr, q; int flag; realtype time, tout, incr; int nout; mem = NULL; yy = yp = NULL; /* Allocate user data. */ data = (UserData) malloc(sizeof(*data)); /* Fill user's data with the appropriate values for coefficients. */ data->k1 = RCONST(18.7); data->k2 = RCONST(0.58); data->k3 = RCONST(0.09); data->k4 = RCONST(0.42); data->K = RCONST(34.4); data->klA = RCONST(3.3); data->Ks = RCONST(115.83); data->pCO2 = RCONST(0.9); data->H = RCONST(737.0); /* Allocate N-vectors. */ yy = N_VNew_Serial(NEQ); if (check_flag((void *)yy, "N_VNew_Serial", 0)) return(1); yp = N_VNew_Serial(NEQ); if (check_flag((void *)yp, "N_VNew_Serial", 0)) return(1); /* Consistent IC for y, y'. */ #define y01 0.444 #define y02 0.00123 #define y03 0.00 #define y04 0.007 #define y05 0.0 Ith(yy,1) = RCONST(y01); Ith(yy,2) = RCONST(y02); Ith(yy,3) = RCONST(y03); Ith(yy,4) = RCONST(y04); Ith(yy,5) = RCONST(y05); Ith(yy,6) = data->Ks * RCONST(y01) * RCONST(y04); /* Get y' = - res(t0, y, 0) */ N_VConst(ZERO, yp); rr = N_VNew_Serial(NEQ); res(T0, yy, yp, rr, data); N_VScale(-ONE, rr, yp); N_VDestroy_Serial(rr); /* Create and initialize q0 for quadratures. */ q = N_VNew_Serial(1); if (check_flag((void *)q, "N_VNew_Serial", 0)) return(1); Ith(q,1) = ZERO; /* Call IDACreate and IDAInit to initialize IDA memory */ mem = IDACreate(); if(check_flag((void *)mem, "IDACreate", 0)) return(1); flag = IDAInit(mem, res, T0, yy, yp); if(check_flag(&flag, "IDAInit", 1)) return(1); /* Set tolerances. */ flag = IDASStolerances(mem, RTOL, ATOL); if(check_flag(&flag, "IDASStolerances", 1)) return(1); /* Attach user data. */ flag = IDASetUserData(mem, data); if(check_flag(&flag, "IDASetUserData", 1)) return(1); /* Attach linear solver. */ flag = IDADense(mem, NEQ); /* Initialize QUADRATURE(S). */ flag = IDAQuadInit(mem, rhsQ, q); if (check_flag(&flag, "IDAQuadInit", 1)) return(1); /* Set tolerances and error control for quadratures. */ flag = IDAQuadSStolerances(mem, RTOLQ, ATOLQ); if (check_flag(&flag, "IDAQuadSStolerances", 1)) return(1); flag = IDASetQuadErrCon(mem, TRUE); if (check_flag(&flag, "IDASetQuadErrCon", 1)) return(1); PrintHeader(RTOL, ATOL, yy); /* Print initial states */ PrintOutput(mem,0.0,yy); tout = T1; nout = 0; incr = RPowerR(TF/T1,ONE/NF); /* FORWARD run. */ while (1) { flag = IDASolve(mem, tout, &time, yy, yp, IDA_NORMAL); if (check_flag(&flag, "IDASolve", 1)) return(1); PrintOutput(mem, time, yy); nout++; tout *= incr; if (nout>NF) break; } flag = IDAGetQuad(mem, &time, q); if (check_flag(&flag, "IDAGetQuad", 1)) return(1); printf("\n--------------------------------------------------------\n"); printf("G: %24.16f \n",Ith(q,1)); printf("--------------------------------------------------------\n\n"); PrintFinalStats(mem); IDAFree(&mem); N_VDestroy_Serial(yy); N_VDestroy_Serial(yp); N_VDestroy_Serial(q); return(0); }
/* * cpLapackDenseProjSetup does the setup operations for the dense * linear solver. * It calls the Jacobian evaluation routine and, depending on ftype, * it performs various factorizations. */ static int cpLapackDenseProjSetup(CPodeMem cp_mem, N_Vector y, N_Vector cy, N_Vector c_tmp1, N_Vector c_tmp2, N_Vector s_tmp1) { int ier; CPDlsProjMem cpdlsP_mem; realtype *col_i, rim1, ri; int i, j, nd, one = 1; int retval; realtype coef_1 = ONE, coef_0 = ZERO; cpdlsP_mem = (CPDlsProjMem) lmemP; nd = ny-nc; /* Call Jacobian function (G will contain the Jacobian transposed) */ retval = djacP(nc, ny, tn, y, cy, G, JP_data, c_tmp1, c_tmp2); njeP++; if (retval < 0) { cpProcessError(cp_mem, CPDIRECT_JACFUNC_UNRECVR, "CPLAPACK", "cpLapackDenseProjSetup", MSGD_JACFUNC_FAILED); return(-1); } else if (retval > 0) { return(1); } /* Save Jacobian before factorization for possible use by lmultP */ dcopy_f77(&(G->ldata), G->data, &one, savedG->data, &one); /* Factorize G, depending on ftype */ switch (ftype) { case CPDIRECT_LU: /* * LU factorization of G^T with partial pivoting * P*G^T = | U1^T | * L^T * | U2^T | * After factorization, P is encoded in pivotsP and * G^T is overwritten with U1 (nc by nc unit upper triangular), * U2 ( nc by ny-nc rectangular), and L (nc by nc lower triangular). * Return ier if factorization failed. */ dgetrf_f77(&ny, &nc, G->data, &ny, pivotsP, &ier); if (ier != 0) return(ier); /* * Build S = U1^{-1} * U2 (in place, S overwrites U2) * For each row j of G, j = nc,...,ny-1, perform * a backward substitution (row version). * After this step, G^T contains U1, S, and L. */ for (j=nc; j<ny; j++) dtrsv_f77("L", "T", "U", &nc, G->data, &ny, (G->data + j), &ny, 1, 1, 1); /* * Build K = D1 + S^T * D2 * S * S^T is stored in g_mat[nc...ny-1][0...nc] * Compute and store only the lower triangular part of K. */ if (pnorm == CP_PROJ_L2NORM) { dsyrk_f77("L", "N", &nd, &nc, &coef_1, (G->data + nc), &ny, &coef_0, K->data, &nd, 1, 1); LapackDenseAddI(K); } else { cplLUcomputeKD(cp_mem, s_tmp1); } /* * Perform Cholesky decomposition of K: K = C*C^T * After factorization, the lower triangular part of K contains C. * Return ier if factorization failed. */ dpotrf_f77("L", &nd, K->data, &nd, &ier, 1); if (ier != 0) return(ier); break; case CPDIRECT_QR: /* * QR factorization of G^T: G^T = Q*R * After factorization, the upper triangular part of G^T * contains the matrix R. The lower trapezoidal part of * G^T, together with the array beta, encodes the orthonormal * columns of Q as elementary reflectors. */ dgeqrf_f77(&ny, &nc, G->data, &ny, beta, wrk, &len_wrk, &ier); if (ier != 0) return(ier); /* If projecting in WRMS norm */ if (pnorm == CP_PROJ_ERRNORM) { /* Build K = Q^T * D^(-1) * Q */ cplQRcomputeKD(cp_mem, s_tmp1); /* Perform Cholesky decomposition of K */ dpotrf_f77("L", &nc, K->data, &nc, &ier, 1); if (ier != 0) return(ier); } break; case CPDIRECT_QRP: /* * QR with pivoting factorization of G^T: G^T * P = Q * R. * After factorization, the upper triangular part of G^T * contains the matrix R. The lower trapezoidal part of * G^T, together with the array beta, encodes the orthonormal * columns of Q as elementary reflectors. * The pivots are stored in 'pivotsP'. */ for (i=0; i<nc; i++) pivotsP[i] = 0; dgeqp3_f77(&ny, &nc, G->data, &ny, pivotsP, beta, wrk, &len_wrk, &ier); if (ier != 0) return(ier); /* * Determine the number of independent constraints. * After the QR factorization, the diagonal elements of R should * be in decreasing order of their absolute values. */ rim1 = ABS(G->data[0]); for (i=1, nr=1; i<nc; i++, nr++) { col_i = G->cols[i]; ri = ABS(col_i[i]); if (ri < 100*uround) break; if (ri/rim1 < RPowerR(uround, THREE/FOUR)) break; } /* If projecting in WRMS norm */ if (pnorm == CP_PROJ_ERRNORM) { /* Build K = Q^T * D^(-1) * Q */ cplQRcomputeKD(cp_mem, s_tmp1); /* Perform Cholesky decomposition of K */ dpotrf_f77("L", &nc, K->data, &nc, &ier, 1); if (ier != 0) return(ier); } break; case CPDIRECT_SC: /* * Build K = G*D^(-1)*G^T * G^T is stored in g_mat[0...ny-1][0...nc] * Compute and store only the lower triangular part of K. */ if (pnorm == CP_PROJ_L2NORM) { dsyrk_f77("L", "T", &nc, &ny, &coef_1, G->data, &ny, &coef_0, K->data, &nc, 1, 1); } else { cplSCcomputeKD(cp_mem, s_tmp1); } /* * Perform Cholesky decomposition of K: K = C*C^T * After factorization, the lower triangular part of K contains C. * Return 1 if factorization failed. */ dpotrf_f77("L", &nc, K->data, &nc, &ier, 1); if (ier != 0) return(ier); break; } return(0); }