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 order of BDF / ADAMS method */ 345 if (cvode->maxord != PETSC_DEFAULT) { 346 flag = CVodeSetMaxOrd(mem,cvode->maxord); 347 if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetMaxOrd() fails, flag %d",flag); 348 } 349 350 /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */ 351 flag = CVodeSetMaxNumSteps(mem,ts->max_steps); 352 if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetMaxNumSteps() fails, flag %d",flag); 353 354 /* call CVSpgmr to use GMRES as the linear solver. */ 355 /* setup the ode integrator with the given preconditioner */ 356 ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 357 ierr = PCGetType(pc,&pctype);CHKERRQ(ierr); 358 ierr = PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr); 359 if (pcnone) { 360 flag = CVSpgmr(mem,PREC_NONE,0); 361 if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 362 } else { 363 flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl); 364 if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 365 366 /* Set preconditioner and solve routines Precond and PSolve, 367 and the pointer to the user-defined block data */ 368 flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials); 369 if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag); 370 } 371 372 flag = CVSpilsSetGSType(mem, MODIFIED_GS); 373 if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag); 374 PetscFunctionReturn(0); 375 } 376 377 /* type of CVODE linear multistep method */ 378 const char *const TSSundialsLmmTypes[] = {"","ADAMS","BDF","TSSundialsLmmType","SUNDIALS_",0}; 379 /* type of G-S orthogonalization used by CVODE linear solver */ 380 const char *const TSSundialsGramSchmidtTypes[] = {"","MODIFIED","CLASSICAL","TSSundialsGramSchmidtType","SUNDIALS_",0}; 381 382 PetscErrorCode TSSetFromOptions_Sundials(PetscOptionItems *PetscOptionsObject,TS ts) 383 { 384 TS_Sundials *cvode = (TS_Sundials*)ts->data; 385 PetscErrorCode ierr; 386 int indx; 387 PetscBool flag; 388 PC pc; 389 390 PetscFunctionBegin; 391 ierr = PetscOptionsHead(PetscOptionsObject,"SUNDIALS ODE solver options");CHKERRQ(ierr); 392 ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr); 393 if (flag) { 394 ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr); 395 } 396 ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr); 397 if (flag) { 398 ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr); 399 } 400 ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,NULL);CHKERRQ(ierr); 401 ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,NULL);CHKERRQ(ierr); 402 ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,NULL);CHKERRQ(ierr); 403 ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,NULL);CHKERRQ(ierr); 404 ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,NULL);CHKERRQ(ierr); 405 ierr = PetscOptionsInt("-ts_sundials_maxord","Max Order for BDF/Adams method","TSSundialsSetMaxOrd",cvode->maxord,&cvode->maxord,NULL);CHKERRQ(ierr); 406 ierr = PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,NULL);CHKERRQ(ierr); 407 ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internal steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,NULL);CHKERRQ(ierr); 408 ierr = PetscOptionsTail();CHKERRQ(ierr); 409 ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 410 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 411 PetscFunctionReturn(0); 412 } 413 414 PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer) 415 { 416 TS_Sundials *cvode = (TS_Sundials*)ts->data; 417 PetscErrorCode ierr; 418 char *type; 419 char atype[] = "Adams"; 420 char btype[] = "BDF: backward differentiation formula"; 421 PetscBool iascii,isstring; 422 long int nsteps,its,nfevals,nlinsetups,nfails,itmp; 423 PetscInt qlast,qcur; 424 PetscReal hinused,hlast,hcur,tcur,tolsfac; 425 PC pc; 426 427 PetscFunctionBegin; 428 if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype; 429 else type = btype; 430 431 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 432 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 433 if (iascii) { 434 ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr); 435 ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr); 436 ierr = PetscViewerASCIIPrintf(viewer,"Sundials maxord %D\n",cvode->maxord);CHKERRQ(ierr); 437 ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",(double)cvode->abstol,(double)cvode->reltol);CHKERRQ(ierr); 438 ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",(double)cvode->linear_tol);CHKERRQ(ierr); 439 ierr = PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);CHKERRQ(ierr); 440 if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 441 ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 442 } else { 443 ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 444 } 445 if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",(double)cvode->mindt);CHKERRQ(ierr);} 446 if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",(double)cvode->maxdt);CHKERRQ(ierr);} 447 448 /* Outputs from CVODE, CVSPILS */ 449 ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr); 450 ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr); 451 ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals, 452 &nlinsetups,&nfails,&qlast,&qcur, 453 &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr); 454 ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr); 455 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr); 456 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr); 457 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr); 458 459 ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr); 460 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr); 461 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr); 462 463 ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */ 464 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr); 465 ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr); 466 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr); 467 468 ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 469 ierr = PCView(pc,viewer);CHKERRQ(ierr); 470 ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr); 471 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr); 472 ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr); 473 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr); 474 475 ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr); 476 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr); 477 ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr); 478 ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr); 479 } else if (isstring) { 480 ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr); 481 } 482 PetscFunctionReturn(0); 483 } 484 485 486 /* --------------------------------------------------------------------------*/ 487 PetscErrorCode TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type) 488 { 489 TS_Sundials *cvode = (TS_Sundials*)ts->data; 490 491 PetscFunctionBegin; 492 cvode->cvode_type = type; 493 PetscFunctionReturn(0); 494 } 495 496 PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl) 497 { 498 TS_Sundials *cvode = (TS_Sundials*)ts->data; 499 500 PetscFunctionBegin; 501 cvode->maxl = maxl; 502 PetscFunctionReturn(0); 503 } 504 505 PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts,double tol) 506 { 507 TS_Sundials *cvode = (TS_Sundials*)ts->data; 508 509 PetscFunctionBegin; 510 cvode->linear_tol = tol; 511 PetscFunctionReturn(0); 512 } 513 514 PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type) 515 { 516 TS_Sundials *cvode = (TS_Sundials*)ts->data; 517 518 PetscFunctionBegin; 519 cvode->gtype = type; 520 PetscFunctionReturn(0); 521 } 522 523 PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel) 524 { 525 TS_Sundials *cvode = (TS_Sundials*)ts->data; 526 527 PetscFunctionBegin; 528 if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 529 if (rel != PETSC_DECIDE) cvode->reltol = rel; 530 PetscFunctionReturn(0); 531 } 532 533 PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt) 534 { 535 TS_Sundials *cvode = (TS_Sundials*)ts->data; 536 537 PetscFunctionBegin; 538 cvode->mindt = mindt; 539 PetscFunctionReturn(0); 540 } 541 542 PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt) 543 { 544 TS_Sundials *cvode = (TS_Sundials*)ts->data; 545 546 PetscFunctionBegin; 547 cvode->maxdt = maxdt; 548 PetscFunctionReturn(0); 549 } 550 PetscErrorCode TSSundialsGetPC_Sundials(TS ts,PC *pc) 551 { 552 SNES snes; 553 KSP ksp; 554 PetscErrorCode ierr; 555 556 PetscFunctionBegin; 557 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 558 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 559 ierr = KSPGetPC(ksp,pc);CHKERRQ(ierr); 560 PetscFunctionReturn(0); 561 } 562 563 PetscErrorCode TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin) 564 { 565 PetscFunctionBegin; 566 if (nonlin) *nonlin = ts->snes_its; 567 if (lin) *lin = ts->ksp_its; 568 PetscFunctionReturn(0); 569 } 570 571 PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s) 572 { 573 TS_Sundials *cvode = (TS_Sundials*)ts->data; 574 575 PetscFunctionBegin; 576 cvode->monitorstep = s; 577 PetscFunctionReturn(0); 578 } 579 /* -------------------------------------------------------------------------------------------*/ 580 581 /*@C 582 TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 583 584 Not Collective 585 586 Input parameters: 587 . ts - the time-step context 588 589 Output Parameters: 590 + nonlin - number of nonlinear iterations 591 - lin - number of linear iterations 592 593 Level: advanced 594 595 Notes: 596 These return the number since the creation of the TS object 597 598 .seealso: TSSundialsSetType(), TSSundialsSetMaxl(), 599 TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 600 TSSundialsGetIterations(), TSSundialsSetType(), 601 TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime() 602 603 @*/ 604 PetscErrorCode TSSundialsGetIterations(TS ts,int *nonlin,int *lin) 605 { 606 PetscErrorCode ierr; 607 608 PetscFunctionBegin; 609 ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr); 610 PetscFunctionReturn(0); 611 } 612 613 /*@ 614 TSSundialsSetType - Sets the method that Sundials will use for integration. 615 616 Logically Collective on TS 617 618 Input parameters: 619 + ts - the time-step context 620 - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 621 622 Level: intermediate 623 624 .seealso: TSSundialsGetIterations(), TSSundialsSetMaxl(), 625 TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 626 TSSundialsGetIterations(), TSSundialsSetType(), 627 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 628 TSSetExactFinalTime() 629 @*/ 630 PetscErrorCode TSSundialsSetType(TS ts,TSSundialsLmmType type) 631 { 632 PetscErrorCode ierr; 633 634 PetscFunctionBegin; 635 ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr); 636 PetscFunctionReturn(0); 637 } 638 639 /*@ 640 TSSundialsSetMaxord - Sets the maximum order for BDF/Adams method used by SUNDIALS. 641 642 Logically Collective on TS 643 644 Input parameters: 645 + ts - the time-step context 646 - maxord - maximum order of BDF / Adams method 647 648 Level: advanced 649 650 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 651 TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 652 TSSundialsGetIterations(), TSSundialsSetType(), 653 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 654 TSSetExactFinalTime() 655 656 @*/ 657 PetscErrorCode TSSundialsSetMaxord(TS ts,PetscInt maxord) 658 { 659 PetscErrorCode ierr; 660 661 PetscFunctionBegin; 662 PetscValidLogicalCollectiveInt(ts,maxord,2); 663 ierr = PetscTryMethod(ts,"TSSundialsSetMaxOrd_C",(TS,PetscInt),(ts,maxord));CHKERRQ(ierr); 664 PetscFunctionReturn(0); 665 } 666 667 /*@ 668 TSSundialsSetMaxl - Sets the dimension of the Krylov space used by 669 GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 670 this is the maximum number of GMRES steps that will be used. 671 672 Logically Collective on TS 673 674 Input parameters: 675 + ts - the time-step context 676 - maxl - number of direction vectors (the dimension of Krylov subspace). 677 678 Level: advanced 679 680 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 681 TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 682 TSSundialsGetIterations(), TSSundialsSetType(), 683 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 684 TSSetExactFinalTime() 685 686 @*/ 687 PetscErrorCode TSSundialsSetMaxl(TS ts,PetscInt maxl) 688 { 689 PetscErrorCode ierr; 690 691 PetscFunctionBegin; 692 PetscValidLogicalCollectiveInt(ts,maxl,2); 693 ierr = PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));CHKERRQ(ierr); 694 PetscFunctionReturn(0); 695 } 696 697 /*@ 698 TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 699 system by SUNDIALS. 700 701 Logically Collective on TS 702 703 Input parameters: 704 + ts - the time-step context 705 - tol - the factor by which the tolerance on the nonlinear solver is 706 multiplied to get the tolerance on the linear solver, .05 by default. 707 708 Level: advanced 709 710 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 711 TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 712 TSSundialsGetIterations(), TSSundialsSetType(), 713 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 714 TSSetExactFinalTime() 715 716 @*/ 717 PetscErrorCode TSSundialsSetLinearTolerance(TS ts,double tol) 718 { 719 PetscErrorCode ierr; 720 721 PetscFunctionBegin; 722 PetscValidLogicalCollectiveReal(ts,tol,2); 723 ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr); 724 PetscFunctionReturn(0); 725 } 726 727 /*@ 728 TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 729 in GMRES method by SUNDIALS linear solver. 730 731 Logically Collective on TS 732 733 Input parameters: 734 + ts - the time-step context 735 - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 736 737 Level: advanced 738 739 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 740 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), 741 TSSundialsGetIterations(), TSSundialsSetType(), 742 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 743 TSSetExactFinalTime() 744 745 @*/ 746 PetscErrorCode TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type) 747 { 748 PetscErrorCode ierr; 749 750 PetscFunctionBegin; 751 ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr); 752 PetscFunctionReturn(0); 753 } 754 755 /*@ 756 TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 757 Sundials for error control. 758 759 Logically Collective on TS 760 761 Input parameters: 762 + ts - the time-step context 763 . aabs - the absolute tolerance 764 - rel - the relative tolerance 765 766 See the Cvode/Sundials users manual for exact details on these parameters. Essentially 767 these regulate the size of the error for a SINGLE timestep. 768 769 Level: intermediate 770 771 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(), 772 TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 773 TSSundialsGetIterations(), TSSundialsSetType(), 774 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 775 TSSetExactFinalTime() 776 777 @*/ 778 PetscErrorCode TSSundialsSetTolerance(TS ts,double aabs,double rel) 779 { 780 PetscErrorCode ierr; 781 782 PetscFunctionBegin; 783 ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr); 784 PetscFunctionReturn(0); 785 } 786 787 /*@ 788 TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 789 790 Input Parameter: 791 . ts - the time-step context 792 793 Output Parameter: 794 . pc - the preconditioner context 795 796 Level: advanced 797 798 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 799 TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 800 TSSundialsGetIterations(), TSSundialsSetType(), 801 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance() 802 @*/ 803 PetscErrorCode TSSundialsGetPC(TS ts,PC *pc) 804 { 805 PetscErrorCode ierr; 806 807 PetscFunctionBegin; 808 ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC*),(ts,pc));CHKERRQ(ierr); 809 PetscFunctionReturn(0); 810 } 811 812 /*@ 813 TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 814 815 Input Parameter: 816 + ts - the time-step context 817 - mindt - lowest time step if positive, negative to deactivate 818 819 Note: 820 Sundials will error if it is not possible to keep the estimated truncation error below 821 the tolerance set with TSSundialsSetTolerance() without going below this step size. 822 823 Level: beginner 824 825 .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 826 @*/ 827 PetscErrorCode TSSundialsSetMinTimeStep(TS ts,PetscReal mindt) 828 { 829 PetscErrorCode ierr; 830 831 PetscFunctionBegin; 832 ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr); 833 PetscFunctionReturn(0); 834 } 835 836 /*@ 837 TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 838 839 Input Parameter: 840 + ts - the time-step context 841 - maxdt - lowest time step if positive, negative to deactivate 842 843 Level: beginner 844 845 .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 846 @*/ 847 PetscErrorCode TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt) 848 { 849 PetscErrorCode ierr; 850 851 PetscFunctionBegin; 852 ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 853 PetscFunctionReturn(0); 854 } 855 856 /*@ 857 TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 858 859 Input Parameter: 860 + ts - the time-step context 861 - ft - PETSC_TRUE if monitor, else PETSC_FALSE 862 863 Level: beginner 864 865 .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 866 TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 867 TSSundialsGetIterations(), TSSundialsSetType(), 868 TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 869 @*/ 870 PetscErrorCode TSSundialsMonitorInternalSteps(TS ts,PetscBool ft) 871 { 872 PetscErrorCode ierr; 873 874 PetscFunctionBegin; 875 ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr); 876 PetscFunctionReturn(0); 877 } 878 /* -------------------------------------------------------------------------------------------*/ 879 /*MC 880 TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 881 882 Options Database: 883 + -ts_sundials_type <bdf,adams> - 884 . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 885 . -ts_sundials_atol <tol> - Absolute tolerance for convergence 886 . -ts_sundials_rtol <tol> - Relative tolerance for convergence 887 . -ts_sundials_linear_tolerance <tol> - 888 . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace 889 - -ts_sundials_monitor_steps - Monitor SUNDIALS internal steps 890 891 892 Notes: 893 This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply, 894 only PETSc PC options. 895 896 Level: beginner 897 898 .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(), 899 TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime() 900 901 M*/ 902 PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts) 903 { 904 TS_Sundials *cvode; 905 PetscErrorCode ierr; 906 PC pc; 907 908 PetscFunctionBegin; 909 ts->ops->reset = TSReset_Sundials; 910 ts->ops->destroy = TSDestroy_Sundials; 911 ts->ops->view = TSView_Sundials; 912 ts->ops->setup = TSSetUp_Sundials; 913 ts->ops->step = TSStep_Sundials; 914 ts->ops->interpolate = TSInterpolate_Sundials; 915 ts->ops->setfromoptions = TSSetFromOptions_Sundials; 916 ts->default_adapt_type = TSADAPTNONE; 917 918 ierr = PetscNewLog(ts,&cvode);CHKERRQ(ierr); 919 920 ts->usessnes = PETSC_TRUE; 921 922 ts->data = (void*)cvode; 923 cvode->cvode_type = SUNDIALS_BDF; 924 cvode->gtype = SUNDIALS_CLASSICAL_GS; 925 cvode->maxl = 5; 926 cvode->maxord = PETSC_DEFAULT; 927 cvode->linear_tol = .05; 928 cvode->monitorstep = PETSC_TRUE; 929 930 ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)ts),&(cvode->comm_sundials));CHKERRQ(ierr); 931 932 cvode->mindt = -1.; 933 cvode->maxdt = -1.; 934 935 /* set tolerance for Sundials */ 936 cvode->reltol = 1e-6; 937 cvode->abstol = 1e-6; 938 939 /* set PCNONE as default pctype */ 940 ierr = TSSundialsGetPC_Sundials(ts,&pc);CHKERRQ(ierr); 941 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 942 943 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",TSSundialsSetType_Sundials);CHKERRQ(ierr); 944 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",TSSundialsSetMaxl_Sundials);CHKERRQ(ierr); 945 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr); 946 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr); 947 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",TSSundialsSetTolerance_Sundials);CHKERRQ(ierr); 948 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr); 949 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr); 950 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",TSSundialsGetPC_Sundials);CHKERRQ(ierr); 951 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",TSSundialsGetIterations_Sundials);CHKERRQ(ierr); 952 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr); 953 PetscFunctionReturn(0); 954 } 955