1 /* 2 Provides a PETSc interface to version 2.5 of SUNDIALS/CVODE solver (a very old version) 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 `TSSUNDIALS`. 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 Note: 610 These return the number since the creation of the `TS` object 611 612 .seealso: [](chapter_ts), `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 613 `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 614 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 615 `TSSundialsSetLinearTolerance()`, `TSSundialsGetPC()`, `TSSetExactFinalTime()` 616 @*/ 617 PetscErrorCode TSSundialsGetIterations(TS ts, int *nonlin, int *lin) 618 { 619 PetscFunctionBegin; 620 PetscUseMethod(ts, "TSSundialsGetIterations_C", (TS, int *, int *), (ts, nonlin, lin)); 621 PetscFunctionReturn(0); 622 } 623 624 /*@ 625 TSSundialsSetType - Sets the method that `TSSUNDIALS` will use for integration. 626 627 Logically Collective on ts 628 629 Input Parameters: 630 + ts - the time-step context 631 - type - one of `SUNDIALS_ADAMS` or `SUNDIALS_BDF` 632 633 Level: intermediate 634 635 .seealso: [](chapter_ts), `TSSundialsGetIterations()`, `TSSundialsSetMaxl()`, 636 `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 637 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 638 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 639 `TSSetExactFinalTime()` 640 @*/ 641 PetscErrorCode TSSundialsSetType(TS ts, TSSundialsLmmType type) 642 { 643 PetscFunctionBegin; 644 PetscTryMethod(ts, "TSSundialsSetType_C", (TS, TSSundialsLmmType), (ts, type)); 645 PetscFunctionReturn(0); 646 } 647 648 /*@ 649 TSSundialsSetMaxord - Sets the maximum order for BDF/Adams method used by `TSSUNDIALS`. 650 651 Logically Collective on ts 652 653 Input Parameters: 654 + ts - the time-step context 655 - maxord - maximum order of BDF / Adams method 656 657 Level: advanced 658 659 .seealso: [](chapter_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, 660 `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 661 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 662 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 663 `TSSetExactFinalTime()` 664 @*/ 665 PetscErrorCode TSSundialsSetMaxord(TS ts, PetscInt maxord) 666 { 667 PetscFunctionBegin; 668 PetscValidLogicalCollectiveInt(ts, maxord, 2); 669 PetscTryMethod(ts, "TSSundialsSetMaxOrd_C", (TS, PetscInt), (ts, maxord)); 670 PetscFunctionReturn(0); 671 } 672 673 /*@ 674 TSSundialsSetMaxl - Sets the dimension of the Krylov space used by 675 GMRES in the linear solver in `TSSUNDIALS`. `TSSUNDIALS` DOES NOT use restarted GMRES so 676 this is the maximum number of GMRES steps that will be used. 677 678 Logically Collective on ts 679 680 Input Parameters: 681 + ts - the time-step context 682 - maxl - number of direction vectors (the dimension of Krylov subspace). 683 684 Level: advanced 685 686 .seealso: [](chapter_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, 687 `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 688 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 689 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 690 `TSSetExactFinalTime()` 691 @*/ 692 PetscErrorCode TSSundialsSetMaxl(TS ts, PetscInt maxl) 693 { 694 PetscFunctionBegin; 695 PetscValidLogicalCollectiveInt(ts, maxl, 2); 696 PetscTryMethod(ts, "TSSundialsSetMaxl_C", (TS, PetscInt), (ts, maxl)); 697 PetscFunctionReturn(0); 698 } 699 700 /*@ 701 TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 702 system by `TSSUNDIALS`. 703 704 Logically Collective on ts 705 706 Input Parameters: 707 + ts - the time-step context 708 - tol - the factor by which the tolerance on the nonlinear solver is 709 multiplied to get the tolerance on the linear solver, .05 by default. 710 711 Level: advanced 712 713 .seealso: [](chapter_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 714 `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 715 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 716 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 717 `TSSetExactFinalTime()` 718 @*/ 719 PetscErrorCode TSSundialsSetLinearTolerance(TS ts, PetscReal tol) 720 { 721 PetscFunctionBegin; 722 PetscValidLogicalCollectiveReal(ts, tol, 2); 723 PetscTryMethod(ts, "TSSundialsSetLinearTolerance_C", (TS, PetscReal), (ts, tol)); 724 PetscFunctionReturn(0); 725 } 726 727 /*@ 728 TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 729 in GMRES method by `TSSUNDIALS` linear solver. 730 731 Logically Collective on ts 732 733 Input Parameters: 734 + ts - the time-step context 735 - type - either `SUNDIALS_MODIFIED_GS` or `SUNDIALS_CLASSICAL_GS` 736 737 Level: advanced 738 739 .seealso: [](chapter_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 740 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, 741 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 742 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 743 `TSSetExactFinalTime()` 744 @*/ 745 PetscErrorCode TSSundialsSetGramSchmidtType(TS ts, TSSundialsGramSchmidtType type) 746 { 747 PetscFunctionBegin; 748 PetscTryMethod(ts, "TSSundialsSetGramSchmidtType_C", (TS, TSSundialsGramSchmidtType), (ts, type)); 749 PetscFunctionReturn(0); 750 } 751 752 /*@ 753 TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 754 `TSSUNDIALS` for error control. 755 756 Logically Collective on ts 757 758 Input Parameters: 759 + ts - the time-step context 760 . aabs - the absolute tolerance 761 - rel - the relative tolerance 762 763 See the CVODE/SUNDIALS users manual for exact details on these parameters. Essentially 764 these regulate the size of the error for a SINGLE timestep. 765 766 Level: intermediate 767 768 .seealso: [](chapter_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetGMRESMaxl()`, 769 `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, 770 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 771 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 772 `TSSetExactFinalTime()` 773 @*/ 774 PetscErrorCode TSSundialsSetTolerance(TS ts, PetscReal aabs, PetscReal rel) 775 { 776 PetscFunctionBegin; 777 PetscTryMethod(ts, "TSSundialsSetTolerance_C", (TS, PetscReal, PetscReal), (ts, aabs, rel)); 778 PetscFunctionReturn(0); 779 } 780 781 /*@ 782 TSSundialsGetPC - Extract the PC context from a time-step context for `TSSUNDIALS`. 783 784 Input Parameter: 785 . ts - the time-step context 786 787 Output Parameter: 788 . pc - the preconditioner context 789 790 Level: advanced 791 792 .seealso: [](chapter_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 793 `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 794 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 795 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()` 796 @*/ 797 PetscErrorCode TSSundialsGetPC(TS ts, PC *pc) 798 { 799 PetscFunctionBegin; 800 PetscUseMethod(ts, "TSSundialsGetPC_C", (TS, PC *), (ts, pc)); 801 PetscFunctionReturn(0); 802 } 803 804 /*@ 805 TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 806 807 Input Parameters: 808 + ts - the time-step context 809 - mindt - lowest time step if positive, negative to deactivate 810 811 Note: 812 `TSSUNDIALS` will error if it is not possible to keep the estimated truncation error below 813 the tolerance set with `TSSundialsSetTolerance()` without going below this step size. 814 815 Level: beginner 816 817 .seealso: [](chapter_ts), `TSSundialsSetType()`, `TSSundialsSetTolerance()`, 818 @*/ 819 PetscErrorCode TSSundialsSetMinTimeStep(TS ts, PetscReal mindt) 820 { 821 PetscFunctionBegin; 822 PetscTryMethod(ts, "TSSundialsSetMinTimeStep_C", (TS, PetscReal), (ts, mindt)); 823 PetscFunctionReturn(0); 824 } 825 826 /*@ 827 TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 828 829 Input Parameters: 830 + ts - the time-step context 831 - maxdt - lowest time step if positive, negative to deactivate 832 833 Level: beginner 834 835 .seealso: [](chapter_ts), `TSSundialsSetType()`, `TSSundialsSetTolerance()`, 836 @*/ 837 PetscErrorCode TSSundialsSetMaxTimeStep(TS ts, PetscReal maxdt) 838 { 839 PetscFunctionBegin; 840 PetscTryMethod(ts, "TSSundialsSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt)); 841 PetscFunctionReturn(0); 842 } 843 844 /*@ 845 TSSundialsMonitorInternalSteps - Monitor `TSSUNDIALS` internal steps (Defaults to false). 846 847 Input Parameters: 848 + ts - the time-step context 849 - ft - `PETSC_TRUE` if monitor, else `PETSC_FALSE` 850 851 Level: beginner 852 853 .seealso: [](chapter_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 854 `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 855 `TSSundialsGetIterations()`, `TSSundialsSetType()`, 856 `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()` 857 @*/ 858 PetscErrorCode TSSundialsMonitorInternalSteps(TS ts, PetscBool ft) 859 { 860 PetscFunctionBegin; 861 PetscTryMethod(ts, "TSSundialsMonitorInternalSteps_C", (TS, PetscBool), (ts, ft)); 862 PetscFunctionReturn(0); 863 } 864 865 /*@ 866 TSSundialsSetUseDense - Set a flag to use a dense linear solver in `TSSUNDIALS` (serial only) 867 868 Logically Collective 869 870 Input Parameters: 871 + ts - the time-step context 872 - use_dense - `PETSC_TRUE` to use the dense solver 873 874 Level: advanced 875 876 .seealso: [](chapter_ts), `TSSUNDIALS` 877 @*/ 878 PetscErrorCode TSSundialsSetUseDense(TS ts, PetscBool use_dense) 879 { 880 PetscFunctionBegin; 881 PetscValidLogicalCollectiveInt(ts, use_dense, 2); 882 PetscTryMethod(ts, "TSSundialsSetUseDense_C", (TS, PetscBool), (ts, use_dense)); 883 PetscFunctionReturn(0); 884 } 885 886 /* -------------------------------------------------------------------------------------------*/ 887 /*MC 888 TSSUNDIALS - ODE solver using a very old version of the LLNL CVODE/SUNDIALS package, version 2.5 (now called SUNDIALS). Requires ./configure --download-sundials 889 890 Options Database Keys: 891 + -ts_sundials_type <bdf,adams> - 892 . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 893 . -ts_sundials_atol <tol> - Absolute tolerance for convergence 894 . -ts_sundials_rtol <tol> - Relative tolerance for convergence 895 . -ts_sundials_linear_tolerance <tol> - 896 . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace 897 . -ts_sundials_monitor_steps - Monitor SUNDIALS internal steps 898 - -ts_sundials_use_dense - Use a dense linear solver within CVODE (serial only) 899 900 Level: beginner 901 902 Note: 903 This uses its own nonlinear solver and Krylov method so PETSc `SNES` and `KSP` options do not apply, 904 only PETSc `PC` options. 905 906 .seealso: [](chapter_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, `TSSundialsSetLinearTolerance()`, `TSType`, 907 `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, `TSSundialsGetIterations()`, `TSSetExactFinalTime()` 908 M*/ 909 PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts) 910 { 911 TS_Sundials *cvode; 912 PC pc; 913 914 PetscFunctionBegin; 915 ts->ops->reset = TSReset_Sundials; 916 ts->ops->destroy = TSDestroy_Sundials; 917 ts->ops->view = TSView_Sundials; 918 ts->ops->setup = TSSetUp_Sundials; 919 ts->ops->step = TSStep_Sundials; 920 ts->ops->interpolate = TSInterpolate_Sundials; 921 ts->ops->setfromoptions = TSSetFromOptions_Sundials; 922 ts->default_adapt_type = TSADAPTNONE; 923 924 PetscCall(PetscNew(&cvode)); 925 926 ts->usessnes = PETSC_TRUE; 927 928 ts->data = (void *)cvode; 929 cvode->cvode_type = SUNDIALS_BDF; 930 cvode->gtype = SUNDIALS_CLASSICAL_GS; 931 cvode->maxl = 5; 932 cvode->maxord = PETSC_DEFAULT; 933 cvode->linear_tol = .05; 934 cvode->monitorstep = PETSC_TRUE; 935 cvode->use_dense = PETSC_FALSE; 936 937 PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)ts), &(cvode->comm_sundials))); 938 939 cvode->mindt = -1.; 940 cvode->maxdt = -1.; 941 942 /* set tolerance for SUNDIALS */ 943 cvode->reltol = 1e-6; 944 cvode->abstol = 1e-6; 945 946 /* set PCNONE as default pctype */ 947 PetscCall(TSSundialsGetPC_Sundials(ts, &pc)); 948 PetscCall(PCSetType(pc, PCNONE)); 949 950 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetType_C", TSSundialsSetType_Sundials)); 951 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxl_C", TSSundialsSetMaxl_Sundials)); 952 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetLinearTolerance_C", TSSundialsSetLinearTolerance_Sundials)); 953 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetGramSchmidtType_C", TSSundialsSetGramSchmidtType_Sundials)); 954 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetTolerance_C", TSSundialsSetTolerance_Sundials)); 955 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMinTimeStep_C", TSSundialsSetMinTimeStep_Sundials)); 956 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxTimeStep_C", TSSundialsSetMaxTimeStep_Sundials)); 957 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetPC_C", TSSundialsGetPC_Sundials)); 958 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetIterations_C", TSSundialsGetIterations_Sundials)); 959 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsMonitorInternalSteps_C", TSSundialsMonitorInternalSteps_Sundials)); 960 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetUseDense_C", TSSundialsSetUseDense_Sundials)); 961 PetscFunctionReturn(0); 962 } 963