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