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