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