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