17cb94ee6SHong Zhang /* 2db268bfaSHong Zhang Provides a PETSc interface to SUNDIALS/CVODE solver. 37cb94ee6SHong Zhang The interface to PVODE (old version of CVODE) was originally contributed 47cb94ee6SHong Zhang by Liyang Xu. It has been redone by Hong Zhang and Dinesh Kaushik. 528adb3f7SHong Zhang 628adb3f7SHong Zhang Reference: sundials-2.4.0/examples/cvode/parallel/cvDiurnal_kry_p.c 77cb94ee6SHong Zhang */ 856a740aaSHong Zhang #include "sundials.h" /*I "petscts.h" I*/ 97cb94ee6SHong Zhang 107cb94ee6SHong Zhang /* 117cb94ee6SHong Zhang TSPrecond_Sundials - function that we provide to SUNDIALS to 127cb94ee6SHong Zhang evaluate the preconditioner. 137cb94ee6SHong Zhang */ 147cb94ee6SHong Zhang #undef __FUNCT__ 157cb94ee6SHong Zhang #define __FUNCT__ "TSPrecond_Sundials" 167cb94ee6SHong Zhang PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy, 177cb94ee6SHong Zhang booleantype jok,booleantype *jcurPtr, 187cb94ee6SHong Zhang realtype _gamma,void *P_data, 197cb94ee6SHong Zhang N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3) 207cb94ee6SHong Zhang { 217cb94ee6SHong Zhang TS ts = (TS) P_data; 227cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 23dcbc6d53SSean Farley PC pc; 247cb94ee6SHong Zhang PetscErrorCode ierr; 250679a0aeSJed Brown Mat J,P; 260679a0aeSJed Brown Vec yy = cvode->w1,yydot = cvode->ydot; 270679a0aeSJed Brown PetscReal gm = (PetscReal)_gamma; 287cb94ee6SHong Zhang MatStructure str = DIFFERENT_NONZERO_PATTERN; 297cb94ee6SHong Zhang PetscScalar *y_data; 307cb94ee6SHong Zhang 317cb94ee6SHong Zhang PetscFunctionBegin; 320679a0aeSJed Brown ierr = TSGetIJacobian(ts,&J,&P,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 337cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(y); 347cb94ee6SHong Zhang ierr = VecPlaceArray(yy,y_data); CHKERRQ(ierr); 350679a0aeSJed Brown ierr = VecZeroEntries(yydot);CHKERRQ(ierr); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */ 360679a0aeSJed Brown /* compute the shifted Jacobian (1/gm)*I + Jrest */ 370679a0aeSJed Brown ierr = TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,&J,&P,&str,PETSC_FALSE);CHKERRQ(ierr); 387cb94ee6SHong Zhang ierr = VecResetArray(yy); CHKERRQ(ierr); 390679a0aeSJed Brown ierr = MatScale(P,gm);CHKERRQ(ierr); /* turn into I-gm*Jrest, J is not used by Sundials */ 407cb94ee6SHong Zhang *jcurPtr = TRUE; 41dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr); 420679a0aeSJed Brown ierr = PCSetOperators(pc,J,P,str);CHKERRQ(ierr); 437cb94ee6SHong Zhang PetscFunctionReturn(0); 447cb94ee6SHong Zhang } 457cb94ee6SHong Zhang 467cb94ee6SHong Zhang /* 477cb94ee6SHong Zhang TSPSolve_Sundials - routine that we provide to Sundials that applies the preconditioner. 487cb94ee6SHong Zhang */ 497cb94ee6SHong Zhang #undef __FUNCT__ 507cb94ee6SHong Zhang #define __FUNCT__ "TSPSolve_Sundials" 514ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z, 524ac7836bSHong Zhang realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp) 537cb94ee6SHong Zhang { 547cb94ee6SHong Zhang TS ts = (TS) P_data; 557cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 56dcbc6d53SSean Farley PC pc; 577cb94ee6SHong Zhang Vec rr = cvode->w1,zz = cvode->w2; 587cb94ee6SHong Zhang PetscErrorCode ierr; 597cb94ee6SHong Zhang PetscScalar *r_data,*z_data; 607cb94ee6SHong Zhang 617cb94ee6SHong Zhang PetscFunctionBegin; 627cb94ee6SHong Zhang /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/ 637cb94ee6SHong Zhang r_data = (PetscScalar *) N_VGetArrayPointer(r); 647cb94ee6SHong Zhang z_data = (PetscScalar *) N_VGetArrayPointer(z); 657cb94ee6SHong Zhang ierr = VecPlaceArray(rr,r_data); CHKERRQ(ierr); 667cb94ee6SHong Zhang ierr = VecPlaceArray(zz,z_data); CHKERRQ(ierr); 674ac7836bSHong Zhang 687cb94ee6SHong Zhang /* Solve the Px=r and put the result in zz */ 69dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr); 707cb94ee6SHong Zhang ierr = PCApply(pc,rr,zz); CHKERRQ(ierr); 717cb94ee6SHong Zhang ierr = VecResetArray(rr); CHKERRQ(ierr); 727cb94ee6SHong Zhang ierr = VecResetArray(zz); CHKERRQ(ierr); 737cb94ee6SHong Zhang PetscFunctionReturn(0); 747cb94ee6SHong Zhang } 757cb94ee6SHong Zhang 767cb94ee6SHong Zhang /* 777cb94ee6SHong Zhang TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side. 787cb94ee6SHong Zhang */ 797cb94ee6SHong Zhang #undef __FUNCT__ 807cb94ee6SHong Zhang #define __FUNCT__ "TSFunction_Sundials" 814ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx) 827cb94ee6SHong Zhang { 837cb94ee6SHong Zhang TS ts = (TS) ctx; 847adad957SLisandro Dalcin MPI_Comm comm = ((PetscObject)ts)->comm; 857cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 860679a0aeSJed Brown Vec yy = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot; 877cb94ee6SHong Zhang PetscScalar *y_data,*ydot_data; 887cb94ee6SHong Zhang PetscErrorCode ierr; 897cb94ee6SHong Zhang 907cb94ee6SHong Zhang PetscFunctionBegin; 915ff03cf7SDinesh Kaushik /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/ 927cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(y); 937cb94ee6SHong Zhang ydot_data = (PetscScalar *) N_VGetArrayPointer(ydot); 944ac7836bSHong Zhang ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr); 954ac7836bSHong Zhang ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr); 964ac7836bSHong Zhang 977cb94ee6SHong Zhang /* now compute the right hand side function */ 98bc0cc02bSJed Brown if (!ts->userops->ifunction) { 99bc0cc02bSJed Brown ierr = TSComputeRHSFunction(ts,t,yy,yyd);CHKERRQ(ierr); 100bc0cc02bSJed Brown } else { /* If rhsfunction is also set, this computes both parts and shifts them to the right */ 101bc0cc02bSJed Brown ierr = VecZeroEntries(yydot);CHKERRQ(ierr); 1020679a0aeSJed Brown ierr = TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE); CHKERRABORT(comm,ierr); 1030679a0aeSJed Brown ierr = VecScale(yyd,-1.);CHKERRQ(ierr); 104bc0cc02bSJed Brown } 1057cb94ee6SHong Zhang ierr = VecResetArray(yy); CHKERRABORT(comm,ierr); 1067cb94ee6SHong Zhang ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr); 1074ac7836bSHong Zhang PetscFunctionReturn(0); 1087cb94ee6SHong Zhang } 1097cb94ee6SHong Zhang 1107cb94ee6SHong Zhang /* 111b4eba00bSSean Farley TSStep_Sundials - Calls Sundials to integrate the ODE. 1127cb94ee6SHong Zhang */ 1137cb94ee6SHong Zhang #undef __FUNCT__ 114b4eba00bSSean Farley #define __FUNCT__ "TSStep_Sundials" 115b4eba00bSSean Farley PetscErrorCode TSStep_Sundials(TS ts) 1167cb94ee6SHong Zhang { 1177cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 1187cb94ee6SHong Zhang PetscErrorCode ierr; 119b4eba00bSSean Farley PetscInt flag; 1209f94935aSHong Zhang long int its,nsteps; 1217cb94ee6SHong Zhang realtype t,tout; 1227cb94ee6SHong Zhang PetscScalar *y_data; 1237cb94ee6SHong Zhang void *mem; 1247cb94ee6SHong Zhang 1257cb94ee6SHong Zhang PetscFunctionBegin; 12616016685SHong Zhang mem = cvode->mem; 1277cb94ee6SHong Zhang tout = ts->max_time; 1287cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr); 1297cb94ee6SHong Zhang N_VSetArrayPointer((realtype *)y_data,cvode->y); 1307cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 131186e87acSLisandro Dalcin 1323f2090d5SJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 133186e87acSLisandro Dalcin 1342bfc04deSHong Zhang if (cvode->monitorstep) { 135*b8123daeSJed Brown /* We would like to call TSPreStep() when starting each step (including rejections) and TSPreStage() before each 136*b8123daeSJed Brown * stage solve, but CVode does not appear to support this. */ 13728adb3f7SHong Zhang flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP); 1382bfc04deSHong Zhang } else { 1392bfc04deSHong Zhang flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL); 1402bfc04deSHong Zhang } 1419f94935aSHong Zhang 1429f94935aSHong Zhang if (flag){ /* display error message */ 1439f94935aSHong Zhang switch (flag){ 1449f94935aSHong Zhang case CV_ILL_INPUT: 1459f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT"); 1469f94935aSHong Zhang break; 1479f94935aSHong Zhang case CV_TOO_CLOSE: 1489f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE"); 1499f94935aSHong Zhang break; 1503c7fefeeSJed Brown case CV_TOO_MUCH_WORK: { 1519f94935aSHong Zhang PetscReal tcur; 1529f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 1539f94935aSHong Zhang ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr); 1549f94935aSHong Zhang 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); 1553c7fefeeSJed Brown } break; 1569f94935aSHong Zhang case CV_TOO_MUCH_ACC: 1579f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC"); 1589f94935aSHong Zhang break; 1599f94935aSHong Zhang case CV_ERR_FAILURE: 1609f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE"); 1619f94935aSHong Zhang break; 1629f94935aSHong Zhang case CV_CONV_FAILURE: 1639f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE"); 1649f94935aSHong Zhang break; 1659f94935aSHong Zhang case CV_LINIT_FAIL: 1669f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL"); 1679f94935aSHong Zhang break; 1689f94935aSHong Zhang case CV_LSETUP_FAIL: 1699f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL"); 1709f94935aSHong Zhang break; 1719f94935aSHong Zhang case CV_LSOLVE_FAIL: 1729f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL"); 1739f94935aSHong Zhang break; 1749f94935aSHong Zhang case CV_RHSFUNC_FAIL: 1759f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL"); 1769f94935aSHong Zhang break; 1779f94935aSHong Zhang case CV_FIRST_RHSFUNC_ERR: 1789f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR"); 1799f94935aSHong Zhang break; 1809f94935aSHong Zhang case CV_REPTD_RHSFUNC_ERR: 1819f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR"); 1829f94935aSHong Zhang break; 1839f94935aSHong Zhang case CV_UNREC_RHSFUNC_ERR: 1849f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR"); 1859f94935aSHong Zhang break; 1869f94935aSHong Zhang case CV_RTFUNC_FAIL: 1879f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL"); 1889f94935aSHong Zhang break; 1899f94935aSHong Zhang default: 1909f94935aSHong Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag); 1919f94935aSHong Zhang } 1929f94935aSHong Zhang } 1939f94935aSHong Zhang 1947cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */ 1957cb94ee6SHong Zhang ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr); 1967cb94ee6SHong Zhang ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr); 1977cb94ee6SHong Zhang ierr = VecResetArray(cvode->w1); CHKERRQ(ierr); 198bc0cc02bSJed Brown ierr = VecCopy(cvode->update,ts->vec_sol);CHKERRQ(ierr); 1997cb94ee6SHong Zhang ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr); 2004ac7836bSHong Zhang ierr = CVSpilsGetNumLinIters(mem, &its); 2015ef26d82SJed Brown ts->snes_its = its; ts->ksp_its = its; 202186e87acSLisandro Dalcin 203186e87acSLisandro Dalcin ts->time_step = t - ts->ptime; 204186e87acSLisandro Dalcin ts->ptime = t; 2057cb94ee6SHong Zhang ts->steps++; 206186e87acSLisandro Dalcin 2079f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 208b4eba00bSSean Farley if (!cvode->monitorstep) ts->steps = nsteps; 209b4eba00bSSean Farley PetscFunctionReturn(0); 210b4eba00bSSean Farley } 211b4eba00bSSean Farley 212b4eba00bSSean Farley #undef __FUNCT__ 213b4eba00bSSean Farley #define __FUNCT__ "TSInterpolate_Sundials" 214b4eba00bSSean Farley static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X) 215b4eba00bSSean Farley { 216b4eba00bSSean Farley TS_Sundials *cvode = (TS_Sundials*)ts->data; 217b4eba00bSSean Farley N_Vector y; 218b4eba00bSSean Farley PetscErrorCode ierr; 219b4eba00bSSean Farley PetscScalar *x_data; 220b4eba00bSSean Farley PetscInt glosize,locsize; 221b4eba00bSSean Farley 222b4eba00bSSean Farley PetscFunctionBegin; 223b4eba00bSSean Farley 224b4eba00bSSean Farley /* get the vector size */ 225b4eba00bSSean Farley ierr = VecGetSize(X,&glosize);CHKERRQ(ierr); 226b4eba00bSSean Farley ierr = VecGetLocalSize(X,&locsize);CHKERRQ(ierr); 227b4eba00bSSean Farley 228b4eba00bSSean Farley /* allocate the memory for N_Vec y */ 229b4eba00bSSean Farley y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 230b4eba00bSSean Farley if (!y) SETERRQ(PETSC_COMM_SELF,1,"Interpolated y is not allocated"); 231b4eba00bSSean Farley 232b4eba00bSSean Farley ierr = VecGetArray(X,&x_data);CHKERRQ(ierr); 233b4eba00bSSean Farley N_VSetArrayPointer((realtype *)x_data,y); 234b4eba00bSSean Farley ierr = CVodeGetDky(cvode->mem,t,0,y);CHKERRQ(ierr); 235b4eba00bSSean Farley ierr = VecRestoreArray(X,&x_data);CHKERRQ(ierr); 236b4eba00bSSean Farley 2377cb94ee6SHong Zhang PetscFunctionReturn(0); 2387cb94ee6SHong Zhang } 2397cb94ee6SHong Zhang 2407cb94ee6SHong Zhang #undef __FUNCT__ 241277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Sundials" 242277b19d0SLisandro Dalcin PetscErrorCode TSReset_Sundials(TS ts) 243277b19d0SLisandro Dalcin { 244277b19d0SLisandro Dalcin TS_Sundials *cvode = (TS_Sundials*)ts->data; 245277b19d0SLisandro Dalcin PetscErrorCode ierr; 246277b19d0SLisandro Dalcin 247277b19d0SLisandro Dalcin PetscFunctionBegin; 2485c6b4a3dSJed Brown ierr = VecDestroy(&cvode->update);CHKERRQ(ierr); 2490679a0aeSJed Brown ierr = VecDestroy(&cvode->ydot);CHKERRQ(ierr); 2505c6b4a3dSJed Brown ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr); 2515c6b4a3dSJed Brown ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr); 252277b19d0SLisandro Dalcin if (cvode->mem) {CVodeFree(&cvode->mem);} 253277b19d0SLisandro Dalcin PetscFunctionReturn(0); 254277b19d0SLisandro Dalcin } 255277b19d0SLisandro Dalcin 256277b19d0SLisandro Dalcin #undef __FUNCT__ 2577cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials" 2587cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts) 2597cb94ee6SHong Zhang { 2607cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2617cb94ee6SHong Zhang PetscErrorCode ierr; 2627cb94ee6SHong Zhang 2637cb94ee6SHong Zhang PetscFunctionBegin; 264277b19d0SLisandro Dalcin ierr = TSReset_Sundials(ts);CHKERRQ(ierr); 265a07356b8SSatish Balay ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr); 266277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 267335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","",PETSC_NULL);CHKERRQ(ierr); 268f61b2b6aSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxl_C","",PETSC_NULL);CHKERRQ(ierr); 269335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C","",PETSC_NULL);CHKERRQ(ierr); 270335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C","",PETSC_NULL);CHKERRQ(ierr); 271335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C","",PETSC_NULL);CHKERRQ(ierr); 272335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C","",PETSC_NULL);CHKERRQ(ierr); 273335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C","",PETSC_NULL);CHKERRQ(ierr); 274335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C","",PETSC_NULL);CHKERRQ(ierr); 275335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C","",PETSC_NULL);CHKERRQ(ierr); 276335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C","",PETSC_NULL);CHKERRQ(ierr); 2777cb94ee6SHong Zhang PetscFunctionReturn(0); 2787cb94ee6SHong Zhang } 2797cb94ee6SHong Zhang 2807cb94ee6SHong Zhang #undef __FUNCT__ 281214bc6a2SJed Brown #define __FUNCT__ "TSSetUp_Sundials" 282214bc6a2SJed Brown PetscErrorCode TSSetUp_Sundials(TS ts) 2837cb94ee6SHong Zhang { 2847cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2857cb94ee6SHong Zhang PetscErrorCode ierr; 28616016685SHong Zhang PetscInt glosize,locsize,i,flag; 2877cb94ee6SHong Zhang PetscScalar *y_data,*parray; 28816016685SHong Zhang void *mem; 289dcbc6d53SSean Farley PC pc; 29016016685SHong Zhang const PCType pctype; 291ace3abfcSBarry Smith PetscBool pcnone; 2927cb94ee6SHong Zhang 2937cb94ee6SHong Zhang PetscFunctionBegin; 2947cb94ee6SHong Zhang /* get the vector size */ 2957cb94ee6SHong Zhang ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr); 2967cb94ee6SHong Zhang ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr); 2977cb94ee6SHong Zhang 2987cb94ee6SHong Zhang /* allocate the memory for N_Vec y */ 299a07356b8SSatish Balay cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 300e32f2f54SBarry Smith if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated"); 3017cb94ee6SHong Zhang 30228adb3f7SHong Zhang /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 3037cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr); 3047cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y); 3057cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 3067cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 30742528757SHong Zhang 3087cb94ee6SHong Zhang ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr); 3090679a0aeSJed Brown ierr = VecDuplicate(ts->vec_sol,&cvode->ydot);CHKERRQ(ierr); 3107cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr); 3110679a0aeSJed Brown ierr = PetscLogObjectParent(ts,cvode->ydot);CHKERRQ(ierr); 3127cb94ee6SHong Zhang 3137cb94ee6SHong Zhang /* 3147cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are 3157cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided 3167cb94ee6SHong Zhang by Sundials and set using VecPlaceArray(). 3177cb94ee6SHong Zhang */ 318778a2246SBarry Smith ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,1,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr); 319778a2246SBarry Smith ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,1,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr); 3207cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr); 3217cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr); 32216016685SHong Zhang 32316016685SHong Zhang /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 32416016685SHong Zhang mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 325e32f2f54SBarry Smith if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails"); 32616016685SHong Zhang cvode->mem = mem; 32716016685SHong Zhang 32816016685SHong Zhang /* Set the pointer to user-defined data */ 32916016685SHong Zhang flag = CVodeSetUserData(mem, ts); 330e32f2f54SBarry Smith if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails"); 33116016685SHong Zhang 332fc6b9e64SJed Brown /* Sundials may choose to use a smaller initial step, but will never use a larger step. */ 333cdbf8f93SLisandro Dalcin flag = CVodeSetInitStep(mem,(realtype)ts->time_step); 334fc6b9e64SJed Brown if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed"); 335f1cd61daSJed Brown if (cvode->mindt > 0) { 336f1cd61daSJed Brown flag = CVodeSetMinStep(mem,(realtype)cvode->mindt); 3379f94935aSHong Zhang if (flag){ 3389f94935aSHong Zhang if (flag == CV_MEM_NULL){ 3399f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL"); 3409f94935aSHong Zhang } else if (flag == CV_ILL_INPUT){ 3419f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size"); 3429f94935aSHong Zhang } else { 3439f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed"); 3449f94935aSHong Zhang } 3459f94935aSHong Zhang } 346f1cd61daSJed Brown } 347f1cd61daSJed Brown if (cvode->maxdt > 0) { 348f1cd61daSJed Brown flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt); 349f1cd61daSJed Brown if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed"); 350f1cd61daSJed Brown } 351f1cd61daSJed Brown 35216016685SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the 35316016685SHong Zhang * user's right hand side function in u'=f(t,u), the inital time T0, and 35416016685SHong Zhang * the initial dependent variable vector cvode->y */ 35516016685SHong Zhang flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y); 35616016685SHong Zhang if (flag){ 357e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag); 35816016685SHong Zhang } 35916016685SHong Zhang 3609f94935aSHong Zhang /* specifies scalar relative and absolute tolerances */ 36116016685SHong Zhang flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol); 36216016685SHong Zhang if (flag){ 363e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag); 36416016685SHong Zhang } 36516016685SHong Zhang 3669f94935aSHong Zhang /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */ 3679f94935aSHong Zhang flag = CVodeSetMaxNumSteps(mem,ts->max_steps); 36816016685SHong Zhang 36916016685SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 37016016685SHong Zhang /* setup the ode integrator with the given preconditioner */ 371dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr); 372dcbc6d53SSean Farley ierr = PCGetType(pc,&pctype);CHKERRQ(ierr); 373251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr); 37416016685SHong Zhang if (pcnone){ 37516016685SHong Zhang flag = CVSpgmr(mem,PREC_NONE,0); 376e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 37716016685SHong Zhang } else { 378f61b2b6aSHong Zhang flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl); 379e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 38016016685SHong Zhang 38116016685SHong Zhang /* Set preconditioner and solve routines Precond and PSolve, 38216016685SHong Zhang and the pointer to the user-defined block data */ 38316016685SHong Zhang flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials); 384e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag); 38516016685SHong Zhang } 38616016685SHong Zhang 38716016685SHong Zhang flag = CVSpilsSetGSType(mem, MODIFIED_GS); 388e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag); 3897cb94ee6SHong Zhang PetscFunctionReturn(0); 3907cb94ee6SHong Zhang } 3917cb94ee6SHong Zhang 3926fadb2cdSHong Zhang /* type of CVODE linear multistep method */ 393dfcd32c8SHong Zhang const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0}; 3946fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */ 395dfcd32c8SHong Zhang const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0}; 396a04cf4d8SBarry Smith 3977cb94ee6SHong Zhang #undef __FUNCT__ 398214bc6a2SJed Brown #define __FUNCT__ "TSSetFromOptions_Sundials" 399214bc6a2SJed Brown PetscErrorCode TSSetFromOptions_Sundials(TS ts) 4007cb94ee6SHong Zhang { 4017cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4027cb94ee6SHong Zhang PetscErrorCode ierr; 4037cb94ee6SHong Zhang int indx; 404ace3abfcSBarry Smith PetscBool flag; 4057bda3a07SJed Brown PC pc; 4067cb94ee6SHong Zhang 4077cb94ee6SHong Zhang PetscFunctionBegin; 4087cb94ee6SHong Zhang ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr); 4096fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr); 4107cb94ee6SHong Zhang if (flag) { 4116fadb2cdSHong Zhang ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr); 4127cb94ee6SHong Zhang } 4136fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr); 4147cb94ee6SHong Zhang if (flag) { 4157cb94ee6SHong Zhang ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr); 4167cb94ee6SHong Zhang } 4177cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr); 4187cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr); 419f1cd61daSJed Brown ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);CHKERRQ(ierr); 420f1cd61daSJed Brown ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);CHKERRQ(ierr); 4217cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr); 422f61b2b6aSHong Zhang ierr = PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,&flag);CHKERRQ(ierr); 423acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr); 4247cb94ee6SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4257bda3a07SJed Brown ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 4267bda3a07SJed Brown ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 4277cb94ee6SHong Zhang PetscFunctionReturn(0); 4287cb94ee6SHong Zhang } 4297cb94ee6SHong Zhang 4307cb94ee6SHong Zhang #undef __FUNCT__ 4317cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials" 4327cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer) 4337cb94ee6SHong Zhang { 4347cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4357cb94ee6SHong Zhang PetscErrorCode ierr; 4367cb94ee6SHong Zhang char *type; 4377cb94ee6SHong Zhang char atype[] = "Adams"; 4387cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 439ace3abfcSBarry Smith PetscBool iascii,isstring; 4402c823083SHong Zhang long int nsteps,its,nfevals,nlinsetups,nfails,itmp; 4412c823083SHong Zhang PetscInt qlast,qcur; 4425d47aee6SHong Zhang PetscReal hinused,hlast,hcur,tcur,tolsfac; 44342528757SHong Zhang PC pc; 4447cb94ee6SHong Zhang 4457cb94ee6SHong Zhang PetscFunctionBegin; 4467cb94ee6SHong Zhang if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;} 4477cb94ee6SHong Zhang else {type = btype;} 4487cb94ee6SHong Zhang 449251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 450251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 4517cb94ee6SHong Zhang if (iascii) { 4527cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr); 4537cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr); 4547cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr); 4557cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr); 456f61b2b6aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);CHKERRQ(ierr); 4577cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 4587cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4597cb94ee6SHong Zhang } else { 4607cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4617cb94ee6SHong Zhang } 462f1cd61daSJed Brown if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);} 463f1cd61daSJed Brown if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);} 4642c823083SHong Zhang 4655d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */ 4665d47aee6SHong Zhang ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr); 4675d47aee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr); 4682c823083SHong Zhang ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals, 4692c823083SHong Zhang &nlinsetups,&nfails,&qlast,&qcur, 4702c823083SHong Zhang &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr); 4712c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr); 4722c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr); 4732c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr); 4742c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr); 4752c823083SHong Zhang 4762c823083SHong Zhang ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr); 4772c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr); 4782c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr); 4792c823083SHong Zhang 4802c823083SHong Zhang ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */ 4812c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr); 4822c823083SHong Zhang ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr); 4832c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr); 48442528757SHong Zhang 48542528757SHong Zhang ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr); 48642528757SHong Zhang ierr = PCView(pc,viewer);CHKERRQ(ierr); 4872c823083SHong Zhang ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4882c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr); 4892c823083SHong Zhang ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr); 4902c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr); 49142528757SHong Zhang 4922c823083SHong Zhang ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4932c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr); 4942c823083SHong Zhang ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4952c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr); 4967cb94ee6SHong Zhang } else if (isstring) { 4977cb94ee6SHong Zhang ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr); 4987cb94ee6SHong Zhang } 4997cb94ee6SHong Zhang PetscFunctionReturn(0); 5007cb94ee6SHong Zhang } 5017cb94ee6SHong Zhang 5027cb94ee6SHong Zhang 5037cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 5047cb94ee6SHong Zhang EXTERN_C_BEGIN 5057cb94ee6SHong Zhang #undef __FUNCT__ 5067cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials" 5077087cfbeSBarry Smith PetscErrorCode TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type) 5087cb94ee6SHong Zhang { 5097cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5107cb94ee6SHong Zhang 5117cb94ee6SHong Zhang PetscFunctionBegin; 5127cb94ee6SHong Zhang cvode->cvode_type = type; 5137cb94ee6SHong Zhang PetscFunctionReturn(0); 5147cb94ee6SHong Zhang } 5157cb94ee6SHong Zhang EXTERN_C_END 5167cb94ee6SHong Zhang 5177cb94ee6SHong Zhang EXTERN_C_BEGIN 5187cb94ee6SHong Zhang #undef __FUNCT__ 519f61b2b6aSHong Zhang #define __FUNCT__ "TSSundialsSetMaxl_Sundials" 520f61b2b6aSHong Zhang PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl) 5217cb94ee6SHong Zhang { 5227cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5237cb94ee6SHong Zhang 5247cb94ee6SHong Zhang PetscFunctionBegin; 525f61b2b6aSHong Zhang cvode->maxl = maxl; 5267cb94ee6SHong Zhang PetscFunctionReturn(0); 5277cb94ee6SHong Zhang } 5287cb94ee6SHong Zhang EXTERN_C_END 5297cb94ee6SHong Zhang 5307cb94ee6SHong Zhang EXTERN_C_BEGIN 5317cb94ee6SHong Zhang #undef __FUNCT__ 5327cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials" 5337087cfbeSBarry Smith PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts,double tol) 5347cb94ee6SHong Zhang { 5357cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5367cb94ee6SHong Zhang 5377cb94ee6SHong Zhang PetscFunctionBegin; 5387cb94ee6SHong Zhang cvode->linear_tol = tol; 5397cb94ee6SHong Zhang PetscFunctionReturn(0); 5407cb94ee6SHong Zhang } 5417cb94ee6SHong Zhang EXTERN_C_END 5427cb94ee6SHong Zhang 5437cb94ee6SHong Zhang EXTERN_C_BEGIN 5447cb94ee6SHong Zhang #undef __FUNCT__ 5457cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials" 5467087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type) 5477cb94ee6SHong Zhang { 5487cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5497cb94ee6SHong Zhang 5507cb94ee6SHong Zhang PetscFunctionBegin; 5517cb94ee6SHong Zhang cvode->gtype = type; 5527cb94ee6SHong Zhang PetscFunctionReturn(0); 5537cb94ee6SHong Zhang } 5547cb94ee6SHong Zhang EXTERN_C_END 5557cb94ee6SHong Zhang 5567cb94ee6SHong Zhang EXTERN_C_BEGIN 5577cb94ee6SHong Zhang #undef __FUNCT__ 5587cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials" 5597087cfbeSBarry Smith PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel) 5607cb94ee6SHong Zhang { 5617cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5627cb94ee6SHong Zhang 5637cb94ee6SHong Zhang PetscFunctionBegin; 5647cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 5657cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 5667cb94ee6SHong Zhang PetscFunctionReturn(0); 5677cb94ee6SHong Zhang } 5687cb94ee6SHong Zhang EXTERN_C_END 5697cb94ee6SHong Zhang 5707cb94ee6SHong Zhang EXTERN_C_BEGIN 5717cb94ee6SHong Zhang #undef __FUNCT__ 572f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials" 5737087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt) 574f1cd61daSJed Brown { 575f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 576f1cd61daSJed Brown 577f1cd61daSJed Brown PetscFunctionBegin; 578f1cd61daSJed Brown cvode->mindt = mindt; 579f1cd61daSJed Brown PetscFunctionReturn(0); 580f1cd61daSJed Brown } 581f1cd61daSJed Brown EXTERN_C_END 582f1cd61daSJed Brown 583f1cd61daSJed Brown EXTERN_C_BEGIN 584f1cd61daSJed Brown #undef __FUNCT__ 585f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials" 5867087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt) 587f1cd61daSJed Brown { 588f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 589f1cd61daSJed Brown 590f1cd61daSJed Brown PetscFunctionBegin; 591f1cd61daSJed Brown cvode->maxdt = maxdt; 592f1cd61daSJed Brown PetscFunctionReturn(0); 593f1cd61daSJed Brown } 594f1cd61daSJed Brown EXTERN_C_END 595f1cd61daSJed Brown 596f1cd61daSJed Brown EXTERN_C_BEGIN 597f1cd61daSJed Brown #undef __FUNCT__ 5987cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials" 5997087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC_Sundials(TS ts,PC *pc) 6007cb94ee6SHong Zhang { 601dcbc6d53SSean Farley SNES snes; 602dcbc6d53SSean Farley KSP ksp; 603dcbc6d53SSean Farley PetscErrorCode ierr; 6047cb94ee6SHong Zhang 6057cb94ee6SHong Zhang PetscFunctionBegin; 606dcbc6d53SSean Farley ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 607dcbc6d53SSean Farley ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 608dcbc6d53SSean Farley ierr = KSPGetPC(ksp,pc); CHKERRQ(ierr); 6097cb94ee6SHong Zhang PetscFunctionReturn(0); 6107cb94ee6SHong Zhang } 6117cb94ee6SHong Zhang EXTERN_C_END 6127cb94ee6SHong Zhang 6137cb94ee6SHong Zhang EXTERN_C_BEGIN 6147cb94ee6SHong Zhang #undef __FUNCT__ 6157cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials" 6167087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin) 6177cb94ee6SHong Zhang { 6187cb94ee6SHong Zhang PetscFunctionBegin; 6195ef26d82SJed Brown if (nonlin) *nonlin = ts->snes_its; 6205ef26d82SJed Brown if (lin) *lin = ts->ksp_its; 6217cb94ee6SHong Zhang PetscFunctionReturn(0); 6227cb94ee6SHong Zhang } 6237cb94ee6SHong Zhang EXTERN_C_END 6247cb94ee6SHong Zhang 6257cb94ee6SHong Zhang EXTERN_C_BEGIN 6267cb94ee6SHong Zhang #undef __FUNCT__ 6272bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials" 6287087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s) 6292bfc04deSHong Zhang { 6302bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 6312bfc04deSHong Zhang 6322bfc04deSHong Zhang PetscFunctionBegin; 6332bfc04deSHong Zhang cvode->monitorstep = s; 6342bfc04deSHong Zhang PetscFunctionReturn(0); 6352bfc04deSHong Zhang } 6362bfc04deSHong Zhang EXTERN_C_END 6377cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 6387cb94ee6SHong Zhang 6397cb94ee6SHong Zhang #undef __FUNCT__ 6407cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations" 6417cb94ee6SHong Zhang /*@C 6427cb94ee6SHong Zhang TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 6437cb94ee6SHong Zhang 6447cb94ee6SHong Zhang Not Collective 6457cb94ee6SHong Zhang 6467cb94ee6SHong Zhang Input parameters: 6477cb94ee6SHong Zhang . ts - the time-step context 6487cb94ee6SHong Zhang 6497cb94ee6SHong Zhang Output Parameters: 6507cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 6517cb94ee6SHong Zhang - lin - number of linear iterations 6527cb94ee6SHong Zhang 6537cb94ee6SHong Zhang Level: advanced 6547cb94ee6SHong Zhang 6557cb94ee6SHong Zhang Notes: 6567cb94ee6SHong Zhang These return the number since the creation of the TS object 6577cb94ee6SHong Zhang 6587cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations 6597cb94ee6SHong Zhang 660f61b2b6aSHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetMaxl(), 6617cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 662f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 663a43b19c4SJed Brown TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime() 6647cb94ee6SHong Zhang 6657cb94ee6SHong Zhang @*/ 6667087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations(TS ts,int *nonlin,int *lin) 6677cb94ee6SHong Zhang { 6684ac538c5SBarry Smith PetscErrorCode ierr; 6697cb94ee6SHong Zhang 6707cb94ee6SHong Zhang PetscFunctionBegin; 6714ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr); 6727cb94ee6SHong Zhang PetscFunctionReturn(0); 6737cb94ee6SHong Zhang } 6747cb94ee6SHong Zhang 6757cb94ee6SHong Zhang #undef __FUNCT__ 6767cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType" 6777cb94ee6SHong Zhang /*@ 6787cb94ee6SHong Zhang TSSundialsSetType - Sets the method that Sundials will use for integration. 6797cb94ee6SHong Zhang 6803f9fe445SBarry Smith Logically Collective on TS 6817cb94ee6SHong Zhang 6827cb94ee6SHong Zhang Input parameters: 6837cb94ee6SHong Zhang + ts - the time-step context 6847cb94ee6SHong Zhang - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 6857cb94ee6SHong Zhang 6867cb94ee6SHong Zhang Level: intermediate 6877cb94ee6SHong Zhang 6887cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula 6897cb94ee6SHong Zhang 690f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetMaxl(), 6917cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 692f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 6937cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 694a43b19c4SJed Brown TSSetExactFinalTime() 6957cb94ee6SHong Zhang @*/ 6967087cfbeSBarry Smith PetscErrorCode TSSundialsSetType(TS ts,TSSundialsLmmType type) 6977cb94ee6SHong Zhang { 6984ac538c5SBarry Smith PetscErrorCode ierr; 6997cb94ee6SHong Zhang 7007cb94ee6SHong Zhang PetscFunctionBegin; 7014ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr); 7027cb94ee6SHong Zhang PetscFunctionReturn(0); 7037cb94ee6SHong Zhang } 7047cb94ee6SHong Zhang 7057cb94ee6SHong Zhang #undef __FUNCT__ 706f61b2b6aSHong Zhang #define __FUNCT__ "TSSundialsSetMaxl" 7077cb94ee6SHong Zhang /*@ 708f61b2b6aSHong Zhang TSSundialsSetMaxl - Sets the dimension of the Krylov space used by 7097cb94ee6SHong Zhang GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 710f61b2b6aSHong Zhang this is the maximum number of GMRES steps that will be used. 7117cb94ee6SHong Zhang 7123f9fe445SBarry Smith Logically Collective on TS 7137cb94ee6SHong Zhang 7147cb94ee6SHong Zhang Input parameters: 7157cb94ee6SHong Zhang + ts - the time-step context 716f61b2b6aSHong Zhang - maxl - number of direction vectors (the dimension of Krylov subspace). 7177cb94ee6SHong Zhang 7187cb94ee6SHong Zhang Level: advanced 7197cb94ee6SHong Zhang 720f61b2b6aSHong Zhang .keywords: GMRES 7217cb94ee6SHong Zhang 7227cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 7237cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 724f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7257cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 726a43b19c4SJed Brown TSSetExactFinalTime() 7277cb94ee6SHong Zhang 7287cb94ee6SHong Zhang @*/ 729f61b2b6aSHong Zhang PetscErrorCode TSSundialsSetMaxl(TS ts,PetscInt maxl) 7307cb94ee6SHong Zhang { 7314ac538c5SBarry Smith PetscErrorCode ierr; 7327cb94ee6SHong Zhang 7337cb94ee6SHong Zhang PetscFunctionBegin; 734f61b2b6aSHong Zhang PetscValidLogicalCollectiveInt(ts,maxl,2); 735f61b2b6aSHong Zhang ierr = PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));CHKERRQ(ierr); 7367cb94ee6SHong Zhang PetscFunctionReturn(0); 7377cb94ee6SHong Zhang } 7387cb94ee6SHong Zhang 7397cb94ee6SHong Zhang #undef __FUNCT__ 7407cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance" 7417cb94ee6SHong Zhang /*@ 7427cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 7437cb94ee6SHong Zhang system by SUNDIALS. 7447cb94ee6SHong Zhang 7453f9fe445SBarry Smith Logically Collective on TS 7467cb94ee6SHong Zhang 7477cb94ee6SHong Zhang Input parameters: 7487cb94ee6SHong Zhang + ts - the time-step context 7497cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 7507cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 7517cb94ee6SHong Zhang 7527cb94ee6SHong Zhang Level: advanced 7537cb94ee6SHong Zhang 7547cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS 7557cb94ee6SHong Zhang 756f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 7577cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 758f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7597cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 760a43b19c4SJed Brown TSSetExactFinalTime() 7617cb94ee6SHong Zhang 7627cb94ee6SHong Zhang @*/ 7637087cfbeSBarry Smith PetscErrorCode TSSundialsSetLinearTolerance(TS ts,double tol) 7647cb94ee6SHong Zhang { 7654ac538c5SBarry Smith PetscErrorCode ierr; 7667cb94ee6SHong Zhang 7677cb94ee6SHong Zhang PetscFunctionBegin; 768c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,tol,2); 7694ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr); 7707cb94ee6SHong Zhang PetscFunctionReturn(0); 7717cb94ee6SHong Zhang } 7727cb94ee6SHong Zhang 7737cb94ee6SHong Zhang #undef __FUNCT__ 7747cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType" 7757cb94ee6SHong Zhang /*@ 7767cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 7777cb94ee6SHong Zhang in GMRES method by SUNDIALS linear solver. 7787cb94ee6SHong Zhang 7793f9fe445SBarry Smith Logically Collective on TS 7807cb94ee6SHong Zhang 7817cb94ee6SHong Zhang Input parameters: 7827cb94ee6SHong Zhang + ts - the time-step context 7837cb94ee6SHong Zhang - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 7847cb94ee6SHong Zhang 7857cb94ee6SHong Zhang Level: advanced 7867cb94ee6SHong Zhang 7877cb94ee6SHong Zhang .keywords: Sundials, orthogonalization 7887cb94ee6SHong Zhang 789f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 7907cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), 791f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7927cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 793a43b19c4SJed Brown TSSetExactFinalTime() 7947cb94ee6SHong Zhang 7957cb94ee6SHong Zhang @*/ 7967087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type) 7977cb94ee6SHong Zhang { 7984ac538c5SBarry Smith PetscErrorCode ierr; 7997cb94ee6SHong Zhang 8007cb94ee6SHong Zhang PetscFunctionBegin; 8014ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr); 8027cb94ee6SHong Zhang PetscFunctionReturn(0); 8037cb94ee6SHong Zhang } 8047cb94ee6SHong Zhang 8057cb94ee6SHong Zhang #undef __FUNCT__ 8067cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance" 8077cb94ee6SHong Zhang /*@ 8087cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 8097cb94ee6SHong Zhang Sundials for error control. 8107cb94ee6SHong Zhang 8113f9fe445SBarry Smith Logically Collective on TS 8127cb94ee6SHong Zhang 8137cb94ee6SHong Zhang Input parameters: 8147cb94ee6SHong Zhang + ts - the time-step context 8157cb94ee6SHong Zhang . aabs - the absolute tolerance 8167cb94ee6SHong Zhang - rel - the relative tolerance 8177cb94ee6SHong Zhang 8187cb94ee6SHong Zhang See the Cvode/Sundials users manual for exact details on these parameters. Essentially 8197cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 8207cb94ee6SHong Zhang 8217cb94ee6SHong Zhang Level: intermediate 8227cb94ee6SHong Zhang 8237cb94ee6SHong Zhang .keywords: Sundials, tolerance 8247cb94ee6SHong Zhang 825f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(), 8267cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 827f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 8287cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 829a43b19c4SJed Brown TSSetExactFinalTime() 8307cb94ee6SHong Zhang 8317cb94ee6SHong Zhang @*/ 8327087cfbeSBarry Smith PetscErrorCode TSSundialsSetTolerance(TS ts,double aabs,double rel) 8337cb94ee6SHong Zhang { 8344ac538c5SBarry Smith PetscErrorCode ierr; 8357cb94ee6SHong Zhang 8367cb94ee6SHong Zhang PetscFunctionBegin; 8374ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr); 8387cb94ee6SHong Zhang PetscFunctionReturn(0); 8397cb94ee6SHong Zhang } 8407cb94ee6SHong Zhang 8417cb94ee6SHong Zhang #undef __FUNCT__ 8427cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC" 8437cb94ee6SHong Zhang /*@ 8447cb94ee6SHong Zhang TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 8457cb94ee6SHong Zhang 8467cb94ee6SHong Zhang Input Parameter: 8477cb94ee6SHong Zhang . ts - the time-step context 8487cb94ee6SHong Zhang 8497cb94ee6SHong Zhang Output Parameter: 8507cb94ee6SHong Zhang . pc - the preconditioner context 8517cb94ee6SHong Zhang 8527cb94ee6SHong Zhang Level: advanced 8537cb94ee6SHong Zhang 854f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 8557cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 856f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 8577cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance() 8587cb94ee6SHong Zhang @*/ 8597087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC(TS ts,PC *pc) 8607cb94ee6SHong Zhang { 8614ac538c5SBarry Smith PetscErrorCode ierr; 8627cb94ee6SHong Zhang 8637cb94ee6SHong Zhang PetscFunctionBegin; 8644ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));CHKERRQ(ierr); 8657cb94ee6SHong Zhang PetscFunctionReturn(0); 8667cb94ee6SHong Zhang } 8677cb94ee6SHong Zhang 8687cb94ee6SHong Zhang #undef __FUNCT__ 869f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep" 870f1cd61daSJed Brown /*@ 871f1cd61daSJed Brown TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 872f1cd61daSJed Brown 873f1cd61daSJed Brown Input Parameter: 874f1cd61daSJed Brown + ts - the time-step context 875f1cd61daSJed Brown - mindt - lowest time step if positive, negative to deactivate 876f1cd61daSJed Brown 877fc6b9e64SJed Brown Note: 878fc6b9e64SJed Brown Sundials will error if it is not possible to keep the estimated truncation error below 879fc6b9e64SJed Brown the tolerance set with TSSundialsSetTolerance() without going below this step size. 880fc6b9e64SJed Brown 881f1cd61daSJed Brown Level: beginner 882f1cd61daSJed Brown 883f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 884f1cd61daSJed Brown @*/ 8857087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep(TS ts,PetscReal mindt) 886f1cd61daSJed Brown { 8874ac538c5SBarry Smith PetscErrorCode ierr; 888f1cd61daSJed Brown 889f1cd61daSJed Brown PetscFunctionBegin; 8904ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr); 891f1cd61daSJed Brown PetscFunctionReturn(0); 892f1cd61daSJed Brown } 893f1cd61daSJed Brown 894f1cd61daSJed Brown #undef __FUNCT__ 895f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep" 896f1cd61daSJed Brown /*@ 897f1cd61daSJed Brown TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 898f1cd61daSJed Brown 899f1cd61daSJed Brown Input Parameter: 900f1cd61daSJed Brown + ts - the time-step context 901f1cd61daSJed Brown - maxdt - lowest time step if positive, negative to deactivate 902f1cd61daSJed Brown 903f1cd61daSJed Brown Level: beginner 904f1cd61daSJed Brown 905f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 906f1cd61daSJed Brown @*/ 9077087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt) 908f1cd61daSJed Brown { 9094ac538c5SBarry Smith PetscErrorCode ierr; 910f1cd61daSJed Brown 911f1cd61daSJed Brown PetscFunctionBegin; 9124ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 913f1cd61daSJed Brown PetscFunctionReturn(0); 914f1cd61daSJed Brown } 915f1cd61daSJed Brown 916f1cd61daSJed Brown #undef __FUNCT__ 9172bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps" 9182bfc04deSHong Zhang /*@ 9192bfc04deSHong Zhang TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 9202bfc04deSHong Zhang 9212bfc04deSHong Zhang Input Parameter: 9222bfc04deSHong Zhang + ts - the time-step context 9232bfc04deSHong Zhang - ft - PETSC_TRUE if monitor, else PETSC_FALSE 9242bfc04deSHong Zhang 9252bfc04deSHong Zhang Level: beginner 9262bfc04deSHong Zhang 927f61b2b6aSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 9282bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 929f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 9302bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 9312bfc04deSHong Zhang @*/ 9327087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps(TS ts,PetscBool ft) 9332bfc04deSHong Zhang { 9344ac538c5SBarry Smith PetscErrorCode ierr; 9352bfc04deSHong Zhang 9362bfc04deSHong Zhang PetscFunctionBegin; 9374ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr); 9382bfc04deSHong Zhang PetscFunctionReturn(0); 9392bfc04deSHong Zhang } 9407cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 9417cb94ee6SHong Zhang /*MC 94296f5712cSJed Brown TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 9437cb94ee6SHong Zhang 9447cb94ee6SHong Zhang Options Database: 9457cb94ee6SHong Zhang + -ts_sundials_type <bdf,adams> 9467cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 9477cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 9487cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 9497cb94ee6SHong Zhang . -ts_sundials_linear_tolerance <tol> 950f61b2b6aSHong Zhang . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace 95116016685SHong Zhang - -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps 95216016685SHong Zhang 9537cb94ee6SHong Zhang 9547cb94ee6SHong Zhang Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply 9557cb94ee6SHong Zhang only PETSc PC options 9567cb94ee6SHong Zhang 9577cb94ee6SHong Zhang Level: beginner 9587cb94ee6SHong Zhang 959f61b2b6aSHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(), 960a43b19c4SJed Brown TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime() 9617cb94ee6SHong Zhang 9627cb94ee6SHong Zhang M*/ 9637cb94ee6SHong Zhang EXTERN_C_BEGIN 9647cb94ee6SHong Zhang #undef __FUNCT__ 9657cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials" 9667087cfbeSBarry Smith PetscErrorCode TSCreate_Sundials(TS ts) 9677cb94ee6SHong Zhang { 9687cb94ee6SHong Zhang TS_Sundials *cvode; 9697cb94ee6SHong Zhang PetscErrorCode ierr; 97042528757SHong Zhang PC pc; 9717cb94ee6SHong Zhang 9727cb94ee6SHong Zhang PetscFunctionBegin; 973277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Sundials; 97428adb3f7SHong Zhang ts->ops->destroy = TSDestroy_Sundials; 97528adb3f7SHong Zhang ts->ops->view = TSView_Sundials; 976214bc6a2SJed Brown ts->ops->setup = TSSetUp_Sundials; 977b4eba00bSSean Farley ts->ops->step = TSStep_Sundials; 978b4eba00bSSean Farley ts->ops->interpolate = TSInterpolate_Sundials; 979214bc6a2SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Sundials; 9807cb94ee6SHong Zhang 98138f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr); 9827cb94ee6SHong Zhang ts->data = (void*)cvode; 9836fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF; 9846fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS; 985f61b2b6aSHong Zhang cvode->maxl = 5; 9867cb94ee6SHong Zhang cvode->linear_tol = .05; 9877cb94ee6SHong Zhang 988b4eba00bSSean Farley cvode->monitorstep = PETSC_TRUE; 9897cb94ee6SHong Zhang 990a07356b8SSatish Balay ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr); 991f1cd61daSJed Brown 992f1cd61daSJed Brown cvode->mindt = -1.; 993f1cd61daSJed Brown cvode->maxdt = -1.; 994f1cd61daSJed Brown 9957cb94ee6SHong Zhang /* set tolerance for Sundials */ 9967cb94ee6SHong Zhang cvode->reltol = 1e-6; 9972c823083SHong Zhang cvode->abstol = 1e-6; 9987cb94ee6SHong Zhang 99942528757SHong Zhang /* set PCNONE as default pctype */ 100042528757SHong Zhang ierr = TSSundialsGetPC_Sundials(ts,&pc);CHKERRQ(ierr); 100142528757SHong Zhang ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 100242528757SHong Zhang 1003b4eba00bSSean Farley if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE; 1004a43b19c4SJed Brown 10057cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials", 10067cb94ee6SHong Zhang TSSundialsSetType_Sundials);CHKERRQ(ierr); 1007f61b2b6aSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxl_C", 1008f61b2b6aSHong Zhang "TSSundialsSetMaxl_Sundials", 1009f61b2b6aSHong Zhang TSSundialsSetMaxl_Sundials);CHKERRQ(ierr); 10107cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C", 10117cb94ee6SHong Zhang "TSSundialsSetLinearTolerance_Sundials", 10127cb94ee6SHong Zhang TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr); 10137cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C", 10147cb94ee6SHong Zhang "TSSundialsSetGramSchmidtType_Sundials", 10157cb94ee6SHong Zhang TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr); 10167cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C", 10177cb94ee6SHong Zhang "TSSundialsSetTolerance_Sundials", 10187cb94ee6SHong Zhang TSSundialsSetTolerance_Sundials);CHKERRQ(ierr); 1019f1cd61daSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C", 1020f1cd61daSJed Brown "TSSundialsSetMinTimeStep_Sundials", 1021f1cd61daSJed Brown TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr); 1022f1cd61daSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C", 1023f1cd61daSJed Brown "TSSundialsSetMaxTimeStep_Sundials", 1024f1cd61daSJed Brown TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr); 10257cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C", 10267cb94ee6SHong Zhang "TSSundialsGetPC_Sundials", 10277cb94ee6SHong Zhang TSSundialsGetPC_Sundials);CHKERRQ(ierr); 10287cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C", 10297cb94ee6SHong Zhang "TSSundialsGetIterations_Sundials", 10307cb94ee6SHong Zhang TSSundialsGetIterations_Sundials);CHKERRQ(ierr); 10312bfc04deSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C", 10322bfc04deSHong Zhang "TSSundialsMonitorInternalSteps_Sundials", 10332bfc04deSHong Zhang TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr); 10342bfc04deSHong Zhang 10357cb94ee6SHong Zhang PetscFunctionReturn(0); 10367cb94ee6SHong Zhang } 10377cb94ee6SHong Zhang EXTERN_C_END 1038