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