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