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