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" 16bbd56ea5SKarl Rupp PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,booleantype jok,booleantype *jcurPtr, 17bbd56ea5SKarl Rupp realtype _gamma,void *P_data,N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3) 187cb94ee6SHong Zhang { 197cb94ee6SHong Zhang TS ts = (TS) P_data; 207cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 21dcbc6d53SSean Farley PC pc; 227cb94ee6SHong Zhang PetscErrorCode ierr; 230679a0aeSJed Brown Mat J,P; 240679a0aeSJed Brown Vec yy = cvode->w1,yydot = cvode->ydot; 250679a0aeSJed Brown PetscReal gm = (PetscReal)_gamma; 267cb94ee6SHong Zhang MatStructure str = DIFFERENT_NONZERO_PATTERN; 277cb94ee6SHong Zhang PetscScalar *y_data; 287cb94ee6SHong Zhang 297cb94ee6SHong Zhang PetscFunctionBegin; 300298fd71SBarry Smith ierr = TSGetIJacobian(ts,&J,&P,NULL,NULL);CHKERRQ(ierr); 317cb94ee6SHong Zhang y_data = (PetscScalar*) N_VGetArrayPointer(y); 327cb94ee6SHong Zhang ierr = VecPlaceArray(yy,y_data);CHKERRQ(ierr); 330679a0aeSJed Brown ierr = VecZeroEntries(yydot);CHKERRQ(ierr); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */ 340679a0aeSJed Brown /* compute the shifted Jacobian (1/gm)*I + Jrest */ 350679a0aeSJed Brown ierr = TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,&J,&P,&str,PETSC_FALSE);CHKERRQ(ierr); 367cb94ee6SHong Zhang ierr = VecResetArray(yy);CHKERRQ(ierr); 370679a0aeSJed Brown ierr = MatScale(P,gm);CHKERRQ(ierr); /* turn into I-gm*Jrest, J is not used by Sundials */ 387cb94ee6SHong Zhang *jcurPtr = TRUE; 39dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 400679a0aeSJed Brown ierr = PCSetOperators(pc,J,P,str);CHKERRQ(ierr); 417cb94ee6SHong Zhang PetscFunctionReturn(0); 427cb94ee6SHong Zhang } 437cb94ee6SHong Zhang 447cb94ee6SHong Zhang /* 457cb94ee6SHong Zhang TSPSolve_Sundials - routine that we provide to Sundials that applies the preconditioner. 467cb94ee6SHong Zhang */ 477cb94ee6SHong Zhang #undef __FUNCT__ 487cb94ee6SHong Zhang #define __FUNCT__ "TSPSolve_Sundials" 494ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z, 504ac7836bSHong Zhang realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp) 517cb94ee6SHong Zhang { 527cb94ee6SHong Zhang TS ts = (TS) P_data; 537cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 54dcbc6d53SSean Farley PC pc; 557cb94ee6SHong Zhang Vec rr = cvode->w1,zz = cvode->w2; 567cb94ee6SHong Zhang PetscErrorCode ierr; 577cb94ee6SHong Zhang PetscScalar *r_data,*z_data; 587cb94ee6SHong Zhang 597cb94ee6SHong Zhang PetscFunctionBegin; 607cb94ee6SHong Zhang /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/ 617cb94ee6SHong Zhang r_data = (PetscScalar*) N_VGetArrayPointer(r); 627cb94ee6SHong Zhang z_data = (PetscScalar*) N_VGetArrayPointer(z); 637cb94ee6SHong Zhang ierr = VecPlaceArray(rr,r_data);CHKERRQ(ierr); 647cb94ee6SHong Zhang ierr = VecPlaceArray(zz,z_data);CHKERRQ(ierr); 654ac7836bSHong Zhang 667cb94ee6SHong Zhang /* Solve the Px=r and put the result in zz */ 67dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 687cb94ee6SHong Zhang ierr = PCApply(pc,rr,zz);CHKERRQ(ierr); 697cb94ee6SHong Zhang ierr = VecResetArray(rr);CHKERRQ(ierr); 707cb94ee6SHong Zhang ierr = VecResetArray(zz);CHKERRQ(ierr); 717cb94ee6SHong Zhang PetscFunctionReturn(0); 727cb94ee6SHong Zhang } 737cb94ee6SHong Zhang 747cb94ee6SHong Zhang /* 757cb94ee6SHong Zhang TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side. 767cb94ee6SHong Zhang */ 777cb94ee6SHong Zhang #undef __FUNCT__ 787cb94ee6SHong Zhang #define __FUNCT__ "TSFunction_Sundials" 794ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx) 807cb94ee6SHong Zhang { 817cb94ee6SHong Zhang TS ts = (TS) ctx; 825fcef5e4SJed Brown DM dm; 83942e3340SBarry Smith DMTS tsdm; 845fcef5e4SJed Brown TSIFunction ifunction; 85ce94432eSBarry Smith MPI_Comm comm; 867cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 870679a0aeSJed Brown Vec yy = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot; 887cb94ee6SHong Zhang PetscScalar *y_data,*ydot_data; 897cb94ee6SHong Zhang PetscErrorCode ierr; 907cb94ee6SHong Zhang 917cb94ee6SHong Zhang PetscFunctionBegin; 92ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 935ff03cf7SDinesh Kaushik /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/ 947cb94ee6SHong Zhang y_data = (PetscScalar*) N_VGetArrayPointer(y); 957cb94ee6SHong Zhang ydot_data = (PetscScalar*) N_VGetArrayPointer(ydot); 964ac7836bSHong Zhang ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr); 974ac7836bSHong Zhang ierr = VecPlaceArray(yyd,ydot_data);CHKERRABORT(comm,ierr); 984ac7836bSHong Zhang 995fcef5e4SJed Brown /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */ 1005fcef5e4SJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 101942e3340SBarry Smith ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 1020298fd71SBarry Smith ierr = DMTSGetIFunction(dm,&ifunction,NULL);CHKERRQ(ierr); 1035fcef5e4SJed Brown if (!ifunction) { 104bc0cc02bSJed Brown ierr = TSComputeRHSFunction(ts,t,yy,yyd);CHKERRQ(ierr); 105bc0cc02bSJed Brown } else { /* If rhsfunction is also set, this computes both parts and shifts them to the right */ 106bc0cc02bSJed Brown ierr = VecZeroEntries(yydot);CHKERRQ(ierr); 1070679a0aeSJed Brown ierr = TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE);CHKERRABORT(comm,ierr); 1080679a0aeSJed Brown ierr = VecScale(yyd,-1.);CHKERRQ(ierr); 109bc0cc02bSJed Brown } 1107cb94ee6SHong Zhang ierr = VecResetArray(yy);CHKERRABORT(comm,ierr); 1117cb94ee6SHong Zhang ierr = VecResetArray(yyd);CHKERRABORT(comm,ierr); 1124ac7836bSHong Zhang PetscFunctionReturn(0); 1137cb94ee6SHong Zhang } 1147cb94ee6SHong Zhang 1157cb94ee6SHong Zhang /* 116b4eba00bSSean Farley TSStep_Sundials - Calls Sundials to integrate the ODE. 1177cb94ee6SHong Zhang */ 1187cb94ee6SHong Zhang #undef __FUNCT__ 119b4eba00bSSean Farley #define __FUNCT__ "TSStep_Sundials" 120b4eba00bSSean Farley PetscErrorCode TSStep_Sundials(TS ts) 1217cb94ee6SHong Zhang { 1227cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 1237cb94ee6SHong Zhang PetscErrorCode ierr; 124b4eba00bSSean Farley PetscInt flag; 1259f94935aSHong Zhang long int its,nsteps; 1267cb94ee6SHong Zhang realtype t,tout; 1277cb94ee6SHong Zhang PetscScalar *y_data; 1287cb94ee6SHong Zhang void *mem; 1297cb94ee6SHong Zhang 1307cb94ee6SHong Zhang PetscFunctionBegin; 13116016685SHong Zhang mem = cvode->mem; 1327cb94ee6SHong Zhang tout = ts->max_time; 1337cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr); 1347cb94ee6SHong Zhang N_VSetArrayPointer((realtype*)y_data,cvode->y); 1350298fd71SBarry Smith ierr = VecRestoreArray(ts->vec_sol,NULL);CHKERRQ(ierr); 136186e87acSLisandro Dalcin 1373f2090d5SJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 138186e87acSLisandro Dalcin 139b8123daeSJed Brown /* We would like to call TSPreStep() when starting each step (including rejections) and TSPreStage() before each 140b8123daeSJed Brown * stage solve, but CVode does not appear to support this. */ 141bbd56ea5SKarl Rupp if (cvode->monitorstep) flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP); 142bbd56ea5SKarl Rupp else flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL); 1439f94935aSHong Zhang 1449f94935aSHong Zhang if (flag) { /* display error message */ 1459f94935aSHong Zhang switch (flag) { 1469f94935aSHong Zhang case CV_ILL_INPUT: 1479f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT"); 1489f94935aSHong Zhang break; 1499f94935aSHong Zhang case CV_TOO_CLOSE: 1509f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE"); 1519f94935aSHong Zhang break; 1523c7fefeeSJed Brown case CV_TOO_MUCH_WORK: { 1539f94935aSHong Zhang PetscReal tcur; 1549f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 1559f94935aSHong Zhang ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr); 1569f94935aSHong 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); 1573c7fefeeSJed Brown } break; 1589f94935aSHong Zhang case CV_TOO_MUCH_ACC: 1599f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC"); 1609f94935aSHong Zhang break; 1619f94935aSHong Zhang case CV_ERR_FAILURE: 1629f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE"); 1639f94935aSHong Zhang break; 1649f94935aSHong Zhang case CV_CONV_FAILURE: 1659f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE"); 1669f94935aSHong Zhang break; 1679f94935aSHong Zhang case CV_LINIT_FAIL: 1689f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL"); 1699f94935aSHong Zhang break; 1709f94935aSHong Zhang case CV_LSETUP_FAIL: 1719f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL"); 1729f94935aSHong Zhang break; 1739f94935aSHong Zhang case CV_LSOLVE_FAIL: 1749f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL"); 1759f94935aSHong Zhang break; 1769f94935aSHong Zhang case CV_RHSFUNC_FAIL: 1779f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL"); 1789f94935aSHong Zhang break; 1799f94935aSHong Zhang case CV_FIRST_RHSFUNC_ERR: 1809f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR"); 1819f94935aSHong Zhang break; 1829f94935aSHong Zhang case CV_REPTD_RHSFUNC_ERR: 1839f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR"); 1849f94935aSHong Zhang break; 1859f94935aSHong Zhang case CV_UNREC_RHSFUNC_ERR: 1869f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR"); 1879f94935aSHong Zhang break; 1889f94935aSHong Zhang case CV_RTFUNC_FAIL: 1899f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL"); 1909f94935aSHong Zhang break; 1919f94935aSHong Zhang default: 1929f94935aSHong Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag); 1939f94935aSHong Zhang } 1949f94935aSHong Zhang } 1959f94935aSHong Zhang 1967cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */ 1977cb94ee6SHong Zhang ierr = VecPlaceArray(cvode->w1,y_data);CHKERRQ(ierr); 1987cb94ee6SHong Zhang ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr); 1997cb94ee6SHong Zhang ierr = VecResetArray(cvode->w1);CHKERRQ(ierr); 200bc0cc02bSJed Brown ierr = VecCopy(cvode->update,ts->vec_sol);CHKERRQ(ierr); 2017cb94ee6SHong Zhang ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr); 2024ac7836bSHong Zhang ierr = CVSpilsGetNumLinIters(mem, &its); 2035ef26d82SJed Brown ts->snes_its = its; ts->ksp_its = its; 204186e87acSLisandro Dalcin 205186e87acSLisandro Dalcin ts->time_step = t - ts->ptime; 206186e87acSLisandro Dalcin ts->ptime = t; 2077cb94ee6SHong Zhang ts->steps++; 208186e87acSLisandro Dalcin 2099f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 210b4eba00bSSean Farley if (!cvode->monitorstep) ts->steps = nsteps; 211b4eba00bSSean Farley PetscFunctionReturn(0); 212b4eba00bSSean Farley } 213b4eba00bSSean Farley 214b4eba00bSSean Farley #undef __FUNCT__ 215b4eba00bSSean Farley #define __FUNCT__ "TSInterpolate_Sundials" 216b4eba00bSSean Farley static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X) 217b4eba00bSSean Farley { 218b4eba00bSSean Farley TS_Sundials *cvode = (TS_Sundials*)ts->data; 219b4eba00bSSean Farley N_Vector y; 220b4eba00bSSean Farley PetscErrorCode ierr; 221b4eba00bSSean Farley PetscScalar *x_data; 222b4eba00bSSean Farley PetscInt glosize,locsize; 223b4eba00bSSean Farley 224b4eba00bSSean Farley PetscFunctionBegin; 225b4eba00bSSean Farley /* get the vector size */ 226b4eba00bSSean Farley ierr = VecGetSize(X,&glosize);CHKERRQ(ierr); 227b4eba00bSSean Farley ierr = VecGetLocalSize(X,&locsize);CHKERRQ(ierr); 228b4eba00bSSean Farley 229b4eba00bSSean Farley /* allocate the memory for N_Vec y */ 230b4eba00bSSean Farley y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 231b4eba00bSSean Farley if (!y) SETERRQ(PETSC_COMM_SELF,1,"Interpolated y is not allocated"); 232b4eba00bSSean Farley 233b4eba00bSSean Farley ierr = VecGetArray(X,&x_data);CHKERRQ(ierr); 234b4eba00bSSean Farley N_VSetArrayPointer((realtype*)x_data,y); 235b4eba00bSSean Farley ierr = CVodeGetDky(cvode->mem,t,0,y);CHKERRQ(ierr); 236b4eba00bSSean Farley ierr = VecRestoreArray(X,&x_data);CHKERRQ(ierr); 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); 252bbd56ea5SKarl Rupp 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); 26700de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C","",NULL);CHKERRQ(ierr); 26800de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C","",NULL);CHKERRQ(ierr); 26900de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C","",NULL);CHKERRQ(ierr); 27000de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C","",NULL);CHKERRQ(ierr); 27100de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C","",NULL);CHKERRQ(ierr); 27200de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C","",NULL);CHKERRQ(ierr); 27300de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C","",NULL);CHKERRQ(ierr); 27400de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C","",NULL);CHKERRQ(ierr); 27500de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C","",NULL);CHKERRQ(ierr); 27600de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C","",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; 29019fd82e9SBarry Smith 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]; 3060298fd71SBarry Smith ierr = VecRestoreArray(ts->vec_sol,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 */ 318ce94432eSBarry Smith ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr); 319ce94432eSBarry Smith ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),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); 334ce94432eSBarry Smith if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetInitStep() failed"); 335f1cd61daSJed Brown if (cvode->mindt > 0) { 336f1cd61daSJed Brown flag = CVodeSetMinStep(mem,(realtype)cvode->mindt); 3379f94935aSHong Zhang if (flag) { 338ce94432eSBarry Smith if (flag == CV_MEM_NULL) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL"); 339ce94432eSBarry Smith else if (flag == CV_ILL_INPUT) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size"); 340ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed"); 3419f94935aSHong Zhang } 342f1cd61daSJed Brown } 343f1cd61daSJed Brown if (cvode->maxdt > 0) { 344f1cd61daSJed Brown flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt); 345ce94432eSBarry Smith if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMaxStep() failed"); 346f1cd61daSJed Brown } 347f1cd61daSJed Brown 34816016685SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the 34916016685SHong Zhang * user's right hand side function in u'=f(t,u), the inital time T0, and 35016016685SHong Zhang * the initial dependent variable vector cvode->y */ 35116016685SHong Zhang flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y); 352bbd56ea5SKarl Rupp if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag); 35316016685SHong Zhang 3549f94935aSHong Zhang /* specifies scalar relative and absolute tolerances */ 35516016685SHong Zhang flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol); 356bbd56ea5SKarl Rupp if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag); 35716016685SHong Zhang 3589f94935aSHong Zhang /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */ 3599f94935aSHong Zhang flag = CVodeSetMaxNumSteps(mem,ts->max_steps); 36016016685SHong Zhang 36116016685SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 36216016685SHong Zhang /* setup the ode integrator with the given preconditioner */ 363dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 364dcbc6d53SSean Farley ierr = PCGetType(pc,&pctype);CHKERRQ(ierr); 365251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr); 36616016685SHong Zhang if (pcnone) { 36716016685SHong Zhang flag = CVSpgmr(mem,PREC_NONE,0); 368e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 36916016685SHong Zhang } else { 370f61b2b6aSHong Zhang flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl); 371e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 37216016685SHong Zhang 37316016685SHong Zhang /* Set preconditioner and solve routines Precond and PSolve, 37416016685SHong Zhang and the pointer to the user-defined block data */ 37516016685SHong Zhang flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials); 376e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag); 37716016685SHong Zhang } 37816016685SHong Zhang 37916016685SHong Zhang flag = CVSpilsSetGSType(mem, MODIFIED_GS); 380e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag); 3817cb94ee6SHong Zhang PetscFunctionReturn(0); 3827cb94ee6SHong Zhang } 3837cb94ee6SHong Zhang 3846fadb2cdSHong Zhang /* type of CVODE linear multistep method */ 3856a6fc655SJed Brown const char *const TSSundialsLmmTypes[] = {"","ADAMS","BDF","TSSundialsLmmType","SUNDIALS_",0}; 3866fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */ 3876a6fc655SJed Brown const char *const TSSundialsGramSchmidtTypes[] = {"","MODIFIED","CLASSICAL","TSSundialsGramSchmidtType","SUNDIALS_",0}; 388a04cf4d8SBarry Smith 3897cb94ee6SHong Zhang #undef __FUNCT__ 390214bc6a2SJed Brown #define __FUNCT__ "TSSetFromOptions_Sundials" 391214bc6a2SJed Brown PetscErrorCode TSSetFromOptions_Sundials(TS ts) 3927cb94ee6SHong Zhang { 3937cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 3947cb94ee6SHong Zhang PetscErrorCode ierr; 3957cb94ee6SHong Zhang int indx; 396ace3abfcSBarry Smith PetscBool flag; 3977bda3a07SJed Brown PC pc; 3987cb94ee6SHong Zhang 3997cb94ee6SHong Zhang PetscFunctionBegin; 4007cb94ee6SHong Zhang ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr); 4016fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr); 4027cb94ee6SHong Zhang if (flag) { 4036fadb2cdSHong Zhang ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr); 4047cb94ee6SHong Zhang } 4056fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr); 4067cb94ee6SHong Zhang if (flag) { 4077cb94ee6SHong Zhang ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr); 4087cb94ee6SHong Zhang } 4090298fd71SBarry Smith ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,NULL);CHKERRQ(ierr); 4100298fd71SBarry Smith ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,NULL);CHKERRQ(ierr); 4110298fd71SBarry Smith ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,NULL);CHKERRQ(ierr); 4120298fd71SBarry Smith ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,NULL);CHKERRQ(ierr); 4137cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr); 414f61b2b6aSHong Zhang ierr = PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,&flag);CHKERRQ(ierr); 4150298fd71SBarry Smith ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,NULL);CHKERRQ(ierr); 4167cb94ee6SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4177bda3a07SJed Brown ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 4187bda3a07SJed Brown ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 4197cb94ee6SHong Zhang PetscFunctionReturn(0); 4207cb94ee6SHong Zhang } 4217cb94ee6SHong Zhang 4227cb94ee6SHong Zhang #undef __FUNCT__ 4237cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials" 4247cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer) 4257cb94ee6SHong Zhang { 4267cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4277cb94ee6SHong Zhang PetscErrorCode ierr; 4287cb94ee6SHong Zhang char *type; 4297cb94ee6SHong Zhang char atype[] = "Adams"; 4307cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 431ace3abfcSBarry Smith PetscBool iascii,isstring; 4322c823083SHong Zhang long int nsteps,its,nfevals,nlinsetups,nfails,itmp; 4332c823083SHong Zhang PetscInt qlast,qcur; 4345d47aee6SHong Zhang PetscReal hinused,hlast,hcur,tcur,tolsfac; 43542528757SHong Zhang PC pc; 4367cb94ee6SHong Zhang 4377cb94ee6SHong Zhang PetscFunctionBegin; 438bbd56ea5SKarl Rupp if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype; 439bbd56ea5SKarl Rupp else type = btype; 4407cb94ee6SHong Zhang 441251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 442251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 4437cb94ee6SHong Zhang if (iascii) { 4447cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr); 4457cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr); 4467cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr); 4477cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr); 448f61b2b6aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);CHKERRQ(ierr); 4497cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 4507cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4517cb94ee6SHong Zhang } else { 4527cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4537cb94ee6SHong Zhang } 454f1cd61daSJed Brown if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);} 455f1cd61daSJed Brown if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);} 4562c823083SHong Zhang 4575d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */ 4585d47aee6SHong Zhang ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr); 4595d47aee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr); 4602c823083SHong Zhang ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals, 4612c823083SHong Zhang &nlinsetups,&nfails,&qlast,&qcur, 4622c823083SHong Zhang &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr); 4632c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr); 4642c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr); 4652c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr); 4662c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr); 4672c823083SHong Zhang 4682c823083SHong Zhang ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr); 4692c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr); 4702c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr); 4712c823083SHong Zhang 4722c823083SHong Zhang ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */ 4732c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr); 4742c823083SHong Zhang ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr); 4752c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr); 47642528757SHong Zhang 47742528757SHong Zhang ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 47842528757SHong Zhang ierr = PCView(pc,viewer);CHKERRQ(ierr); 4792c823083SHong Zhang ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4802c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr); 4812c823083SHong Zhang ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr); 4822c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr); 48342528757SHong Zhang 4842c823083SHong Zhang ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4852c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr); 4862c823083SHong Zhang ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4872c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr); 4887cb94ee6SHong Zhang } else if (isstring) { 4897cb94ee6SHong Zhang ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr); 4907cb94ee6SHong Zhang } 4917cb94ee6SHong Zhang PetscFunctionReturn(0); 4927cb94ee6SHong Zhang } 4937cb94ee6SHong Zhang 4947cb94ee6SHong Zhang 4957cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 4967cb94ee6SHong Zhang #undef __FUNCT__ 4977cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials" 4987087cfbeSBarry Smith PetscErrorCode TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type) 4997cb94ee6SHong Zhang { 5007cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5017cb94ee6SHong Zhang 5027cb94ee6SHong Zhang PetscFunctionBegin; 5037cb94ee6SHong Zhang cvode->cvode_type = type; 5047cb94ee6SHong Zhang PetscFunctionReturn(0); 5057cb94ee6SHong Zhang } 5067cb94ee6SHong Zhang 5077cb94ee6SHong Zhang #undef __FUNCT__ 508f61b2b6aSHong Zhang #define __FUNCT__ "TSSundialsSetMaxl_Sundials" 509f61b2b6aSHong Zhang PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl) 5107cb94ee6SHong Zhang { 5117cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5127cb94ee6SHong Zhang 5137cb94ee6SHong Zhang PetscFunctionBegin; 514f61b2b6aSHong Zhang cvode->maxl = maxl; 5157cb94ee6SHong Zhang PetscFunctionReturn(0); 5167cb94ee6SHong Zhang } 5177cb94ee6SHong Zhang 5187cb94ee6SHong Zhang #undef __FUNCT__ 5197cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials" 5207087cfbeSBarry Smith PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts,double tol) 5217cb94ee6SHong Zhang { 5227cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5237cb94ee6SHong Zhang 5247cb94ee6SHong Zhang PetscFunctionBegin; 5257cb94ee6SHong Zhang cvode->linear_tol = tol; 5267cb94ee6SHong Zhang PetscFunctionReturn(0); 5277cb94ee6SHong Zhang } 5287cb94ee6SHong Zhang 5297cb94ee6SHong Zhang #undef __FUNCT__ 5307cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials" 5317087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type) 5327cb94ee6SHong Zhang { 5337cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5347cb94ee6SHong Zhang 5357cb94ee6SHong Zhang PetscFunctionBegin; 5367cb94ee6SHong Zhang cvode->gtype = type; 5377cb94ee6SHong Zhang PetscFunctionReturn(0); 5387cb94ee6SHong Zhang } 5397cb94ee6SHong Zhang 5407cb94ee6SHong Zhang #undef __FUNCT__ 5417cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials" 5427087cfbeSBarry Smith PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel) 5437cb94ee6SHong Zhang { 5447cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5457cb94ee6SHong Zhang 5467cb94ee6SHong Zhang PetscFunctionBegin; 5477cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 5487cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 5497cb94ee6SHong Zhang PetscFunctionReturn(0); 5507cb94ee6SHong Zhang } 5517cb94ee6SHong Zhang 5527cb94ee6SHong Zhang #undef __FUNCT__ 553f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials" 5547087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt) 555f1cd61daSJed Brown { 556f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 557f1cd61daSJed Brown 558f1cd61daSJed Brown PetscFunctionBegin; 559f1cd61daSJed Brown cvode->mindt = mindt; 560f1cd61daSJed Brown PetscFunctionReturn(0); 561f1cd61daSJed Brown } 562f1cd61daSJed Brown 563f1cd61daSJed Brown #undef __FUNCT__ 564f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials" 5657087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt) 566f1cd61daSJed Brown { 567f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 568f1cd61daSJed Brown 569f1cd61daSJed Brown PetscFunctionBegin; 570f1cd61daSJed Brown cvode->maxdt = maxdt; 571f1cd61daSJed Brown PetscFunctionReturn(0); 572f1cd61daSJed Brown } 573f1cd61daSJed Brown #undef __FUNCT__ 5747cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials" 5757087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC_Sundials(TS ts,PC *pc) 5767cb94ee6SHong Zhang { 577dcbc6d53SSean Farley SNES snes; 578dcbc6d53SSean Farley KSP ksp; 579dcbc6d53SSean Farley PetscErrorCode ierr; 5807cb94ee6SHong Zhang 5817cb94ee6SHong Zhang PetscFunctionBegin; 582dcbc6d53SSean Farley ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 583dcbc6d53SSean Farley ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 584dcbc6d53SSean Farley ierr = KSPGetPC(ksp,pc);CHKERRQ(ierr); 5857cb94ee6SHong Zhang PetscFunctionReturn(0); 5867cb94ee6SHong Zhang } 5877cb94ee6SHong Zhang 5887cb94ee6SHong Zhang #undef __FUNCT__ 5897cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials" 5907087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin) 5917cb94ee6SHong Zhang { 5927cb94ee6SHong Zhang PetscFunctionBegin; 5935ef26d82SJed Brown if (nonlin) *nonlin = ts->snes_its; 5945ef26d82SJed Brown if (lin) *lin = ts->ksp_its; 5957cb94ee6SHong Zhang PetscFunctionReturn(0); 5967cb94ee6SHong Zhang } 5977cb94ee6SHong Zhang 5987cb94ee6SHong Zhang #undef __FUNCT__ 5992bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials" 6007087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s) 6012bfc04deSHong Zhang { 6022bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 6032bfc04deSHong Zhang 6042bfc04deSHong Zhang PetscFunctionBegin; 6052bfc04deSHong Zhang cvode->monitorstep = s; 6062bfc04deSHong Zhang PetscFunctionReturn(0); 6072bfc04deSHong Zhang } 6087cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 6097cb94ee6SHong Zhang 6107cb94ee6SHong Zhang #undef __FUNCT__ 6117cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations" 6127cb94ee6SHong Zhang /*@C 6137cb94ee6SHong Zhang TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 6147cb94ee6SHong Zhang 6157cb94ee6SHong Zhang Not Collective 6167cb94ee6SHong Zhang 6177cb94ee6SHong Zhang Input parameters: 6187cb94ee6SHong Zhang . ts - the time-step context 6197cb94ee6SHong Zhang 6207cb94ee6SHong Zhang Output Parameters: 6217cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 6227cb94ee6SHong Zhang - lin - number of linear iterations 6237cb94ee6SHong Zhang 6247cb94ee6SHong Zhang Level: advanced 6257cb94ee6SHong Zhang 6267cb94ee6SHong Zhang Notes: 6277cb94ee6SHong Zhang These return the number since the creation of the TS object 6287cb94ee6SHong Zhang 6297cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations 6307cb94ee6SHong Zhang 631f61b2b6aSHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetMaxl(), 6327cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 633f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 634a43b19c4SJed Brown TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime() 6357cb94ee6SHong Zhang 6367cb94ee6SHong Zhang @*/ 6377087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations(TS ts,int *nonlin,int *lin) 6387cb94ee6SHong Zhang { 6394ac538c5SBarry Smith PetscErrorCode ierr; 6407cb94ee6SHong Zhang 6417cb94ee6SHong Zhang PetscFunctionBegin; 6424ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr); 6437cb94ee6SHong Zhang PetscFunctionReturn(0); 6447cb94ee6SHong Zhang } 6457cb94ee6SHong Zhang 6467cb94ee6SHong Zhang #undef __FUNCT__ 6477cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType" 6487cb94ee6SHong Zhang /*@ 6497cb94ee6SHong Zhang TSSundialsSetType - Sets the method that Sundials will use for integration. 6507cb94ee6SHong Zhang 6513f9fe445SBarry Smith Logically Collective on TS 6527cb94ee6SHong Zhang 6537cb94ee6SHong Zhang Input parameters: 6547cb94ee6SHong Zhang + ts - the time-step context 6557cb94ee6SHong Zhang - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 6567cb94ee6SHong Zhang 6577cb94ee6SHong Zhang Level: intermediate 6587cb94ee6SHong Zhang 6597cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula 6607cb94ee6SHong Zhang 661f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetMaxl(), 6627cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 663f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 6647cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 665a43b19c4SJed Brown TSSetExactFinalTime() 6667cb94ee6SHong Zhang @*/ 6677087cfbeSBarry Smith PetscErrorCode TSSundialsSetType(TS ts,TSSundialsLmmType type) 6687cb94ee6SHong Zhang { 6694ac538c5SBarry Smith PetscErrorCode ierr; 6707cb94ee6SHong Zhang 6717cb94ee6SHong Zhang PetscFunctionBegin; 6724ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr); 6737cb94ee6SHong Zhang PetscFunctionReturn(0); 6747cb94ee6SHong Zhang } 6757cb94ee6SHong Zhang 6767cb94ee6SHong Zhang #undef __FUNCT__ 677f61b2b6aSHong Zhang #define __FUNCT__ "TSSundialsSetMaxl" 6787cb94ee6SHong Zhang /*@ 679f61b2b6aSHong Zhang TSSundialsSetMaxl - Sets the dimension of the Krylov space used by 6807cb94ee6SHong Zhang GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 681f61b2b6aSHong Zhang this is the maximum number of GMRES steps that will be used. 6827cb94ee6SHong Zhang 6833f9fe445SBarry Smith Logically Collective on TS 6847cb94ee6SHong Zhang 6857cb94ee6SHong Zhang Input parameters: 6867cb94ee6SHong Zhang + ts - the time-step context 687f61b2b6aSHong Zhang - maxl - number of direction vectors (the dimension of Krylov subspace). 6887cb94ee6SHong Zhang 6897cb94ee6SHong Zhang Level: advanced 6907cb94ee6SHong Zhang 691f61b2b6aSHong Zhang .keywords: GMRES 6927cb94ee6SHong Zhang 6937cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 6947cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 695f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 6967cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 697a43b19c4SJed Brown TSSetExactFinalTime() 6987cb94ee6SHong Zhang 6997cb94ee6SHong Zhang @*/ 700f61b2b6aSHong Zhang PetscErrorCode TSSundialsSetMaxl(TS ts,PetscInt maxl) 7017cb94ee6SHong Zhang { 7024ac538c5SBarry Smith PetscErrorCode ierr; 7037cb94ee6SHong Zhang 7047cb94ee6SHong Zhang PetscFunctionBegin; 705f61b2b6aSHong Zhang PetscValidLogicalCollectiveInt(ts,maxl,2); 706f61b2b6aSHong Zhang ierr = PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));CHKERRQ(ierr); 7077cb94ee6SHong Zhang PetscFunctionReturn(0); 7087cb94ee6SHong Zhang } 7097cb94ee6SHong Zhang 7107cb94ee6SHong Zhang #undef __FUNCT__ 7117cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance" 7127cb94ee6SHong Zhang /*@ 7137cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 7147cb94ee6SHong Zhang system by SUNDIALS. 7157cb94ee6SHong Zhang 7163f9fe445SBarry Smith Logically Collective on TS 7177cb94ee6SHong Zhang 7187cb94ee6SHong Zhang Input parameters: 7197cb94ee6SHong Zhang + ts - the time-step context 7207cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 7217cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 7227cb94ee6SHong Zhang 7237cb94ee6SHong Zhang Level: advanced 7247cb94ee6SHong Zhang 7257cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS 7267cb94ee6SHong Zhang 727f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 7287cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 729f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7307cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 731a43b19c4SJed Brown TSSetExactFinalTime() 7327cb94ee6SHong Zhang 7337cb94ee6SHong Zhang @*/ 7347087cfbeSBarry Smith PetscErrorCode TSSundialsSetLinearTolerance(TS ts,double tol) 7357cb94ee6SHong Zhang { 7364ac538c5SBarry Smith PetscErrorCode ierr; 7377cb94ee6SHong Zhang 7387cb94ee6SHong Zhang PetscFunctionBegin; 739c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,tol,2); 7404ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr); 7417cb94ee6SHong Zhang PetscFunctionReturn(0); 7427cb94ee6SHong Zhang } 7437cb94ee6SHong Zhang 7447cb94ee6SHong Zhang #undef __FUNCT__ 7457cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType" 7467cb94ee6SHong Zhang /*@ 7477cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 7487cb94ee6SHong Zhang in GMRES method by SUNDIALS linear solver. 7497cb94ee6SHong Zhang 7503f9fe445SBarry Smith Logically Collective on TS 7517cb94ee6SHong Zhang 7527cb94ee6SHong Zhang Input parameters: 7537cb94ee6SHong Zhang + ts - the time-step context 7547cb94ee6SHong Zhang - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 7557cb94ee6SHong Zhang 7567cb94ee6SHong Zhang Level: advanced 7577cb94ee6SHong Zhang 7587cb94ee6SHong Zhang .keywords: Sundials, orthogonalization 7597cb94ee6SHong Zhang 760f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 7617cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), 762f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7637cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 764a43b19c4SJed Brown TSSetExactFinalTime() 7657cb94ee6SHong Zhang 7667cb94ee6SHong Zhang @*/ 7677087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type) 7687cb94ee6SHong Zhang { 7694ac538c5SBarry Smith PetscErrorCode ierr; 7707cb94ee6SHong Zhang 7717cb94ee6SHong Zhang PetscFunctionBegin; 7724ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr); 7737cb94ee6SHong Zhang PetscFunctionReturn(0); 7747cb94ee6SHong Zhang } 7757cb94ee6SHong Zhang 7767cb94ee6SHong Zhang #undef __FUNCT__ 7777cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance" 7787cb94ee6SHong Zhang /*@ 7797cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 7807cb94ee6SHong Zhang Sundials for error control. 7817cb94ee6SHong Zhang 7823f9fe445SBarry Smith Logically Collective on TS 7837cb94ee6SHong Zhang 7847cb94ee6SHong Zhang Input parameters: 7857cb94ee6SHong Zhang + ts - the time-step context 7867cb94ee6SHong Zhang . aabs - the absolute tolerance 7877cb94ee6SHong Zhang - rel - the relative tolerance 7887cb94ee6SHong Zhang 7897cb94ee6SHong Zhang See the Cvode/Sundials users manual for exact details on these parameters. Essentially 7907cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 7917cb94ee6SHong Zhang 7927cb94ee6SHong Zhang Level: intermediate 7937cb94ee6SHong Zhang 7947cb94ee6SHong Zhang .keywords: Sundials, tolerance 7957cb94ee6SHong Zhang 796f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(), 7977cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 798f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7997cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 800a43b19c4SJed Brown TSSetExactFinalTime() 8017cb94ee6SHong Zhang 8027cb94ee6SHong Zhang @*/ 8037087cfbeSBarry Smith PetscErrorCode TSSundialsSetTolerance(TS ts,double aabs,double rel) 8047cb94ee6SHong Zhang { 8054ac538c5SBarry Smith PetscErrorCode ierr; 8067cb94ee6SHong Zhang 8077cb94ee6SHong Zhang PetscFunctionBegin; 8084ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr); 8097cb94ee6SHong Zhang PetscFunctionReturn(0); 8107cb94ee6SHong Zhang } 8117cb94ee6SHong Zhang 8127cb94ee6SHong Zhang #undef __FUNCT__ 8137cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC" 8147cb94ee6SHong Zhang /*@ 8157cb94ee6SHong Zhang TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 8167cb94ee6SHong Zhang 8177cb94ee6SHong Zhang Input Parameter: 8187cb94ee6SHong Zhang . ts - the time-step context 8197cb94ee6SHong Zhang 8207cb94ee6SHong Zhang Output Parameter: 8217cb94ee6SHong Zhang . pc - the preconditioner context 8227cb94ee6SHong Zhang 8237cb94ee6SHong Zhang Level: advanced 8247cb94ee6SHong Zhang 825f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 8267cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 827f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 8287cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance() 8297cb94ee6SHong Zhang @*/ 8307087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC(TS ts,PC *pc) 8317cb94ee6SHong Zhang { 8324ac538c5SBarry Smith PetscErrorCode ierr; 8337cb94ee6SHong Zhang 8347cb94ee6SHong Zhang PetscFunctionBegin; 8354ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC*),(ts,pc));CHKERRQ(ierr); 8367cb94ee6SHong Zhang PetscFunctionReturn(0); 8377cb94ee6SHong Zhang } 8387cb94ee6SHong Zhang 8397cb94ee6SHong Zhang #undef __FUNCT__ 840f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep" 841f1cd61daSJed Brown /*@ 842f1cd61daSJed Brown TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 843f1cd61daSJed Brown 844f1cd61daSJed Brown Input Parameter: 845f1cd61daSJed Brown + ts - the time-step context 846f1cd61daSJed Brown - mindt - lowest time step if positive, negative to deactivate 847f1cd61daSJed Brown 848fc6b9e64SJed Brown Note: 849fc6b9e64SJed Brown Sundials will error if it is not possible to keep the estimated truncation error below 850fc6b9e64SJed Brown the tolerance set with TSSundialsSetTolerance() without going below this step size. 851fc6b9e64SJed Brown 852f1cd61daSJed Brown Level: beginner 853f1cd61daSJed Brown 854f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 855f1cd61daSJed Brown @*/ 8567087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep(TS ts,PetscReal mindt) 857f1cd61daSJed Brown { 8584ac538c5SBarry Smith PetscErrorCode ierr; 859f1cd61daSJed Brown 860f1cd61daSJed Brown PetscFunctionBegin; 8614ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr); 862f1cd61daSJed Brown PetscFunctionReturn(0); 863f1cd61daSJed Brown } 864f1cd61daSJed Brown 865f1cd61daSJed Brown #undef __FUNCT__ 866f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep" 867f1cd61daSJed Brown /*@ 868f1cd61daSJed Brown TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 869f1cd61daSJed Brown 870f1cd61daSJed Brown Input Parameter: 871f1cd61daSJed Brown + ts - the time-step context 872f1cd61daSJed Brown - maxdt - lowest time step if positive, negative to deactivate 873f1cd61daSJed Brown 874f1cd61daSJed Brown Level: beginner 875f1cd61daSJed Brown 876f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 877f1cd61daSJed Brown @*/ 8787087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt) 879f1cd61daSJed Brown { 8804ac538c5SBarry Smith PetscErrorCode ierr; 881f1cd61daSJed Brown 882f1cd61daSJed Brown PetscFunctionBegin; 8834ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 884f1cd61daSJed Brown PetscFunctionReturn(0); 885f1cd61daSJed Brown } 886f1cd61daSJed Brown 887f1cd61daSJed Brown #undef __FUNCT__ 8882bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps" 8892bfc04deSHong Zhang /*@ 8902bfc04deSHong Zhang TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 8912bfc04deSHong Zhang 8922bfc04deSHong Zhang Input Parameter: 8932bfc04deSHong Zhang + ts - the time-step context 8942bfc04deSHong Zhang - ft - PETSC_TRUE if monitor, else PETSC_FALSE 8952bfc04deSHong Zhang 8962bfc04deSHong Zhang Level: beginner 8972bfc04deSHong Zhang 898f61b2b6aSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 8992bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 900f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 9012bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 9022bfc04deSHong Zhang @*/ 9037087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps(TS ts,PetscBool ft) 9042bfc04deSHong Zhang { 9054ac538c5SBarry Smith PetscErrorCode ierr; 9062bfc04deSHong Zhang 9072bfc04deSHong Zhang PetscFunctionBegin; 9084ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr); 9092bfc04deSHong Zhang PetscFunctionReturn(0); 9102bfc04deSHong Zhang } 9117cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 9127cb94ee6SHong Zhang /*MC 91396f5712cSJed Brown TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 9147cb94ee6SHong Zhang 9157cb94ee6SHong Zhang Options Database: 9167cb94ee6SHong Zhang + -ts_sundials_type <bdf,adams> 9177cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 9187cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 9197cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 9207cb94ee6SHong Zhang . -ts_sundials_linear_tolerance <tol> 921f61b2b6aSHong Zhang . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace 92216016685SHong Zhang - -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps 92316016685SHong Zhang 9247cb94ee6SHong Zhang 9257cb94ee6SHong Zhang Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply 9267cb94ee6SHong Zhang only PETSc PC options 9277cb94ee6SHong Zhang 9287cb94ee6SHong Zhang Level: beginner 9297cb94ee6SHong Zhang 930f61b2b6aSHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(), 931a43b19c4SJed Brown TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime() 9327cb94ee6SHong Zhang 9337cb94ee6SHong Zhang M*/ 9347cb94ee6SHong Zhang #undef __FUNCT__ 9357cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials" 936*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts) 9377cb94ee6SHong Zhang { 9387cb94ee6SHong Zhang TS_Sundials *cvode; 9397cb94ee6SHong Zhang PetscErrorCode ierr; 94042528757SHong Zhang PC pc; 9417cb94ee6SHong Zhang 9427cb94ee6SHong Zhang PetscFunctionBegin; 943277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Sundials; 94428adb3f7SHong Zhang ts->ops->destroy = TSDestroy_Sundials; 94528adb3f7SHong Zhang ts->ops->view = TSView_Sundials; 946214bc6a2SJed Brown ts->ops->setup = TSSetUp_Sundials; 947b4eba00bSSean Farley ts->ops->step = TSStep_Sundials; 948b4eba00bSSean Farley ts->ops->interpolate = TSInterpolate_Sundials; 949214bc6a2SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Sundials; 9507cb94ee6SHong Zhang 95138f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr); 952bbd56ea5SKarl Rupp 9537cb94ee6SHong Zhang ts->data = (void*)cvode; 9546fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF; 9556fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS; 956f61b2b6aSHong Zhang cvode->maxl = 5; 9577cb94ee6SHong Zhang cvode->linear_tol = .05; 958b4eba00bSSean Farley cvode->monitorstep = PETSC_TRUE; 9597cb94ee6SHong Zhang 960ce94432eSBarry Smith ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)ts),&(cvode->comm_sundials));CHKERRQ(ierr); 961f1cd61daSJed Brown 962f1cd61daSJed Brown cvode->mindt = -1.; 963f1cd61daSJed Brown cvode->maxdt = -1.; 964f1cd61daSJed Brown 9657cb94ee6SHong Zhang /* set tolerance for Sundials */ 9667cb94ee6SHong Zhang cvode->reltol = 1e-6; 9672c823083SHong Zhang cvode->abstol = 1e-6; 9687cb94ee6SHong Zhang 96942528757SHong Zhang /* set PCNONE as default pctype */ 97042528757SHong Zhang ierr = TSSundialsGetPC_Sundials(ts,&pc);CHKERRQ(ierr); 97142528757SHong Zhang ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 97242528757SHong Zhang 97349354f04SShri Abhyankar if (ts->exact_final_time != TS_EXACTFINALTIME_STEPOVER) ts->exact_final_time = TS_EXACTFINALTIME_STEPOVER; 974a43b19c4SJed Brown 97500de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",TSSundialsSetType_Sundials);CHKERRQ(ierr); 97600de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C","TSSundialsSetMaxl_Sundials",TSSundialsSetMaxl_Sundials);CHKERRQ(ierr); 97700de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C","TSSundialsSetLinearTolerance_Sundials",TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr); 97800de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C","TSSundialsSetGramSchmidtType_Sundials",TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr); 97900de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C","TSSundialsSetTolerance_Sundials",TSSundialsSetTolerance_Sundials);CHKERRQ(ierr); 98000de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C","TSSundialsSetMinTimeStep_Sundials",TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr); 98100de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C","TSSundialsSetMaxTimeStep_Sundials",TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr); 98200de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C","TSSundialsGetPC_Sundials",TSSundialsGetPC_Sundials);CHKERRQ(ierr); 98300de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C","TSSundialsGetIterations_Sundials",TSSundialsGetIterations_Sundials);CHKERRQ(ierr); 98400de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C","TSSundialsMonitorInternalSteps_Sundials",TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr); 9857cb94ee6SHong Zhang PetscFunctionReturn(0); 9867cb94ee6SHong Zhang } 987