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 219 /* allocate the memory for N_Vec y */ 220 y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 221 if (!y) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Interpolated y is not allocated"); 222 223 ierr = VecGetArray(X,&x_data);CHKERRQ(ierr); 224 N_VSetArrayPointer((realtype*)x_data,y); 225 ierr = CVodeGetDky(cvode->mem,t,0,y);CHKERRQ(ierr); 226 ierr = VecRestoreArray(X,&x_data);CHKERRQ(ierr); 227 PetscFunctionReturn(0); 228 } 229 230 PetscErrorCode TSReset_Sundials(TS ts) 231 { 232 TS_Sundials *cvode = (TS_Sundials*)ts->data; 233 PetscErrorCode ierr; 234 235 PetscFunctionBegin; 236 ierr = VecDestroy(&cvode->update);CHKERRQ(ierr); 237 ierr = VecDestroy(&cvode->ydot);CHKERRQ(ierr); 238 ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr); 239 ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr); 240 if (cvode->mem) CVodeFree(&cvode->mem); 241 PetscFunctionReturn(0); 242 } 243 244 PetscErrorCode TSDestroy_Sundials(TS ts) 245 { 246 TS_Sundials *cvode = (TS_Sundials*)ts->data; 247 PetscErrorCode ierr; 248 249 PetscFunctionBegin; 250 ierr = TSReset_Sundials(ts);CHKERRQ(ierr); 251 ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr); 252 ierr = PetscFree(ts->data);CHKERRQ(ierr); 253 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",NULL);CHKERRQ(ierr); 254 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",NULL);CHKERRQ(ierr); 255 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",NULL);CHKERRQ(ierr); 256 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",NULL);CHKERRQ(ierr); 257 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",NULL);CHKERRQ(ierr); 258 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",NULL);CHKERRQ(ierr); 259 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",NULL);CHKERRQ(ierr); 260 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",NULL);CHKERRQ(ierr); 261 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",NULL);CHKERRQ(ierr); 262 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",NULL);CHKERRQ(ierr); 263 PetscFunctionReturn(0); 264 } 265 266 PetscErrorCode TSSetUp_Sundials(TS ts) 267 { 268 TS_Sundials *cvode = (TS_Sundials*)ts->data; 269 PetscErrorCode ierr; 270 PetscInt glosize,locsize,i,flag; 271 PetscScalar *y_data,*parray; 272 void *mem; 273 PC pc; 274 PCType pctype; 275 PetscBool pcnone; 276 277 PetscFunctionBegin; 278 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"); 279 280 /* get the vector size */ 281 ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr); 282 ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr); 283 284 /* allocate the memory for N_Vec y */ 285 cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 286 if (!cvode->y) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cvode->y is not allocated"); 287 288 /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 289 ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr); 290 y_data = (PetscScalar*) N_VGetArrayPointer(cvode->y); 291 for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 292 ierr = VecRestoreArray(ts->vec_sol,NULL);CHKERRQ(ierr); 293 294 ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr); 295 ierr = VecDuplicate(ts->vec_sol,&cvode->ydot);CHKERRQ(ierr); 296 ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->update);CHKERRQ(ierr); 297 ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->ydot);CHKERRQ(ierr); 298 299 /* 300 Create work vectors for the TSPSolve_Sundials() routine. Note these are 301 allocated with zero space arrays because the actual array space is provided 302 by Sundials and set using VecPlaceArray(). 303 */ 304 ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr); 305 ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr); 306 ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w1);CHKERRQ(ierr); 307 ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w2);CHKERRQ(ierr); 308 309 /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 310 mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 311 if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails"); 312 cvode->mem = mem; 313 314 /* Set the pointer to user-defined data */ 315 flag = CVodeSetUserData(mem, ts); 316 if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails"); 317 318 /* Sundials may choose to use a smaller initial step, but will never use a larger step. */ 319 flag = CVodeSetInitStep(mem,(realtype)ts->time_step); 320 if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetInitStep() failed"); 321 if (cvode->mindt > 0) { 322 flag = CVodeSetMinStep(mem,(realtype)cvode->mindt); 323 if (flag) { 324 if (flag == CV_MEM_NULL) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL"); 325 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"); 326 else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed"); 327 } 328 } 329 if (cvode->maxdt > 0) { 330 flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt); 331 if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMaxStep() failed"); 332 } 333 334 /* Call CVodeInit to initialize the integrator memory and specify the 335 * user's right hand side function in u'=f(t,u), the inital time T0, and 336 * the initial dependent variable vector cvode->y */ 337 flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y); 338 if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag); 339 340 /* specifies scalar relative and absolute tolerances */ 341 flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol); 342 if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag); 343 344 /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */ 345 flag = CVodeSetMaxNumSteps(mem,ts->max_steps); 346 if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetMaxNumSteps() fails, flag %d",flag); 347 348 /* call CVSpgmr to use GMRES as the linear solver. */ 349 /* setup the ode integrator with the given preconditioner */ 350 ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 351 ierr = PCGetType(pc,&pctype);CHKERRQ(ierr); 352 ierr = PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr); 353 if (pcnone) { 354 flag = CVSpgmr(mem,PREC_NONE,0); 355 if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 356 } else { 357 flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl); 358 if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 359 360 /* Set preconditioner and solve routines Precond and PSolve, 361 and the pointer to the user-defined block data */ 362 flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials); 363 if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag); 364 } 365 366 flag = CVSpilsSetGSType(mem, MODIFIED_GS); 367 if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag); 368 PetscFunctionReturn(0); 369 } 370 371 /* type of CVODE linear multistep method */ 372 const char *const TSSundialsLmmTypes[] = {"","ADAMS","BDF","TSSundialsLmmType","SUNDIALS_",0}; 373 /* type of G-S orthogonalization used by CVODE linear solver */ 374 const char *const TSSundialsGramSchmidtTypes[] = {"","MODIFIED","CLASSICAL","TSSundialsGramSchmidtType","SUNDIALS_",0}; 375 376 PetscErrorCode TSSetFromOptions_Sundials(PetscOptionItems *PetscOptionsObject,TS ts) 377 { 378 TS_Sundials *cvode = (TS_Sundials*)ts->data; 379 PetscErrorCode ierr; 380 int indx; 381 PetscBool flag; 382 PC pc; 383 384 PetscFunctionBegin; 385 ierr = PetscOptionsHead(PetscOptionsObject,"SUNDIALS ODE solver options");CHKERRQ(ierr); 386 ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr); 387 if (flag) { 388 ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr); 389 } 390 ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr); 391 if (flag) { 392 ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr); 393 } 394 ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,NULL);CHKERRQ(ierr); 395 ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,NULL);CHKERRQ(ierr); 396 ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,NULL);CHKERRQ(ierr); 397 ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,NULL);CHKERRQ(ierr); 398 ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,NULL);CHKERRQ(ierr); 399 ierr = PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,NULL);CHKERRQ(ierr); 400 ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internal steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,NULL);CHKERRQ(ierr); 401 ierr = PetscOptionsTail();CHKERRQ(ierr); 402 ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 403 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 404 PetscFunctionReturn(0); 405 } 406 407 PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer) 408 { 409 TS_Sundials *cvode = (TS_Sundials*)ts->data; 410 PetscErrorCode ierr; 411 char *type; 412 char atype[] = "Adams"; 413 char btype[] = "BDF: backward differentiation formula"; 414 PetscBool iascii,isstring; 415 long int nsteps,its,nfevals,nlinsetups,nfails,itmp; 416 PetscInt qlast,qcur; 417 PetscReal hinused,hlast,hcur,tcur,tolsfac; 418 PC pc; 419 420 PetscFunctionBegin; 421 if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype; 422 else type = btype; 423 424 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 425 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 426 if (iascii) { 427 ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr); 428 ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr); 429 ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr); 430 ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr); 431 ierr = PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);CHKERRQ(ierr); 432 if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 433 ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 434 } else { 435 ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 436 } 437 if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);} 438 if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);} 439 440 /* Outputs from CVODE, CVSPILS */ 441 ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr); 442 ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr); 443 ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals, 444 &nlinsetups,&nfails,&qlast,&qcur, 445 &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr); 446 ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr); 447 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr); 448 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr); 449 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr); 450 451 ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr); 452 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr); 453 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr); 454 455 ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */ 456 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr); 457 ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr); 458 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr); 459 460 ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 461 ierr = PCView(pc,viewer);CHKERRQ(ierr); 462 ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr); 463 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr); 464 ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr); 465 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr); 466 467 ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr); 468 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr); 469 ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr); 470 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr); 471 } else if (isstring) { 472 ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr); 473 } 474 PetscFunctionReturn(0); 475 } 476 477 478 /* --------------------------------------------------------------------------*/ 479 PetscErrorCode TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type) 480 { 481 TS_Sundials *cvode = (TS_Sundials*)ts->data; 482 483 PetscFunctionBegin; 484 cvode->cvode_type = type; 485 PetscFunctionReturn(0); 486 } 487 488 PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl) 489 { 490 TS_Sundials *cvode = (TS_Sundials*)ts->data; 491 492 PetscFunctionBegin; 493 cvode->maxl = maxl; 494 PetscFunctionReturn(0); 495 } 496 497 PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts,double tol) 498 { 499 TS_Sundials *cvode = (TS_Sundials*)ts->data; 500 501 PetscFunctionBegin; 502 cvode->linear_tol = tol; 503 PetscFunctionReturn(0); 504 } 505 506 PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type) 507 { 508 TS_Sundials *cvode = (TS_Sundials*)ts->data; 509 510 PetscFunctionBegin; 511 cvode->gtype = type; 512 PetscFunctionReturn(0); 513 } 514 515 PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel) 516 { 517 TS_Sundials *cvode = (TS_Sundials*)ts->data; 518 519 PetscFunctionBegin; 520 if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 521 if (rel != PETSC_DECIDE) cvode->reltol = rel; 522 PetscFunctionReturn(0); 523 } 524 525 PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt) 526 { 527 TS_Sundials *cvode = (TS_Sundials*)ts->data; 528 529 PetscFunctionBegin; 530 cvode->mindt = mindt; 531 PetscFunctionReturn(0); 532 } 533 534 PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt) 535 { 536 TS_Sundials *cvode = (TS_Sundials*)ts->data; 537 538 PetscFunctionBegin; 539 cvode->maxdt = maxdt; 540 PetscFunctionReturn(0); 541 } 542 PetscErrorCode TSSundialsGetPC_Sundials(TS ts,PC *pc) 543 { 544 SNES snes; 545 KSP ksp; 546 PetscErrorCode ierr; 547 548 PetscFunctionBegin; 549 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 550 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 551 ierr = KSPGetPC(ksp,pc);CHKERRQ(ierr); 552 PetscFunctionReturn(0); 553 } 554 555 PetscErrorCode TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin) 556 { 557 PetscFunctionBegin; 558 if (nonlin) *nonlin = ts->snes_its; 559 if (lin) *lin = ts->ksp_its; 560 PetscFunctionReturn(0); 561 } 562 563 PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s) 564 { 565 TS_Sundials *cvode = (TS_Sundials*)ts->data; 566 567 PetscFunctionBegin; 568 cvode->monitorstep = s; 569 PetscFunctionReturn(0); 570 } 571 /* -------------------------------------------------------------------------------------------*/ 572 573 /*@C 574 TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 575 576 Not Collective 577 578 Input parameters: 579 . ts - the time-step context 580 581 Output Parameters: 582 + nonlin - number of nonlinear iterations 583 - lin - number of linear iterations 584 585 Level: advanced 586 587 Notes: 588 These return the number since the creation of the TS object 589 590 .seealso: TSSundialsSetType(), TSSundialsSetMaxl(), 591 TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 592 TSSundialsGetIterations(), TSSundialsSetType(), 593 TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime() 594 595 @*/ 596 PetscErrorCode TSSundialsGetIterations(TS ts,int *nonlin,int *lin) 597 { 598 PetscErrorCode ierr; 599 600 PetscFunctionBegin; 601 ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr); 602 PetscFunctionReturn(0); 603 } 604 605 /*@ 606 TSSundialsSetType - Sets the method that Sundials will use for integration. 607 608 Logically Collective on TS 609 610 Input parameters: 611 + ts - the time-step context 612 - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 613 614 Level: intermediate 615 616 .seealso: TSSundialsGetIterations(), TSSundialsSetMaxl(), 617 TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 618 TSSundialsGetIterations(), TSSundialsSetType(), 619 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 620 TSSetExactFinalTime() 621 @*/ 622 PetscErrorCode TSSundialsSetType(TS ts,TSSundialsLmmType type) 623 { 624 PetscErrorCode ierr; 625 626 PetscFunctionBegin; 627 ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr); 628 PetscFunctionReturn(0); 629 } 630 631 /*@ 632 TSSundialsSetMaxl - Sets the dimension of the Krylov space used by 633 GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 634 this is the maximum number of GMRES steps that will be used. 635 636 Logically Collective on TS 637 638 Input parameters: 639 + ts - the time-step context 640 - maxl - number of direction vectors (the dimension of Krylov subspace). 641 642 Level: advanced 643 644 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 645 TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 646 TSSundialsGetIterations(), TSSundialsSetType(), 647 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 648 TSSetExactFinalTime() 649 650 @*/ 651 PetscErrorCode TSSundialsSetMaxl(TS ts,PetscInt maxl) 652 { 653 PetscErrorCode ierr; 654 655 PetscFunctionBegin; 656 PetscValidLogicalCollectiveInt(ts,maxl,2); 657 ierr = PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));CHKERRQ(ierr); 658 PetscFunctionReturn(0); 659 } 660 661 /*@ 662 TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 663 system by SUNDIALS. 664 665 Logically Collective on TS 666 667 Input parameters: 668 + ts - the time-step context 669 - tol - the factor by which the tolerance on the nonlinear solver is 670 multiplied to get the tolerance on the linear solver, .05 by default. 671 672 Level: advanced 673 674 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 675 TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 676 TSSundialsGetIterations(), TSSundialsSetType(), 677 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 678 TSSetExactFinalTime() 679 680 @*/ 681 PetscErrorCode TSSundialsSetLinearTolerance(TS ts,double tol) 682 { 683 PetscErrorCode ierr; 684 685 PetscFunctionBegin; 686 PetscValidLogicalCollectiveReal(ts,tol,2); 687 ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr); 688 PetscFunctionReturn(0); 689 } 690 691 /*@ 692 TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 693 in GMRES method by SUNDIALS linear solver. 694 695 Logically Collective on TS 696 697 Input parameters: 698 + ts - the time-step context 699 - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 700 701 Level: advanced 702 703 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 704 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), 705 TSSundialsGetIterations(), TSSundialsSetType(), 706 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 707 TSSetExactFinalTime() 708 709 @*/ 710 PetscErrorCode TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type) 711 { 712 PetscErrorCode ierr; 713 714 PetscFunctionBegin; 715 ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr); 716 PetscFunctionReturn(0); 717 } 718 719 /*@ 720 TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 721 Sundials for error control. 722 723 Logically Collective on TS 724 725 Input parameters: 726 + ts - the time-step context 727 . aabs - the absolute tolerance 728 - rel - the relative tolerance 729 730 See the Cvode/Sundials users manual for exact details on these parameters. Essentially 731 these regulate the size of the error for a SINGLE timestep. 732 733 Level: intermediate 734 735 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(), 736 TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 737 TSSundialsGetIterations(), TSSundialsSetType(), 738 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 739 TSSetExactFinalTime() 740 741 @*/ 742 PetscErrorCode TSSundialsSetTolerance(TS ts,double aabs,double rel) 743 { 744 PetscErrorCode ierr; 745 746 PetscFunctionBegin; 747 ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr); 748 PetscFunctionReturn(0); 749 } 750 751 /*@ 752 TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 753 754 Input Parameter: 755 . ts - the time-step context 756 757 Output Parameter: 758 . pc - the preconditioner context 759 760 Level: advanced 761 762 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 763 TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 764 TSSundialsGetIterations(), TSSundialsSetType(), 765 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance() 766 @*/ 767 PetscErrorCode TSSundialsGetPC(TS ts,PC *pc) 768 { 769 PetscErrorCode ierr; 770 771 PetscFunctionBegin; 772 ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC*),(ts,pc));CHKERRQ(ierr); 773 PetscFunctionReturn(0); 774 } 775 776 /*@ 777 TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 778 779 Input Parameter: 780 + ts - the time-step context 781 - mindt - lowest time step if positive, negative to deactivate 782 783 Note: 784 Sundials will error if it is not possible to keep the estimated truncation error below 785 the tolerance set with TSSundialsSetTolerance() without going below this step size. 786 787 Level: beginner 788 789 .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 790 @*/ 791 PetscErrorCode TSSundialsSetMinTimeStep(TS ts,PetscReal mindt) 792 { 793 PetscErrorCode ierr; 794 795 PetscFunctionBegin; 796 ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr); 797 PetscFunctionReturn(0); 798 } 799 800 /*@ 801 TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 802 803 Input Parameter: 804 + ts - the time-step context 805 - maxdt - lowest time step if positive, negative to deactivate 806 807 Level: beginner 808 809 .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 810 @*/ 811 PetscErrorCode TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt) 812 { 813 PetscErrorCode ierr; 814 815 PetscFunctionBegin; 816 ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 817 PetscFunctionReturn(0); 818 } 819 820 /*@ 821 TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 822 823 Input Parameter: 824 + ts - the time-step context 825 - ft - PETSC_TRUE if monitor, else PETSC_FALSE 826 827 Level: beginner 828 829 .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 830 TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 831 TSSundialsGetIterations(), TSSundialsSetType(), 832 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 833 @*/ 834 PetscErrorCode TSSundialsMonitorInternalSteps(TS ts,PetscBool ft) 835 { 836 PetscErrorCode ierr; 837 838 PetscFunctionBegin; 839 ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr); 840 PetscFunctionReturn(0); 841 } 842 /* -------------------------------------------------------------------------------------------*/ 843 /*MC 844 TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 845 846 Options Database: 847 + -ts_sundials_type <bdf,adams> - 848 . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 849 . -ts_sundials_atol <tol> - Absolute tolerance for convergence 850 . -ts_sundials_rtol <tol> - Relative tolerance for convergence 851 . -ts_sundials_linear_tolerance <tol> - 852 . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace 853 - -ts_sundials_monitor_steps - Monitor SUNDIALS internal steps 854 855 856 Notes: 857 This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply, 858 only PETSc PC options. 859 860 Level: beginner 861 862 .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(), 863 TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime() 864 865 M*/ 866 PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts) 867 { 868 TS_Sundials *cvode; 869 PetscErrorCode ierr; 870 PC pc; 871 872 PetscFunctionBegin; 873 ts->ops->reset = TSReset_Sundials; 874 ts->ops->destroy = TSDestroy_Sundials; 875 ts->ops->view = TSView_Sundials; 876 ts->ops->setup = TSSetUp_Sundials; 877 ts->ops->step = TSStep_Sundials; 878 ts->ops->interpolate = TSInterpolate_Sundials; 879 ts->ops->setfromoptions = TSSetFromOptions_Sundials; 880 ts->default_adapt_type = TSADAPTNONE; 881 882 ierr = PetscNewLog(ts,&cvode);CHKERRQ(ierr); 883 884 ts->usessnes = PETSC_TRUE; 885 886 ts->data = (void*)cvode; 887 cvode->cvode_type = SUNDIALS_BDF; 888 cvode->gtype = SUNDIALS_CLASSICAL_GS; 889 cvode->maxl = 5; 890 cvode->linear_tol = .05; 891 cvode->monitorstep = PETSC_TRUE; 892 893 ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)ts),&(cvode->comm_sundials));CHKERRQ(ierr); 894 895 cvode->mindt = -1.; 896 cvode->maxdt = -1.; 897 898 /* set tolerance for Sundials */ 899 cvode->reltol = 1e-6; 900 cvode->abstol = 1e-6; 901 902 /* set PCNONE as default pctype */ 903 ierr = TSSundialsGetPC_Sundials(ts,&pc);CHKERRQ(ierr); 904 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 905 906 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",TSSundialsSetType_Sundials);CHKERRQ(ierr); 907 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",TSSundialsSetMaxl_Sundials);CHKERRQ(ierr); 908 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr); 909 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr); 910 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",TSSundialsSetTolerance_Sundials);CHKERRQ(ierr); 911 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr); 912 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr); 913 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",TSSundialsGetPC_Sundials);CHKERRQ(ierr); 914 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",TSSundialsGetIterations_Sundials);CHKERRQ(ierr); 915 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr); 916 PetscFunctionReturn(0); 917 } 918