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