1 /* 2 Provides a PETSc interface to SUNDIALS/CVODE solver. 3 The interface to PVODE (old version of CVODE) was originally contributed 4 by Liyang Xu. It has been redone by Hong Zhang and Dinesh Kaushik. 5 6 Reference: sundials-2.4.0/examples/cvode/parallel/cvDiurnal_kry_p.c 7 */ 8 #include <../src/ts/impls/implicit/sundials/sundials.h> /*I "petscts.h" I*/ 9 10 /* 11 TSPrecond_Sundials - function that we provide to SUNDIALS to 12 evaluate the preconditioner. 13 */ 14 PetscErrorCode TSPrecond_Sundials(realtype tn, N_Vector y, N_Vector fy, booleantype jok, booleantype *jcurPtr, realtype _gamma, void *P_data, N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3) 15 { 16 TS ts = (TS)P_data; 17 TS_Sundials *cvode = (TS_Sundials *)ts->data; 18 PC pc; 19 Mat J, P; 20 Vec yy = cvode->w1, yydot = cvode->ydot; 21 PetscReal gm = (PetscReal)_gamma; 22 PetscScalar *y_data; 23 24 PetscFunctionBegin; 25 PetscCall(TSGetIJacobian(ts, &J, &P, NULL, NULL)); 26 y_data = (PetscScalar *)N_VGetArrayPointer(y); 27 PetscCall(VecPlaceArray(yy, y_data)); 28 PetscCall(VecZeroEntries(yydot)); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */ 29 /* compute the shifted Jacobian (1/gm)*I + Jrest */ 30 PetscCall(TSComputeIJacobian(ts, ts->ptime, yy, yydot, 1 / gm, J, P, PETSC_FALSE)); 31 PetscCall(VecResetArray(yy)); 32 PetscCall(MatScale(P, gm)); /* turn into I-gm*Jrest, J is not used by Sundials */ 33 *jcurPtr = TRUE; 34 PetscCall(TSSundialsGetPC(ts, &pc)); 35 PetscCall(PCSetOperators(pc, J, P)); 36 PetscFunctionReturn(0); 37 } 38 39 /* 40 TSPSolve_Sundials - routine that we provide to Sundials that applies the preconditioner. 41 */ 42 PetscErrorCode TSPSolve_Sundials(realtype tn, N_Vector y, N_Vector fy, N_Vector r, N_Vector z, realtype _gamma, realtype delta, int lr, void *P_data, N_Vector vtemp) 43 { 44 TS ts = (TS)P_data; 45 TS_Sundials *cvode = (TS_Sundials *)ts->data; 46 PC pc; 47 Vec rr = cvode->w1, zz = cvode->w2; 48 PetscScalar *r_data, *z_data; 49 50 PetscFunctionBegin; 51 /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/ 52 r_data = (PetscScalar *)N_VGetArrayPointer(r); 53 z_data = (PetscScalar *)N_VGetArrayPointer(z); 54 PetscCall(VecPlaceArray(rr, r_data)); 55 PetscCall(VecPlaceArray(zz, z_data)); 56 57 /* Solve the Px=r and put the result in zz */ 58 PetscCall(TSSundialsGetPC(ts, &pc)); 59 PetscCall(PCApply(pc, rr, zz)); 60 PetscCall(VecResetArray(rr)); 61 PetscCall(VecResetArray(zz)); 62 PetscFunctionReturn(0); 63 } 64 65 /* 66 TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side. 67 */ 68 int TSFunction_Sundials(realtype t, N_Vector y, N_Vector ydot, void *ctx) 69 { 70 TS ts = (TS)ctx; 71 DM dm; 72 DMTS tsdm; 73 TSIFunction ifunction; 74 MPI_Comm comm; 75 TS_Sundials *cvode = (TS_Sundials *)ts->data; 76 Vec yy = cvode->w1, yyd = cvode->w2, yydot = cvode->ydot; 77 PetscScalar *y_data, *ydot_data; 78 79 PetscFunctionBegin; 80 PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 81 /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/ 82 y_data = (PetscScalar *)N_VGetArrayPointer(y); 83 ydot_data = (PetscScalar *)N_VGetArrayPointer(ydot); 84 PetscCallAbort(comm, VecPlaceArray(yy, y_data)); 85 PetscCallAbort(comm, VecPlaceArray(yyd, ydot_data)); 86 87 /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */ 88 PetscCall(TSGetDM(ts, &dm)); 89 PetscCall(DMGetDMTS(dm, &tsdm)); 90 PetscCall(DMTSGetIFunction(dm, &ifunction, NULL)); 91 if (!ifunction) { 92 PetscCall(TSComputeRHSFunction(ts, t, yy, yyd)); 93 } else { /* If rhsfunction is also set, this computes both parts and shifts them to the right */ 94 PetscCall(VecZeroEntries(yydot)); 95 PetscCallAbort(comm, TSComputeIFunction(ts, t, yy, yydot, yyd, PETSC_FALSE)); 96 PetscCall(VecScale(yyd, -1.)); 97 } 98 PetscCallAbort(comm, VecResetArray(yy)); 99 PetscCallAbort(comm, VecResetArray(yyd)); 100 PetscFunctionReturn(0); 101 } 102 103 /* 104 TSStep_Sundials - Calls Sundials to integrate the ODE. 105 */ 106 PetscErrorCode TSStep_Sundials(TS ts) 107 { 108 TS_Sundials *cvode = (TS_Sundials *)ts->data; 109 PetscInt flag; 110 long int nits, lits, nsteps; 111 realtype t, tout; 112 PetscScalar *y_data; 113 void *mem; 114 115 PetscFunctionBegin; 116 mem = cvode->mem; 117 tout = ts->max_time; 118 PetscCall(VecGetArray(ts->vec_sol, &y_data)); 119 N_VSetArrayPointer((realtype *)y_data, cvode->y); 120 PetscCall(VecRestoreArray(ts->vec_sol, NULL)); 121 122 /* We would like to TSPreStage() and TSPostStage() 123 * before each stage solve but CVode does not appear to support this. */ 124 if (cvode->monitorstep) flag = CVode(mem, tout, cvode->y, &t, CV_ONE_STEP); 125 else flag = CVode(mem, tout, cvode->y, &t, CV_NORMAL); 126 127 if (flag) { /* display error message */ 128 switch (flag) { 129 case CV_ILL_INPUT: 130 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_ILL_INPUT"); 131 break; 132 case CV_TOO_CLOSE: 133 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_TOO_CLOSE"); 134 break; 135 case CV_TOO_MUCH_WORK: { 136 PetscReal tcur; 137 PetscCall(CVodeGetNumSteps(mem, &nsteps)); 138 PetscCall(CVodeGetCurrentTime(mem, &tcur)); 139 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_TOO_MUCH_WORK. At t=%g, nsteps %ld exceeds maxstep %" PetscInt_FMT ". Increase '-ts_max_steps <>' or modify TSSetMaxSteps()", (double)tcur, nsteps, ts->max_steps); 140 } break; 141 case CV_TOO_MUCH_ACC: 142 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_TOO_MUCH_ACC"); 143 break; 144 case CV_ERR_FAILURE: 145 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_ERR_FAILURE"); 146 break; 147 case CV_CONV_FAILURE: 148 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_CONV_FAILURE"); 149 break; 150 case CV_LINIT_FAIL: 151 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LINIT_FAIL"); 152 break; 153 case CV_LSETUP_FAIL: 154 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LSETUP_FAIL"); 155 break; 156 case CV_LSOLVE_FAIL: 157 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LSOLVE_FAIL"); 158 break; 159 case CV_RHSFUNC_FAIL: 160 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_RHSFUNC_FAIL"); 161 break; 162 case CV_FIRST_RHSFUNC_ERR: 163 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_FIRST_RHSFUNC_ERR"); 164 break; 165 case CV_REPTD_RHSFUNC_ERR: 166 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_REPTD_RHSFUNC_ERR"); 167 break; 168 case CV_UNREC_RHSFUNC_ERR: 169 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_UNREC_RHSFUNC_ERR"); 170 break; 171 case CV_RTFUNC_FAIL: 172 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_RTFUNC_FAIL"); 173 break; 174 default: 175 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, flag %d", flag); 176 } 177 } 178 179 /* log inner nonlinear and linear iterations */ 180 PetscCall(CVodeGetNumNonlinSolvIters(mem, &nits)); 181 PetscCall(CVSpilsGetNumLinIters(mem, &lits)); 182 ts->snes_its += nits; 183 ts->ksp_its = lits; 184 185 /* copy the solution from cvode->y to cvode->update and sol */ 186 PetscCall(VecPlaceArray(cvode->w1, y_data)); 187 PetscCall(VecCopy(cvode->w1, cvode->update)); 188 PetscCall(VecResetArray(cvode->w1)); 189 PetscCall(VecCopy(cvode->update, ts->vec_sol)); 190 191 ts->time_step = t - ts->ptime; 192 ts->ptime = t; 193 194 PetscCall(CVodeGetNumSteps(mem, &nsteps)); 195 if (!cvode->monitorstep) ts->steps += nsteps - 1; /* TSStep() increments the step counter by one */ 196 PetscFunctionReturn(0); 197 } 198 199 static PetscErrorCode TSInterpolate_Sundials(TS ts, PetscReal t, Vec X) 200 { 201 TS_Sundials *cvode = (TS_Sundials *)ts->data; 202 N_Vector y; 203 PetscScalar *x_data; 204 PetscInt glosize, locsize; 205 206 PetscFunctionBegin; 207 /* get the vector size */ 208 PetscCall(VecGetSize(X, &glosize)); 209 PetscCall(VecGetLocalSize(X, &locsize)); 210 PetscCall(VecGetArray(X, &x_data)); 211 212 /* Initialize N_Vec y with x_data */ 213 if (cvode->use_dense) { 214 PetscMPIInt size; 215 216 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 217 PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "TSSUNDIALS only supports a dense solve in the serial case"); 218 y = N_VMake_Serial(locsize, (realtype *)x_data); 219 } else { 220 y = N_VMake_Parallel(cvode->comm_sundials, locsize, glosize, (realtype *)x_data); 221 } 222 223 PetscCheck(y, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Interpolated y is not allocated"); 224 225 PetscCall(CVodeGetDky(cvode->mem, t, 0, y)); 226 PetscCall(VecRestoreArray(X, &x_data)); 227 PetscFunctionReturn(0); 228 } 229 230 PetscErrorCode TSReset_Sundials(TS ts) 231 { 232 TS_Sundials *cvode = (TS_Sundials *)ts->data; 233 234 PetscFunctionBegin; 235 PetscCall(VecDestroy(&cvode->update)); 236 PetscCall(VecDestroy(&cvode->ydot)); 237 PetscCall(VecDestroy(&cvode->w1)); 238 PetscCall(VecDestroy(&cvode->w2)); 239 if (cvode->mem) CVodeFree(&cvode->mem); 240 PetscFunctionReturn(0); 241 } 242 243 PetscErrorCode TSDestroy_Sundials(TS ts) 244 { 245 TS_Sundials *cvode = (TS_Sundials *)ts->data; 246 247 PetscFunctionBegin; 248 PetscCall(TSReset_Sundials(ts)); 249 PetscCallMPI(MPI_Comm_free(&(cvode->comm_sundials))); 250 PetscCall(PetscFree(ts->data)); 251 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetType_C", NULL)); 252 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxl_C", NULL)); 253 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetLinearTolerance_C", NULL)); 254 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetGramSchmidtType_C", NULL)); 255 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetTolerance_C", NULL)); 256 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMinTimeStep_C", NULL)); 257 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxTimeStep_C", NULL)); 258 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetPC_C", NULL)); 259 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetIterations_C", NULL)); 260 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsMonitorInternalSteps_C", NULL)); 261 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetUseDense_C", NULL)); 262 PetscFunctionReturn(0); 263 } 264 265 PetscErrorCode TSSetUp_Sundials(TS ts) 266 { 267 TS_Sundials *cvode = (TS_Sundials *)ts->data; 268 PetscInt glosize, locsize, i, flag; 269 PetscScalar *y_data, *parray; 270 void *mem; 271 PC pc; 272 PCType pctype; 273 PetscBool pcnone; 274 275 PetscFunctionBegin; 276 PetscCheck(ts->exact_final_time != TS_EXACTFINALTIME_MATCHSTEP, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for exact final time option 'MATCHSTEP' when using Sundials"); 277 278 /* get the vector size */ 279 PetscCall(VecGetSize(ts->vec_sol, &glosize)); 280 PetscCall(VecGetLocalSize(ts->vec_sol, &locsize)); 281 282 /* allocate the memory for N_Vec y */ 283 if (cvode->use_dense) { 284 PetscMPIInt size; 285 286 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 287 PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "TSSUNDIALS only supports a dense solve in the serial case"); 288 cvode->y = N_VNew_Serial(locsize); 289 } else { 290 cvode->y = N_VNew_Parallel(cvode->comm_sundials, locsize, glosize); 291 } 292 PetscCheck(cvode->y, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "cvode->y is not allocated"); 293 294 /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 295 PetscCall(VecGetArray(ts->vec_sol, &parray)); 296 y_data = (PetscScalar *)N_VGetArrayPointer(cvode->y); 297 for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 298 PetscCall(VecRestoreArray(ts->vec_sol, NULL)); 299 300 PetscCall(VecDuplicate(ts->vec_sol, &cvode->update)); 301 PetscCall(VecDuplicate(ts->vec_sol, &cvode->ydot)); 302 303 /* 304 Create work vectors for the TSPSolve_Sundials() routine. Note these are 305 allocated with zero space arrays because the actual array space is provided 306 by Sundials and set using VecPlaceArray(). 307 */ 308 PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts), 1, locsize, PETSC_DECIDE, NULL, &cvode->w1)); 309 PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts), 1, locsize, PETSC_DECIDE, NULL, &cvode->w2)); 310 311 /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 312 mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 313 PetscCheck(mem, PETSC_COMM_SELF, PETSC_ERR_MEM, "CVodeCreate() fails"); 314 cvode->mem = mem; 315 316 /* Set the pointer to user-defined data */ 317 flag = CVodeSetUserData(mem, ts); 318 PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeSetUserData() fails"); 319 320 /* Sundials may choose to use a smaller initial step, but will never use a larger step. */ 321 flag = CVodeSetInitStep(mem, (realtype)ts->time_step); 322 PetscCheck(!flag, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetInitStep() failed"); 323 if (cvode->mindt > 0) { 324 flag = CVodeSetMinStep(mem, (realtype)cvode->mindt); 325 if (flag) { 326 PetscCheck(flag != CV_MEM_NULL, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMinStep() failed, cvode_mem pointer is NULL"); 327 PetscCheck(flag != CV_ILL_INPUT, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size"); 328 SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMinStep() failed"); 329 } 330 } 331 if (cvode->maxdt > 0) { 332 flag = CVodeSetMaxStep(mem, (realtype)cvode->maxdt); 333 PetscCheck(!flag, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMaxStep() failed"); 334 } 335 336 /* Call CVodeInit to initialize the integrator memory and specify the 337 * user's right hand side function in u'=f(t,u), the initial time T0, and 338 * the initial dependent variable vector cvode->y */ 339 flag = CVodeInit(mem, TSFunction_Sundials, ts->ptime, cvode->y); 340 PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeInit() fails, flag %d", flag); 341 342 /* specifies scalar relative and absolute tolerances */ 343 flag = CVodeSStolerances(mem, cvode->reltol, cvode->abstol); 344 PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeSStolerances() fails, flag %d", flag); 345 346 /* Specify max order of BDF / ADAMS method */ 347 if (cvode->maxord != PETSC_DEFAULT) { 348 flag = CVodeSetMaxOrd(mem, cvode->maxord); 349 PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeSetMaxOrd() fails, flag %d", flag); 350 } 351 352 /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */ 353 flag = CVodeSetMaxNumSteps(mem, ts->max_steps); 354 PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeSetMaxNumSteps() fails, flag %d", flag); 355 356 if (cvode->use_dense) { 357 /* call CVDense to use a dense linear solver. */ 358 PetscMPIInt size; 359 360 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 361 PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "TSSUNDIALS only supports a dense solve in the serial case"); 362 flag = CVDense(mem, locsize); 363 PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVDense() fails, flag %d", flag); 364 } else { 365 /* call CVSpgmr to use GMRES as the linear solver. */ 366 /* setup the ode integrator with the given preconditioner */ 367 PetscCall(TSSundialsGetPC(ts, &pc)); 368 PetscCall(PCGetType(pc, &pctype)); 369 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCNONE, &pcnone)); 370 if (pcnone) { 371 flag = CVSpgmr(mem, PREC_NONE, 0); 372 PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVSpgmr() fails, flag %d", flag); 373 } else { 374 flag = CVSpgmr(mem, PREC_LEFT, cvode->maxl); 375 PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVSpgmr() fails, flag %d", flag); 376 377 /* Set preconditioner and solve routines Precond and PSolve, 378 and the pointer to the user-defined block data */ 379 flag = CVSpilsSetPreconditioner(mem, TSPrecond_Sundials, TSPSolve_Sundials); 380 PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVSpilsSetPreconditioner() fails, flag %d", flag); 381 } 382 } 383 PetscFunctionReturn(0); 384 } 385 386 /* type of CVODE linear multistep method */ 387 const char *const TSSundialsLmmTypes[] = {"", "ADAMS", "BDF", "TSSundialsLmmType", "SUNDIALS_", NULL}; 388 /* type of G-S orthogonalization used by CVODE linear solver */ 389 const char *const TSSundialsGramSchmidtTypes[] = {"", "MODIFIED", "CLASSICAL", "TSSundialsGramSchmidtType", "SUNDIALS_", NULL}; 390 391 PetscErrorCode TSSetFromOptions_Sundials(TS ts, PetscOptionItems *PetscOptionsObject) 392 { 393 TS_Sundials *cvode = (TS_Sundials *)ts->data; 394 int indx; 395 PetscBool flag; 396 PC pc; 397 398 PetscFunctionBegin; 399 PetscOptionsHeadBegin(PetscOptionsObject, "SUNDIALS ODE solver options"); 400 PetscCall(PetscOptionsEList("-ts_sundials_type", "Scheme", "TSSundialsSetType", TSSundialsLmmTypes, 3, TSSundialsLmmTypes[cvode->cvode_type], &indx, &flag)); 401 if (flag) PetscCall(TSSundialsSetType(ts, (TSSundialsLmmType)indx)); 402 PetscCall(PetscOptionsEList("-ts_sundials_gramschmidt_type", "Type of orthogonalization", "TSSundialsSetGramSchmidtType", TSSundialsGramSchmidtTypes, 3, TSSundialsGramSchmidtTypes[cvode->gtype], &indx, &flag)); 403 if (flag) PetscCall(TSSundialsSetGramSchmidtType(ts, (TSSundialsGramSchmidtType)indx)); 404 PetscCall(PetscOptionsReal("-ts_sundials_atol", "Absolute tolerance for convergence", "TSSundialsSetTolerance", cvode->abstol, &cvode->abstol, NULL)); 405 PetscCall(PetscOptionsReal("-ts_sundials_rtol", "Relative tolerance for convergence", "TSSundialsSetTolerance", cvode->reltol, &cvode->reltol, NULL)); 406 PetscCall(PetscOptionsReal("-ts_sundials_mindt", "Minimum step size", "TSSundialsSetMinTimeStep", cvode->mindt, &cvode->mindt, NULL)); 407 PetscCall(PetscOptionsReal("-ts_sundials_maxdt", "Maximum step size", "TSSundialsSetMaxTimeStep", cvode->maxdt, &cvode->maxdt, NULL)); 408 PetscCall(PetscOptionsReal("-ts_sundials_linear_tolerance", "Convergence tolerance for linear solve", "TSSundialsSetLinearTolerance", cvode->linear_tol, &cvode->linear_tol, NULL)); 409 PetscCall(PetscOptionsInt("-ts_sundials_maxord", "Max Order for BDF/Adams method", "TSSundialsSetMaxOrd", cvode->maxord, &cvode->maxord, NULL)); 410 PetscCall(PetscOptionsInt("-ts_sundials_maxl", "Max dimension of the Krylov subspace", "TSSundialsSetMaxl", cvode->maxl, &cvode->maxl, NULL)); 411 PetscCall(PetscOptionsBool("-ts_sundials_monitor_steps", "Monitor SUNDIALS internal steps", "TSSundialsMonitorInternalSteps", cvode->monitorstep, &cvode->monitorstep, NULL)); 412 PetscCall(PetscOptionsBool("-ts_sundials_use_dense", "Use dense internal solver in SUNDIALS (serial only)", "TSSundialsSetUseDense", cvode->use_dense, &cvode->use_dense, NULL)); 413 PetscOptionsHeadEnd(); 414 PetscCall(TSSundialsGetPC(ts, &pc)); 415 PetscCall(PCSetFromOptions(pc)); 416 PetscFunctionReturn(0); 417 } 418 419 PetscErrorCode TSView_Sundials(TS ts, PetscViewer viewer) 420 { 421 TS_Sundials *cvode = (TS_Sundials *)ts->data; 422 char *type; 423 char atype[] = "Adams"; 424 char btype[] = "BDF: backward differentiation formula"; 425 PetscBool iascii, isstring; 426 long int nsteps, its, nfevals, nlinsetups, nfails, itmp; 427 PetscInt qlast, qcur; 428 PetscReal hinused, hlast, hcur, tcur, tolsfac; 429 PC pc; 430 431 PetscFunctionBegin; 432 if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype; 433 else type = btype; 434 435 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 436 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 437 if (iascii) { 438 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials integrator does not use SNES!\n")); 439 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials integrator type %s\n", type)); 440 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials maxord %" PetscInt_FMT "\n", cvode->maxord)); 441 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials abs tol %g rel tol %g\n", (double)cvode->abstol, (double)cvode->reltol)); 442 if (cvode->use_dense) { 443 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials integrator using a dense linear solve\n")); 444 } else { 445 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials linear solver tolerance factor %g\n", (double)cvode->linear_tol)); 446 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials max dimension of Krylov subspace %" PetscInt_FMT "\n", cvode->maxl)); 447 if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 448 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n")); 449 } else { 450 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n")); 451 } 452 } 453 if (cvode->mindt > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials minimum time step %g\n", (double)cvode->mindt)); 454 if (cvode->maxdt > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials maximum time step %g\n", (double)cvode->maxdt)); 455 456 /* Outputs from CVODE, CVSPILS */ 457 PetscCall(CVodeGetTolScaleFactor(cvode->mem, &tolsfac)); 458 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials suggested factor for tolerance scaling %g\n", tolsfac)); 459 PetscCall(CVodeGetIntegratorStats(cvode->mem, &nsteps, &nfevals, &nlinsetups, &nfails, &qlast, &qcur, &hinused, &hlast, &hcur, &tcur)); 460 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials cumulative number of internal steps %ld\n", nsteps)); 461 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of calls to rhs function %ld\n", nfevals)); 462 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of calls to linear solver setup function %ld\n", nlinsetups)); 463 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of error test failures %ld\n", nfails)); 464 465 PetscCall(CVodeGetNonlinSolvStats(cvode->mem, &its, &nfails)); 466 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of nonlinear solver iterations %ld\n", its)); 467 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of nonlinear convergence failure %ld\n", nfails)); 468 if (!cvode->use_dense) { 469 PetscCall(CVSpilsGetNumLinIters(cvode->mem, &its)); /* its = no. of calls to TSPrecond_Sundials() */ 470 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of linear iterations %ld\n", its)); 471 PetscCall(CVSpilsGetNumConvFails(cvode->mem, &itmp)); 472 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of linear convergence failures %ld\n", itmp)); 473 474 PetscCall(TSSundialsGetPC(ts, &pc)); 475 PetscCall(PCView(pc, viewer)); 476 PetscCall(CVSpilsGetNumPrecEvals(cvode->mem, &itmp)); 477 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of preconditioner evaluations %ld\n", itmp)); 478 PetscCall(CVSpilsGetNumPrecSolves(cvode->mem, &itmp)); 479 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of preconditioner solves %ld\n", itmp)); 480 } 481 PetscCall(CVSpilsGetNumJtimesEvals(cvode->mem, &itmp)); 482 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of Jacobian-vector product evaluations %ld\n", itmp)); 483 PetscCall(CVSpilsGetNumRhsEvals(cvode->mem, &itmp)); 484 PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of rhs calls for finite diff. Jacobian-vector evals %ld\n", itmp)); 485 } else if (isstring) { 486 PetscCall(PetscViewerStringSPrintf(viewer, "Sundials type %s", type)); 487 } 488 PetscFunctionReturn(0); 489 } 490 491 /* --------------------------------------------------------------------------*/ 492 PetscErrorCode TSSundialsSetType_Sundials(TS ts, TSSundialsLmmType type) 493 { 494 TS_Sundials *cvode = (TS_Sundials *)ts->data; 495 496 PetscFunctionBegin; 497 cvode->cvode_type = type; 498 PetscFunctionReturn(0); 499 } 500 501 PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts, PetscInt maxl) 502 { 503 TS_Sundials *cvode = (TS_Sundials *)ts->data; 504 505 PetscFunctionBegin; 506 cvode->maxl = maxl; 507 PetscFunctionReturn(0); 508 } 509 510 PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts, PetscReal tol) 511 { 512 TS_Sundials *cvode = (TS_Sundials *)ts->data; 513 514 PetscFunctionBegin; 515 cvode->linear_tol = tol; 516 PetscFunctionReturn(0); 517 } 518 519 PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts, TSSundialsGramSchmidtType type) 520 { 521 TS_Sundials *cvode = (TS_Sundials *)ts->data; 522 523 PetscFunctionBegin; 524 cvode->gtype = type; 525 PetscFunctionReturn(0); 526 } 527 528 PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts, PetscReal aabs, PetscReal rel) 529 { 530 TS_Sundials *cvode = (TS_Sundials *)ts->data; 531 532 PetscFunctionBegin; 533 if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 534 if (rel != PETSC_DECIDE) cvode->reltol = rel; 535 PetscFunctionReturn(0); 536 } 537 538 PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts, PetscReal mindt) 539 { 540 TS_Sundials *cvode = (TS_Sundials *)ts->data; 541 542 PetscFunctionBegin; 543 cvode->mindt = mindt; 544 PetscFunctionReturn(0); 545 } 546 547 PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts, PetscReal maxdt) 548 { 549 TS_Sundials *cvode = (TS_Sundials *)ts->data; 550 551 PetscFunctionBegin; 552 cvode->maxdt = maxdt; 553 PetscFunctionReturn(0); 554 } 555 556 PetscErrorCode TSSundialsSetUseDense_Sundials(TS ts, PetscBool use_dense) 557 { 558 TS_Sundials *cvode = (TS_Sundials *)ts->data; 559 560 PetscFunctionBegin; 561 cvode->use_dense = use_dense; 562 PetscFunctionReturn(0); 563 } 564 565 PetscErrorCode TSSundialsGetPC_Sundials(TS ts, PC *pc) 566 { 567 SNES snes; 568 KSP ksp; 569 570 PetscFunctionBegin; 571 PetscCall(TSGetSNES(ts, &snes)); 572 PetscCall(SNESGetKSP(snes, &ksp)); 573 PetscCall(KSPGetPC(ksp, pc)); 574 PetscFunctionReturn(0); 575 } 576 577 PetscErrorCode TSSundialsGetIterations_Sundials(TS ts, int *nonlin, int *lin) 578 { 579 PetscFunctionBegin; 580 if (nonlin) *nonlin = ts->snes_its; 581 if (lin) *lin = ts->ksp_its; 582 PetscFunctionReturn(0); 583 } 584 585 PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts, PetscBool s) 586 { 587 TS_Sundials *cvode = (TS_Sundials *)ts->data; 588 589 PetscFunctionBegin; 590 cvode->monitorstep = s; 591 PetscFunctionReturn(0); 592 } 593 /* -------------------------------------------------------------------------------------------*/ 594 595 /*@C 596 TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 597 598 Not Collective 599 600 Input Parameter: 601 . ts - the time-step context 602 603 Output Parameters: 604 + nonlin - number of nonlinear iterations 605 - lin - number of linear iterations 606 607 Level: advanced 608 609 Notes: 610 These return the number since the creation of the TS object 611 612 .seealso: `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 613 `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 614 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 615 `TSSundialsSetLinearTolerance()`, `TSSundialsGetPC()`, `TSSetExactFinalTime()` 616 617 @*/ 618 PetscErrorCode TSSundialsGetIterations(TS ts, int *nonlin, int *lin) 619 { 620 PetscFunctionBegin; 621 PetscUseMethod(ts, "TSSundialsGetIterations_C", (TS, int *, int *), (ts, nonlin, lin)); 622 PetscFunctionReturn(0); 623 } 624 625 /*@ 626 TSSundialsSetType - Sets the method that Sundials will use for integration. 627 628 Logically Collective on TS 629 630 Input Parameters: 631 + ts - the time-step context 632 - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 633 634 Level: intermediate 635 636 .seealso: `TSSundialsGetIterations()`, `TSSundialsSetMaxl()`, 637 `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 638 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 639 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 640 `TSSetExactFinalTime()` 641 @*/ 642 PetscErrorCode TSSundialsSetType(TS ts, TSSundialsLmmType type) 643 { 644 PetscFunctionBegin; 645 PetscTryMethod(ts, "TSSundialsSetType_C", (TS, TSSundialsLmmType), (ts, type)); 646 PetscFunctionReturn(0); 647 } 648 649 /*@ 650 TSSundialsSetMaxord - Sets the maximum order for BDF/Adams method used by SUNDIALS. 651 652 Logically Collective on TS 653 654 Input Parameters: 655 + ts - the time-step context 656 - maxord - maximum order of BDF / Adams method 657 658 Level: advanced 659 660 .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, 661 `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 662 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 663 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 664 `TSSetExactFinalTime()` 665 666 @*/ 667 PetscErrorCode TSSundialsSetMaxord(TS ts, PetscInt maxord) 668 { 669 PetscFunctionBegin; 670 PetscValidLogicalCollectiveInt(ts, maxord, 2); 671 PetscTryMethod(ts, "TSSundialsSetMaxOrd_C", (TS, PetscInt), (ts, maxord)); 672 PetscFunctionReturn(0); 673 } 674 675 /*@ 676 TSSundialsSetMaxl - Sets the dimension of the Krylov space used by 677 GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 678 this is the maximum number of GMRES steps that will be used. 679 680 Logically Collective on TS 681 682 Input Parameters: 683 + ts - the time-step context 684 - maxl - number of direction vectors (the dimension of Krylov subspace). 685 686 Level: advanced 687 688 .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, 689 `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 690 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 691 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 692 `TSSetExactFinalTime()` 693 694 @*/ 695 PetscErrorCode TSSundialsSetMaxl(TS ts, PetscInt maxl) 696 { 697 PetscFunctionBegin; 698 PetscValidLogicalCollectiveInt(ts, maxl, 2); 699 PetscTryMethod(ts, "TSSundialsSetMaxl_C", (TS, PetscInt), (ts, maxl)); 700 PetscFunctionReturn(0); 701 } 702 703 /*@ 704 TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 705 system by SUNDIALS. 706 707 Logically Collective on TS 708 709 Input Parameters: 710 + ts - the time-step context 711 - tol - the factor by which the tolerance on the nonlinear solver is 712 multiplied to get the tolerance on the linear solver, .05 by default. 713 714 Level: advanced 715 716 .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 717 `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 718 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 719 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 720 `TSSetExactFinalTime()` 721 722 @*/ 723 PetscErrorCode TSSundialsSetLinearTolerance(TS ts, PetscReal tol) 724 { 725 PetscFunctionBegin; 726 PetscValidLogicalCollectiveReal(ts, tol, 2); 727 PetscTryMethod(ts, "TSSundialsSetLinearTolerance_C", (TS, PetscReal), (ts, tol)); 728 PetscFunctionReturn(0); 729 } 730 731 /*@ 732 TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 733 in GMRES method by SUNDIALS linear solver. 734 735 Logically Collective on TS 736 737 Input Parameters: 738 + ts - the time-step context 739 - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 740 741 Level: advanced 742 743 .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 744 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, 745 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 746 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 747 `TSSetExactFinalTime()` 748 749 @*/ 750 PetscErrorCode TSSundialsSetGramSchmidtType(TS ts, TSSundialsGramSchmidtType type) 751 { 752 PetscFunctionBegin; 753 PetscTryMethod(ts, "TSSundialsSetGramSchmidtType_C", (TS, TSSundialsGramSchmidtType), (ts, type)); 754 PetscFunctionReturn(0); 755 } 756 757 /*@ 758 TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 759 Sundials for error control. 760 761 Logically Collective on TS 762 763 Input Parameters: 764 + ts - the time-step context 765 . aabs - the absolute tolerance 766 - rel - the relative tolerance 767 768 See the Cvode/Sundials users manual for exact details on these parameters. Essentially 769 these regulate the size of the error for a SINGLE timestep. 770 771 Level: intermediate 772 773 .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetGMRESMaxl()`, 774 `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, 775 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 776 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 777 `TSSetExactFinalTime()` 778 779 @*/ 780 PetscErrorCode TSSundialsSetTolerance(TS ts, PetscReal aabs, PetscReal rel) 781 { 782 PetscFunctionBegin; 783 PetscTryMethod(ts, "TSSundialsSetTolerance_C", (TS, PetscReal, PetscReal), (ts, aabs, rel)); 784 PetscFunctionReturn(0); 785 } 786 787 /*@ 788 TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 789 790 Input Parameter: 791 . ts - the time-step context 792 793 Output Parameter: 794 . pc - the preconditioner context 795 796 Level: advanced 797 798 .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 799 `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 800 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 801 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()` 802 @*/ 803 PetscErrorCode TSSundialsGetPC(TS ts, PC *pc) 804 { 805 PetscFunctionBegin; 806 PetscUseMethod(ts, "TSSundialsGetPC_C", (TS, PC *), (ts, pc)); 807 PetscFunctionReturn(0); 808 } 809 810 /*@ 811 TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 812 813 Input Parameters: 814 + ts - the time-step context 815 - mindt - lowest time step if positive, negative to deactivate 816 817 Note: 818 Sundials will error if it is not possible to keep the estimated truncation error below 819 the tolerance set with TSSundialsSetTolerance() without going below this step size. 820 821 Level: beginner 822 823 .seealso: `TSSundialsSetType()`, `TSSundialsSetTolerance()`, 824 @*/ 825 PetscErrorCode TSSundialsSetMinTimeStep(TS ts, PetscReal mindt) 826 { 827 PetscFunctionBegin; 828 PetscTryMethod(ts, "TSSundialsSetMinTimeStep_C", (TS, PetscReal), (ts, mindt)); 829 PetscFunctionReturn(0); 830 } 831 832 /*@ 833 TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 834 835 Input Parameters: 836 + ts - the time-step context 837 - maxdt - lowest time step if positive, negative to deactivate 838 839 Level: beginner 840 841 .seealso: `TSSundialsSetType()`, `TSSundialsSetTolerance()`, 842 @*/ 843 PetscErrorCode TSSundialsSetMaxTimeStep(TS ts, PetscReal maxdt) 844 { 845 PetscFunctionBegin; 846 PetscTryMethod(ts, "TSSundialsSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt)); 847 PetscFunctionReturn(0); 848 } 849 850 /*@ 851 TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 852 853 Input Parameters: 854 + ts - the time-step context 855 - ft - PETSC_TRUE if monitor, else PETSC_FALSE 856 857 Level: beginner 858 859 .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 860 TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 861 TSSundialsGetIterations(), TSSundialsSetType(), 862 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 863 @*/ 864 PetscErrorCode TSSundialsMonitorInternalSteps(TS ts, PetscBool ft) 865 { 866 PetscFunctionBegin; 867 PetscTryMethod(ts, "TSSundialsMonitorInternalSteps_C", (TS, PetscBool), (ts, ft)); 868 PetscFunctionReturn(0); 869 } 870 871 /*@ 872 TSSundialsSetUseDense - Set a flag to use a dense linear solver in SUNDIALS (serial only) 873 874 Logically Collective 875 876 Input Parameters: 877 + ts - the time-step context 878 - use_dense - PETSC_TRUE to use the dense solver 879 880 Level: advanced 881 882 .seealso: `TSSUNDIALS` 883 884 @*/ 885 PetscErrorCode TSSundialsSetUseDense(TS ts, PetscBool use_dense) 886 { 887 PetscFunctionBegin; 888 PetscValidLogicalCollectiveInt(ts, use_dense, 2); 889 PetscTryMethod(ts, "TSSundialsSetUseDense_C", (TS, PetscBool), (ts, use_dense)); 890 PetscFunctionReturn(0); 891 } 892 893 /* -------------------------------------------------------------------------------------------*/ 894 /*MC 895 TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 896 897 Options Database: 898 + -ts_sundials_type <bdf,adams> - 899 . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 900 . -ts_sundials_atol <tol> - Absolute tolerance for convergence 901 . -ts_sundials_rtol <tol> - Relative tolerance for convergence 902 . -ts_sundials_linear_tolerance <tol> - 903 . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace 904 . -ts_sundials_monitor_steps - Monitor SUNDIALS internal steps 905 - -ts_sundials_use_dense - Use a dense linear solver within CVODE (serial only) 906 907 Notes: 908 This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply, 909 only PETSc PC options. 910 911 Level: beginner 912 913 .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, `TSSundialsSetLinearTolerance()`, 914 `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, `TSSundialsGetIterations()`, `TSSetExactFinalTime()` 915 916 M*/ 917 PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts) 918 { 919 TS_Sundials *cvode; 920 PC pc; 921 922 PetscFunctionBegin; 923 ts->ops->reset = TSReset_Sundials; 924 ts->ops->destroy = TSDestroy_Sundials; 925 ts->ops->view = TSView_Sundials; 926 ts->ops->setup = TSSetUp_Sundials; 927 ts->ops->step = TSStep_Sundials; 928 ts->ops->interpolate = TSInterpolate_Sundials; 929 ts->ops->setfromoptions = TSSetFromOptions_Sundials; 930 ts->default_adapt_type = TSADAPTNONE; 931 932 PetscCall(PetscNew(&cvode)); 933 934 ts->usessnes = PETSC_TRUE; 935 936 ts->data = (void *)cvode; 937 cvode->cvode_type = SUNDIALS_BDF; 938 cvode->gtype = SUNDIALS_CLASSICAL_GS; 939 cvode->maxl = 5; 940 cvode->maxord = PETSC_DEFAULT; 941 cvode->linear_tol = .05; 942 cvode->monitorstep = PETSC_TRUE; 943 cvode->use_dense = PETSC_FALSE; 944 945 PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)ts), &(cvode->comm_sundials))); 946 947 cvode->mindt = -1.; 948 cvode->maxdt = -1.; 949 950 /* set tolerance for Sundials */ 951 cvode->reltol = 1e-6; 952 cvode->abstol = 1e-6; 953 954 /* set PCNONE as default pctype */ 955 PetscCall(TSSundialsGetPC_Sundials(ts, &pc)); 956 PetscCall(PCSetType(pc, PCNONE)); 957 958 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetType_C", TSSundialsSetType_Sundials)); 959 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxl_C", TSSundialsSetMaxl_Sundials)); 960 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetLinearTolerance_C", TSSundialsSetLinearTolerance_Sundials)); 961 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetGramSchmidtType_C", TSSundialsSetGramSchmidtType_Sundials)); 962 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetTolerance_C", TSSundialsSetTolerance_Sundials)); 963 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMinTimeStep_C", TSSundialsSetMinTimeStep_Sundials)); 964 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxTimeStep_C", TSSundialsSetMaxTimeStep_Sundials)); 965 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetPC_C", TSSundialsGetPC_Sundials)); 966 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetIterations_C", TSSundialsGetIterations_Sundials)); 967 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsMonitorInternalSteps_C", TSSundialsMonitorInternalSteps_Sundials)); 968 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetUseDense_C", TSSundialsSetUseDense_Sundials)); 969 PetscFunctionReturn(0); 970 } 971