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