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