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