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 */ 81c7d2463SJed Brown #include <../src/ts/impls/implicit/sundials/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 PetscScalar *y_data; 277cb94ee6SHong Zhang 287cb94ee6SHong Zhang PetscFunctionBegin; 290298fd71SBarry Smith ierr = TSGetIJacobian(ts,&J,&P,NULL,NULL);CHKERRQ(ierr); 307cb94ee6SHong Zhang y_data = (PetscScalar*) N_VGetArrayPointer(y); 317cb94ee6SHong Zhang ierr = VecPlaceArray(yy,y_data);CHKERRQ(ierr); 320679a0aeSJed Brown ierr = VecZeroEntries(yydot);CHKERRQ(ierr); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */ 330679a0aeSJed Brown /* compute the shifted Jacobian (1/gm)*I + Jrest */ 3409e82e2bSBarry Smith ierr = TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,J,P,PETSC_FALSE);CHKERRQ(ierr); 357cb94ee6SHong Zhang ierr = VecResetArray(yy);CHKERRQ(ierr); 360679a0aeSJed Brown ierr = MatScale(P,gm);CHKERRQ(ierr); /* turn into I-gm*Jrest, J is not used by Sundials */ 377cb94ee6SHong Zhang *jcurPtr = TRUE; 38dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 3909e82e2bSBarry Smith ierr = PCSetOperators(pc,J,P);CHKERRQ(ierr); 407cb94ee6SHong Zhang PetscFunctionReturn(0); 417cb94ee6SHong Zhang } 427cb94ee6SHong Zhang 437cb94ee6SHong Zhang /* 447cb94ee6SHong Zhang TSPSolve_Sundials - routine that we provide to Sundials that applies the preconditioner. 457cb94ee6SHong Zhang */ 467cb94ee6SHong Zhang #undef __FUNCT__ 477cb94ee6SHong Zhang #define __FUNCT__ "TSPSolve_Sundials" 484ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z, 494ac7836bSHong Zhang realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp) 507cb94ee6SHong Zhang { 517cb94ee6SHong Zhang TS ts = (TS) P_data; 527cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 53dcbc6d53SSean Farley PC pc; 547cb94ee6SHong Zhang Vec rr = cvode->w1,zz = cvode->w2; 557cb94ee6SHong Zhang PetscErrorCode ierr; 567cb94ee6SHong Zhang PetscScalar *r_data,*z_data; 577cb94ee6SHong Zhang 587cb94ee6SHong Zhang PetscFunctionBegin; 597cb94ee6SHong Zhang /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/ 607cb94ee6SHong Zhang r_data = (PetscScalar*) N_VGetArrayPointer(r); 617cb94ee6SHong Zhang z_data = (PetscScalar*) N_VGetArrayPointer(z); 627cb94ee6SHong Zhang ierr = VecPlaceArray(rr,r_data);CHKERRQ(ierr); 637cb94ee6SHong Zhang ierr = VecPlaceArray(zz,z_data);CHKERRQ(ierr); 644ac7836bSHong Zhang 657cb94ee6SHong Zhang /* Solve the Px=r and put the result in zz */ 66dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 677cb94ee6SHong Zhang ierr = PCApply(pc,rr,zz);CHKERRQ(ierr); 687cb94ee6SHong Zhang ierr = VecResetArray(rr);CHKERRQ(ierr); 697cb94ee6SHong Zhang ierr = VecResetArray(zz);CHKERRQ(ierr); 707cb94ee6SHong Zhang PetscFunctionReturn(0); 717cb94ee6SHong Zhang } 727cb94ee6SHong Zhang 737cb94ee6SHong Zhang /* 747cb94ee6SHong Zhang TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side. 757cb94ee6SHong Zhang */ 767cb94ee6SHong Zhang #undef __FUNCT__ 777cb94ee6SHong Zhang #define __FUNCT__ "TSFunction_Sundials" 784ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx) 797cb94ee6SHong Zhang { 807cb94ee6SHong Zhang TS ts = (TS) ctx; 815fcef5e4SJed Brown DM dm; 82942e3340SBarry Smith DMTS tsdm; 835fcef5e4SJed Brown TSIFunction ifunction; 84ce94432eSBarry Smith MPI_Comm 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; 91ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 925ff03cf7SDinesh Kaushik /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/ 937cb94ee6SHong Zhang y_data = (PetscScalar*) N_VGetArrayPointer(y); 947cb94ee6SHong Zhang ydot_data = (PetscScalar*) N_VGetArrayPointer(ydot); 954ac7836bSHong Zhang ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr); 964ac7836bSHong Zhang ierr = VecPlaceArray(yyd,ydot_data);CHKERRABORT(comm,ierr); 974ac7836bSHong Zhang 985fcef5e4SJed Brown /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */ 995fcef5e4SJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 100942e3340SBarry Smith ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 1010298fd71SBarry Smith ierr = DMTSGetIFunction(dm,&ifunction,NULL);CHKERRQ(ierr); 1025fcef5e4SJed Brown if (!ifunction) { 103bc0cc02bSJed Brown ierr = TSComputeRHSFunction(ts,t,yy,yyd);CHKERRQ(ierr); 104bc0cc02bSJed Brown } else { /* If rhsfunction is also set, this computes both parts and shifts them to the right */ 105bc0cc02bSJed Brown ierr = VecZeroEntries(yydot);CHKERRQ(ierr); 1060679a0aeSJed Brown ierr = TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE);CHKERRABORT(comm,ierr); 1070679a0aeSJed Brown ierr = VecScale(yyd,-1.);CHKERRQ(ierr); 108bc0cc02bSJed Brown } 1097cb94ee6SHong Zhang ierr = VecResetArray(yy);CHKERRABORT(comm,ierr); 1107cb94ee6SHong Zhang ierr = VecResetArray(yyd);CHKERRABORT(comm,ierr); 1114ac7836bSHong Zhang PetscFunctionReturn(0); 1127cb94ee6SHong Zhang } 1137cb94ee6SHong Zhang 1147cb94ee6SHong Zhang /* 115b4eba00bSSean Farley TSStep_Sundials - Calls Sundials to integrate the ODE. 1167cb94ee6SHong Zhang */ 1177cb94ee6SHong Zhang #undef __FUNCT__ 118b4eba00bSSean Farley #define __FUNCT__ "TSStep_Sundials" 119b4eba00bSSean Farley PetscErrorCode TSStep_Sundials(TS ts) 1207cb94ee6SHong Zhang { 1217cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 1227cb94ee6SHong Zhang PetscErrorCode ierr; 123b4eba00bSSean Farley PetscInt flag; 1249f94935aSHong Zhang long int its,nsteps; 1257cb94ee6SHong Zhang realtype t,tout; 1267cb94ee6SHong Zhang PetscScalar *y_data; 1277cb94ee6SHong Zhang void *mem; 1287cb94ee6SHong Zhang 1297cb94ee6SHong Zhang PetscFunctionBegin; 13016016685SHong Zhang mem = cvode->mem; 1317cb94ee6SHong Zhang tout = ts->max_time; 1327cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr); 1337cb94ee6SHong Zhang N_VSetArrayPointer((realtype*)y_data,cvode->y); 1340298fd71SBarry Smith ierr = VecRestoreArray(ts->vec_sol,NULL);CHKERRQ(ierr); 135186e87acSLisandro Dalcin 1363f2090d5SJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 137186e87acSLisandro Dalcin 1389be3e283SDebojyoti Ghosh /* We would like to call TSPreStep() when starting each step (including rejections), TSPreStage(), 1399be3e283SDebojyoti Ghosh * and TSPostStage() before each stage solve, but CVode does not appear to support this. */ 140bbd56ea5SKarl Rupp if (cvode->monitorstep) flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP); 141bbd56ea5SKarl Rupp else flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL); 1429f94935aSHong Zhang 1439f94935aSHong Zhang if (flag) { /* display error message */ 1449f94935aSHong Zhang switch (flag) { 1459f94935aSHong Zhang case CV_ILL_INPUT: 1469f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT"); 1479f94935aSHong Zhang break; 1489f94935aSHong Zhang case CV_TOO_CLOSE: 1499f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE"); 1509f94935aSHong Zhang break; 1513c7fefeeSJed Brown case CV_TOO_MUCH_WORK: { 1529f94935aSHong Zhang PetscReal tcur; 1539f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 1549f94935aSHong Zhang ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr); 1557c8652ddSBarry Smith 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()",(double)tcur,nsteps,ts->max_steps); 1563c7fefeeSJed Brown } break; 1579f94935aSHong Zhang case CV_TOO_MUCH_ACC: 1589f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC"); 1599f94935aSHong Zhang break; 1609f94935aSHong Zhang case CV_ERR_FAILURE: 1619f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE"); 1629f94935aSHong Zhang break; 1639f94935aSHong Zhang case CV_CONV_FAILURE: 1649f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE"); 1659f94935aSHong Zhang break; 1669f94935aSHong Zhang case CV_LINIT_FAIL: 1679f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL"); 1689f94935aSHong Zhang break; 1699f94935aSHong Zhang case CV_LSETUP_FAIL: 1709f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL"); 1719f94935aSHong Zhang break; 1729f94935aSHong Zhang case CV_LSOLVE_FAIL: 1739f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL"); 1749f94935aSHong Zhang break; 1759f94935aSHong Zhang case CV_RHSFUNC_FAIL: 1769f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL"); 1779f94935aSHong Zhang break; 1789f94935aSHong Zhang case CV_FIRST_RHSFUNC_ERR: 1799f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR"); 1809f94935aSHong Zhang break; 1819f94935aSHong Zhang case CV_REPTD_RHSFUNC_ERR: 1829f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR"); 1839f94935aSHong Zhang break; 1849f94935aSHong Zhang case CV_UNREC_RHSFUNC_ERR: 1859f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR"); 1869f94935aSHong Zhang break; 1879f94935aSHong Zhang case CV_RTFUNC_FAIL: 1889f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL"); 1899f94935aSHong Zhang break; 1909f94935aSHong Zhang default: 1919f94935aSHong Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag); 1929f94935aSHong Zhang } 1939f94935aSHong Zhang } 1949f94935aSHong Zhang 1957cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */ 1967cb94ee6SHong Zhang ierr = VecPlaceArray(cvode->w1,y_data);CHKERRQ(ierr); 1977cb94ee6SHong Zhang ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr); 1987cb94ee6SHong Zhang ierr = VecResetArray(cvode->w1);CHKERRQ(ierr); 199bc0cc02bSJed Brown ierr = VecCopy(cvode->update,ts->vec_sol);CHKERRQ(ierr); 2007cb94ee6SHong Zhang ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr); 2014ac7836bSHong Zhang ierr = CVSpilsGetNumLinIters(mem, &its); 2025ef26d82SJed Brown ts->snes_its = its; ts->ksp_its = its; 203186e87acSLisandro Dalcin 204186e87acSLisandro Dalcin ts->time_step = t - ts->ptime; 205186e87acSLisandro Dalcin ts->ptime = t; 2067cb94ee6SHong Zhang ts->steps++; 207186e87acSLisandro Dalcin 2089f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 209b4eba00bSSean Farley if (!cvode->monitorstep) ts->steps = nsteps; 210b4eba00bSSean Farley PetscFunctionReturn(0); 211b4eba00bSSean Farley } 212b4eba00bSSean Farley 213b4eba00bSSean Farley #undef __FUNCT__ 214b4eba00bSSean Farley #define __FUNCT__ "TSInterpolate_Sundials" 215b4eba00bSSean Farley static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X) 216b4eba00bSSean Farley { 217b4eba00bSSean Farley TS_Sundials *cvode = (TS_Sundials*)ts->data; 218b4eba00bSSean Farley N_Vector y; 219b4eba00bSSean Farley PetscErrorCode ierr; 220b4eba00bSSean Farley PetscScalar *x_data; 221b4eba00bSSean Farley PetscInt glosize,locsize; 222b4eba00bSSean Farley 223b4eba00bSSean Farley PetscFunctionBegin; 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); 2367cb94ee6SHong Zhang PetscFunctionReturn(0); 2377cb94ee6SHong Zhang } 2387cb94ee6SHong Zhang 2397cb94ee6SHong Zhang #undef __FUNCT__ 240277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Sundials" 241277b19d0SLisandro Dalcin PetscErrorCode TSReset_Sundials(TS ts) 242277b19d0SLisandro Dalcin { 243277b19d0SLisandro Dalcin TS_Sundials *cvode = (TS_Sundials*)ts->data; 244277b19d0SLisandro Dalcin PetscErrorCode ierr; 245277b19d0SLisandro Dalcin 246277b19d0SLisandro Dalcin PetscFunctionBegin; 2475c6b4a3dSJed Brown ierr = VecDestroy(&cvode->update);CHKERRQ(ierr); 2480679a0aeSJed Brown ierr = VecDestroy(&cvode->ydot);CHKERRQ(ierr); 2495c6b4a3dSJed Brown ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr); 2505c6b4a3dSJed Brown ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr); 251bbd56ea5SKarl Rupp if (cvode->mem) CVodeFree(&cvode->mem); 252277b19d0SLisandro Dalcin PetscFunctionReturn(0); 253277b19d0SLisandro Dalcin } 254277b19d0SLisandro Dalcin 255277b19d0SLisandro Dalcin #undef __FUNCT__ 2567cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials" 2577cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts) 2587cb94ee6SHong Zhang { 2597cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2607cb94ee6SHong Zhang PetscErrorCode ierr; 2617cb94ee6SHong Zhang 2627cb94ee6SHong Zhang PetscFunctionBegin; 263277b19d0SLisandro Dalcin ierr = TSReset_Sundials(ts);CHKERRQ(ierr); 264a07356b8SSatish Balay ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr); 265277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 266bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",NULL);CHKERRQ(ierr); 267bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",NULL);CHKERRQ(ierr); 268bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",NULL);CHKERRQ(ierr); 269bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",NULL);CHKERRQ(ierr); 270bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",NULL);CHKERRQ(ierr); 271bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",NULL);CHKERRQ(ierr); 272bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",NULL);CHKERRQ(ierr); 273bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",NULL);CHKERRQ(ierr); 274bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",NULL);CHKERRQ(ierr); 275bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",NULL);CHKERRQ(ierr); 2767cb94ee6SHong Zhang PetscFunctionReturn(0); 2777cb94ee6SHong Zhang } 2787cb94ee6SHong Zhang 2797cb94ee6SHong Zhang #undef __FUNCT__ 280214bc6a2SJed Brown #define __FUNCT__ "TSSetUp_Sundials" 281214bc6a2SJed Brown PetscErrorCode TSSetUp_Sundials(TS ts) 2827cb94ee6SHong Zhang { 2837cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2847cb94ee6SHong Zhang PetscErrorCode ierr; 28516016685SHong Zhang PetscInt glosize,locsize,i,flag; 2867cb94ee6SHong Zhang PetscScalar *y_data,*parray; 28716016685SHong Zhang void *mem; 288dcbc6d53SSean Farley PC pc; 28919fd82e9SBarry Smith PCType pctype; 290ace3abfcSBarry Smith PetscBool pcnone; 2917cb94ee6SHong Zhang 2927cb94ee6SHong Zhang PetscFunctionBegin; 293*236c45dcSShri Abhyankar if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for exact final time option 'MATCHSTEP' when using Sundials"); 294*236c45dcSShri Abhyankar 2957cb94ee6SHong Zhang /* get the vector size */ 2967cb94ee6SHong Zhang ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr); 2977cb94ee6SHong Zhang ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr); 2987cb94ee6SHong Zhang 2997cb94ee6SHong Zhang /* allocate the memory for N_Vec y */ 300a07356b8SSatish Balay cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 301e32f2f54SBarry Smith if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated"); 3027cb94ee6SHong Zhang 30328adb3f7SHong Zhang /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 3047cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr); 3057cb94ee6SHong Zhang y_data = (PetscScalar*) N_VGetArrayPointer(cvode->y); 3067cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 3070298fd71SBarry Smith ierr = VecRestoreArray(ts->vec_sol,NULL);CHKERRQ(ierr); 30842528757SHong Zhang 3097cb94ee6SHong Zhang ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr); 3100679a0aeSJed Brown ierr = VecDuplicate(ts->vec_sol,&cvode->ydot);CHKERRQ(ierr); 3113bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->update);CHKERRQ(ierr); 3123bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->ydot);CHKERRQ(ierr); 3137cb94ee6SHong Zhang 3147cb94ee6SHong Zhang /* 3157cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are 3167cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided 3177cb94ee6SHong Zhang by Sundials and set using VecPlaceArray(). 3187cb94ee6SHong Zhang */ 319ce94432eSBarry Smith ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr); 320ce94432eSBarry Smith ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr); 3213bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w1);CHKERRQ(ierr); 3223bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w2);CHKERRQ(ierr); 32316016685SHong Zhang 32416016685SHong Zhang /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 32516016685SHong Zhang mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 326e32f2f54SBarry Smith if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails"); 32716016685SHong Zhang cvode->mem = mem; 32816016685SHong Zhang 32916016685SHong Zhang /* Set the pointer to user-defined data */ 33016016685SHong Zhang flag = CVodeSetUserData(mem, ts); 331e32f2f54SBarry Smith if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails"); 33216016685SHong Zhang 333fc6b9e64SJed Brown /* Sundials may choose to use a smaller initial step, but will never use a larger step. */ 334cdbf8f93SLisandro Dalcin flag = CVodeSetInitStep(mem,(realtype)ts->time_step); 335ce94432eSBarry Smith if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetInitStep() failed"); 336f1cd61daSJed Brown if (cvode->mindt > 0) { 337f1cd61daSJed Brown flag = CVodeSetMinStep(mem,(realtype)cvode->mindt); 3389f94935aSHong Zhang if (flag) { 339ce94432eSBarry Smith if (flag == CV_MEM_NULL) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL"); 340ce94432eSBarry 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"); 341ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed"); 3429f94935aSHong Zhang } 343f1cd61daSJed Brown } 344f1cd61daSJed Brown if (cvode->maxdt > 0) { 345f1cd61daSJed Brown flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt); 346ce94432eSBarry Smith if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMaxStep() failed"); 347f1cd61daSJed Brown } 348f1cd61daSJed Brown 34916016685SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the 35016016685SHong Zhang * user's right hand side function in u'=f(t,u), the inital time T0, and 35116016685SHong Zhang * the initial dependent variable vector cvode->y */ 35216016685SHong Zhang flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y); 353bbd56ea5SKarl Rupp if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag); 35416016685SHong Zhang 3559f94935aSHong Zhang /* specifies scalar relative and absolute tolerances */ 35616016685SHong Zhang flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol); 357bbd56ea5SKarl Rupp if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag); 35816016685SHong Zhang 3599f94935aSHong Zhang /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */ 3609f94935aSHong Zhang flag = CVodeSetMaxNumSteps(mem,ts->max_steps); 36116016685SHong Zhang 36216016685SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 36316016685SHong Zhang /* setup the ode integrator with the given preconditioner */ 364dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 365dcbc6d53SSean Farley ierr = PCGetType(pc,&pctype);CHKERRQ(ierr); 366251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr); 36716016685SHong Zhang if (pcnone) { 36816016685SHong Zhang flag = CVSpgmr(mem,PREC_NONE,0); 369e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 37016016685SHong Zhang } else { 371f61b2b6aSHong Zhang flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl); 372e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 37316016685SHong Zhang 37416016685SHong Zhang /* Set preconditioner and solve routines Precond and PSolve, 37516016685SHong Zhang and the pointer to the user-defined block data */ 37616016685SHong Zhang flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials); 377e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag); 37816016685SHong Zhang } 37916016685SHong Zhang 38016016685SHong Zhang flag = CVSpilsSetGSType(mem, MODIFIED_GS); 381e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag); 3827cb94ee6SHong Zhang PetscFunctionReturn(0); 3837cb94ee6SHong Zhang } 3847cb94ee6SHong Zhang 3856fadb2cdSHong Zhang /* type of CVODE linear multistep method */ 3866a6fc655SJed Brown const char *const TSSundialsLmmTypes[] = {"","ADAMS","BDF","TSSundialsLmmType","SUNDIALS_",0}; 3876fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */ 3886a6fc655SJed Brown const char *const TSSundialsGramSchmidtTypes[] = {"","MODIFIED","CLASSICAL","TSSundialsGramSchmidtType","SUNDIALS_",0}; 389a04cf4d8SBarry Smith 3907cb94ee6SHong Zhang #undef __FUNCT__ 391214bc6a2SJed Brown #define __FUNCT__ "TSSetFromOptions_Sundials" 392214bc6a2SJed Brown PetscErrorCode TSSetFromOptions_Sundials(TS ts) 3937cb94ee6SHong Zhang { 3947cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 3957cb94ee6SHong Zhang PetscErrorCode ierr; 3967cb94ee6SHong Zhang int indx; 397ace3abfcSBarry Smith PetscBool flag; 3987bda3a07SJed Brown PC pc; 3997cb94ee6SHong Zhang 4007cb94ee6SHong Zhang PetscFunctionBegin; 4017cb94ee6SHong Zhang ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr); 4026fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr); 4037cb94ee6SHong Zhang if (flag) { 4046fadb2cdSHong Zhang ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr); 4057cb94ee6SHong Zhang } 4066fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr); 4077cb94ee6SHong Zhang if (flag) { 4087cb94ee6SHong Zhang ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr); 4097cb94ee6SHong Zhang } 4100298fd71SBarry Smith ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,NULL);CHKERRQ(ierr); 4110298fd71SBarry Smith ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,NULL);CHKERRQ(ierr); 4120298fd71SBarry Smith ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,NULL);CHKERRQ(ierr); 4130298fd71SBarry Smith ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,NULL);CHKERRQ(ierr); 41494ae4db5SBarry Smith ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,NULL);CHKERRQ(ierr); 41594ae4db5SBarry Smith ierr = PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,NULL);CHKERRQ(ierr); 4160298fd71SBarry Smith ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,NULL);CHKERRQ(ierr); 4177cb94ee6SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4187bda3a07SJed Brown ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 4197bda3a07SJed Brown ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 4207cb94ee6SHong Zhang PetscFunctionReturn(0); 4217cb94ee6SHong Zhang } 4227cb94ee6SHong Zhang 4237cb94ee6SHong Zhang #undef __FUNCT__ 4247cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials" 4257cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer) 4267cb94ee6SHong Zhang { 4277cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4287cb94ee6SHong Zhang PetscErrorCode ierr; 4297cb94ee6SHong Zhang char *type; 4307cb94ee6SHong Zhang char atype[] = "Adams"; 4317cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 432ace3abfcSBarry Smith PetscBool iascii,isstring; 4332c823083SHong Zhang long int nsteps,its,nfevals,nlinsetups,nfails,itmp; 4342c823083SHong Zhang PetscInt qlast,qcur; 4355d47aee6SHong Zhang PetscReal hinused,hlast,hcur,tcur,tolsfac; 43642528757SHong Zhang PC pc; 4377cb94ee6SHong Zhang 4387cb94ee6SHong Zhang PetscFunctionBegin; 439bbd56ea5SKarl Rupp if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype; 440bbd56ea5SKarl Rupp else type = btype; 4417cb94ee6SHong Zhang 442251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 443251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 4447cb94ee6SHong Zhang if (iascii) { 4457cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr); 4467cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr); 4477cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr); 4487cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr); 449f61b2b6aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);CHKERRQ(ierr); 4507cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 4517cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4527cb94ee6SHong Zhang } else { 4537cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4547cb94ee6SHong Zhang } 455f1cd61daSJed Brown if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);} 456f1cd61daSJed Brown if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);} 4572c823083SHong Zhang 4585d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */ 4595d47aee6SHong Zhang ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr); 4605d47aee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr); 4612c823083SHong Zhang ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals, 4622c823083SHong Zhang &nlinsetups,&nfails,&qlast,&qcur, 4632c823083SHong Zhang &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr); 4642c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr); 4652c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr); 4662c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr); 4672c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr); 4682c823083SHong Zhang 4692c823083SHong Zhang ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr); 4702c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr); 4712c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr); 4722c823083SHong Zhang 4732c823083SHong Zhang ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */ 4742c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr); 4752c823083SHong Zhang ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr); 4762c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr); 47742528757SHong Zhang 47842528757SHong Zhang ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 47942528757SHong Zhang ierr = PCView(pc,viewer);CHKERRQ(ierr); 4802c823083SHong Zhang ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4812c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr); 4822c823083SHong Zhang ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr); 4832c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr); 48442528757SHong Zhang 4852c823083SHong Zhang ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4862c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr); 4872c823083SHong Zhang ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4882c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr); 4897cb94ee6SHong Zhang } else if (isstring) { 4907cb94ee6SHong Zhang ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr); 4917cb94ee6SHong Zhang } 4927cb94ee6SHong Zhang PetscFunctionReturn(0); 4937cb94ee6SHong Zhang } 4947cb94ee6SHong Zhang 4957cb94ee6SHong Zhang 4967cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 4977cb94ee6SHong Zhang #undef __FUNCT__ 4987cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials" 4997087cfbeSBarry Smith PetscErrorCode TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type) 5007cb94ee6SHong Zhang { 5017cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5027cb94ee6SHong Zhang 5037cb94ee6SHong Zhang PetscFunctionBegin; 5047cb94ee6SHong Zhang cvode->cvode_type = type; 5057cb94ee6SHong Zhang PetscFunctionReturn(0); 5067cb94ee6SHong Zhang } 5077cb94ee6SHong Zhang 5087cb94ee6SHong Zhang #undef __FUNCT__ 509f61b2b6aSHong Zhang #define __FUNCT__ "TSSundialsSetMaxl_Sundials" 510f61b2b6aSHong Zhang PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl) 5117cb94ee6SHong Zhang { 5127cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5137cb94ee6SHong Zhang 5147cb94ee6SHong Zhang PetscFunctionBegin; 515f61b2b6aSHong Zhang cvode->maxl = maxl; 5167cb94ee6SHong Zhang PetscFunctionReturn(0); 5177cb94ee6SHong Zhang } 5187cb94ee6SHong Zhang 5197cb94ee6SHong Zhang #undef __FUNCT__ 5207cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials" 5217087cfbeSBarry Smith PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts,double tol) 5227cb94ee6SHong Zhang { 5237cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5247cb94ee6SHong Zhang 5257cb94ee6SHong Zhang PetscFunctionBegin; 5267cb94ee6SHong Zhang cvode->linear_tol = tol; 5277cb94ee6SHong Zhang PetscFunctionReturn(0); 5287cb94ee6SHong Zhang } 5297cb94ee6SHong Zhang 5307cb94ee6SHong Zhang #undef __FUNCT__ 5317cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials" 5327087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type) 5337cb94ee6SHong Zhang { 5347cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5357cb94ee6SHong Zhang 5367cb94ee6SHong Zhang PetscFunctionBegin; 5377cb94ee6SHong Zhang cvode->gtype = type; 5387cb94ee6SHong Zhang PetscFunctionReturn(0); 5397cb94ee6SHong Zhang } 5407cb94ee6SHong Zhang 5417cb94ee6SHong Zhang #undef __FUNCT__ 5427cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials" 5437087cfbeSBarry Smith PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel) 5447cb94ee6SHong Zhang { 5457cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5467cb94ee6SHong Zhang 5477cb94ee6SHong Zhang PetscFunctionBegin; 5487cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 5497cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 5507cb94ee6SHong Zhang PetscFunctionReturn(0); 5517cb94ee6SHong Zhang } 5527cb94ee6SHong Zhang 5537cb94ee6SHong Zhang #undef __FUNCT__ 554f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials" 5557087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt) 556f1cd61daSJed Brown { 557f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 558f1cd61daSJed Brown 559f1cd61daSJed Brown PetscFunctionBegin; 560f1cd61daSJed Brown cvode->mindt = mindt; 561f1cd61daSJed Brown PetscFunctionReturn(0); 562f1cd61daSJed Brown } 563f1cd61daSJed Brown 564f1cd61daSJed Brown #undef __FUNCT__ 565f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials" 5667087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt) 567f1cd61daSJed Brown { 568f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 569f1cd61daSJed Brown 570f1cd61daSJed Brown PetscFunctionBegin; 571f1cd61daSJed Brown cvode->maxdt = maxdt; 572f1cd61daSJed Brown PetscFunctionReturn(0); 573f1cd61daSJed Brown } 574f1cd61daSJed Brown #undef __FUNCT__ 5757cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials" 5767087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC_Sundials(TS ts,PC *pc) 5777cb94ee6SHong Zhang { 578dcbc6d53SSean Farley SNES snes; 579dcbc6d53SSean Farley KSP ksp; 580dcbc6d53SSean Farley PetscErrorCode ierr; 5817cb94ee6SHong Zhang 5827cb94ee6SHong Zhang PetscFunctionBegin; 583dcbc6d53SSean Farley ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 584dcbc6d53SSean Farley ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 585dcbc6d53SSean Farley ierr = KSPGetPC(ksp,pc);CHKERRQ(ierr); 5867cb94ee6SHong Zhang PetscFunctionReturn(0); 5877cb94ee6SHong Zhang } 5887cb94ee6SHong Zhang 5897cb94ee6SHong Zhang #undef __FUNCT__ 5907cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials" 5917087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin) 5927cb94ee6SHong Zhang { 5937cb94ee6SHong Zhang PetscFunctionBegin; 5945ef26d82SJed Brown if (nonlin) *nonlin = ts->snes_its; 5955ef26d82SJed Brown if (lin) *lin = ts->ksp_its; 5967cb94ee6SHong Zhang PetscFunctionReturn(0); 5977cb94ee6SHong Zhang } 5987cb94ee6SHong Zhang 5997cb94ee6SHong Zhang #undef __FUNCT__ 6002bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials" 6017087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s) 6022bfc04deSHong Zhang { 6032bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 6042bfc04deSHong Zhang 6052bfc04deSHong Zhang PetscFunctionBegin; 6062bfc04deSHong Zhang cvode->monitorstep = s; 6072bfc04deSHong Zhang PetscFunctionReturn(0); 6082bfc04deSHong Zhang } 6097cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 6107cb94ee6SHong Zhang 6117cb94ee6SHong Zhang #undef __FUNCT__ 6127cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations" 6137cb94ee6SHong Zhang /*@C 6147cb94ee6SHong Zhang TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 6157cb94ee6SHong Zhang 6167cb94ee6SHong Zhang Not Collective 6177cb94ee6SHong Zhang 6187cb94ee6SHong Zhang Input parameters: 6197cb94ee6SHong Zhang . ts - the time-step context 6207cb94ee6SHong Zhang 6217cb94ee6SHong Zhang Output Parameters: 6227cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 6237cb94ee6SHong Zhang - lin - number of linear iterations 6247cb94ee6SHong Zhang 6257cb94ee6SHong Zhang Level: advanced 6267cb94ee6SHong Zhang 6277cb94ee6SHong Zhang Notes: 6287cb94ee6SHong Zhang These return the number since the creation of the TS object 6297cb94ee6SHong Zhang 6307cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations 6317cb94ee6SHong Zhang 632f61b2b6aSHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetMaxl(), 6337cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 634f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 635a43b19c4SJed Brown TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime() 6367cb94ee6SHong Zhang 6377cb94ee6SHong Zhang @*/ 6387087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations(TS ts,int *nonlin,int *lin) 6397cb94ee6SHong Zhang { 6404ac538c5SBarry Smith PetscErrorCode ierr; 6417cb94ee6SHong Zhang 6427cb94ee6SHong Zhang PetscFunctionBegin; 6434ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr); 6447cb94ee6SHong Zhang PetscFunctionReturn(0); 6457cb94ee6SHong Zhang } 6467cb94ee6SHong Zhang 6477cb94ee6SHong Zhang #undef __FUNCT__ 6487cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType" 6497cb94ee6SHong Zhang /*@ 6507cb94ee6SHong Zhang TSSundialsSetType - Sets the method that Sundials will use for integration. 6517cb94ee6SHong Zhang 6523f9fe445SBarry Smith Logically Collective on TS 6537cb94ee6SHong Zhang 6547cb94ee6SHong Zhang Input parameters: 6557cb94ee6SHong Zhang + ts - the time-step context 6567cb94ee6SHong Zhang - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 6577cb94ee6SHong Zhang 6587cb94ee6SHong Zhang Level: intermediate 6597cb94ee6SHong Zhang 6607cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula 6617cb94ee6SHong Zhang 662f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetMaxl(), 6637cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 664f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 6657cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 666a43b19c4SJed Brown TSSetExactFinalTime() 6677cb94ee6SHong Zhang @*/ 6687087cfbeSBarry Smith PetscErrorCode TSSundialsSetType(TS ts,TSSundialsLmmType type) 6697cb94ee6SHong Zhang { 6704ac538c5SBarry Smith PetscErrorCode ierr; 6717cb94ee6SHong Zhang 6727cb94ee6SHong Zhang PetscFunctionBegin; 6734ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr); 6747cb94ee6SHong Zhang PetscFunctionReturn(0); 6757cb94ee6SHong Zhang } 6767cb94ee6SHong Zhang 6777cb94ee6SHong Zhang #undef __FUNCT__ 678f61b2b6aSHong Zhang #define __FUNCT__ "TSSundialsSetMaxl" 6797cb94ee6SHong Zhang /*@ 680f61b2b6aSHong Zhang TSSundialsSetMaxl - Sets the dimension of the Krylov space used by 6817cb94ee6SHong Zhang GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 682f61b2b6aSHong Zhang this is the maximum number of GMRES steps that will be used. 6837cb94ee6SHong Zhang 6843f9fe445SBarry Smith Logically Collective on TS 6857cb94ee6SHong Zhang 6867cb94ee6SHong Zhang Input parameters: 6877cb94ee6SHong Zhang + ts - the time-step context 688f61b2b6aSHong Zhang - maxl - number of direction vectors (the dimension of Krylov subspace). 6897cb94ee6SHong Zhang 6907cb94ee6SHong Zhang Level: advanced 6917cb94ee6SHong Zhang 692f61b2b6aSHong Zhang .keywords: GMRES 6937cb94ee6SHong Zhang 6947cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 6957cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 696f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 6977cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 698a43b19c4SJed Brown TSSetExactFinalTime() 6997cb94ee6SHong Zhang 7007cb94ee6SHong Zhang @*/ 701f61b2b6aSHong Zhang PetscErrorCode TSSundialsSetMaxl(TS ts,PetscInt maxl) 7027cb94ee6SHong Zhang { 7034ac538c5SBarry Smith PetscErrorCode ierr; 7047cb94ee6SHong Zhang 7057cb94ee6SHong Zhang PetscFunctionBegin; 706f61b2b6aSHong Zhang PetscValidLogicalCollectiveInt(ts,maxl,2); 707f61b2b6aSHong Zhang ierr = PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));CHKERRQ(ierr); 7087cb94ee6SHong Zhang PetscFunctionReturn(0); 7097cb94ee6SHong Zhang } 7107cb94ee6SHong Zhang 7117cb94ee6SHong Zhang #undef __FUNCT__ 7127cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance" 7137cb94ee6SHong Zhang /*@ 7147cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 7157cb94ee6SHong Zhang system by SUNDIALS. 7167cb94ee6SHong Zhang 7173f9fe445SBarry Smith Logically Collective on TS 7187cb94ee6SHong Zhang 7197cb94ee6SHong Zhang Input parameters: 7207cb94ee6SHong Zhang + ts - the time-step context 7217cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 7227cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 7237cb94ee6SHong Zhang 7247cb94ee6SHong Zhang Level: advanced 7257cb94ee6SHong Zhang 7267cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS 7277cb94ee6SHong Zhang 728f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 7297cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 730f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7317cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 732a43b19c4SJed Brown TSSetExactFinalTime() 7337cb94ee6SHong Zhang 7347cb94ee6SHong Zhang @*/ 7357087cfbeSBarry Smith PetscErrorCode TSSundialsSetLinearTolerance(TS ts,double tol) 7367cb94ee6SHong Zhang { 7374ac538c5SBarry Smith PetscErrorCode ierr; 7387cb94ee6SHong Zhang 7397cb94ee6SHong Zhang PetscFunctionBegin; 740c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,tol,2); 7414ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr); 7427cb94ee6SHong Zhang PetscFunctionReturn(0); 7437cb94ee6SHong Zhang } 7447cb94ee6SHong Zhang 7457cb94ee6SHong Zhang #undef __FUNCT__ 7467cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType" 7477cb94ee6SHong Zhang /*@ 7487cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 7497cb94ee6SHong Zhang in GMRES method by SUNDIALS linear solver. 7507cb94ee6SHong Zhang 7513f9fe445SBarry Smith Logically Collective on TS 7527cb94ee6SHong Zhang 7537cb94ee6SHong Zhang Input parameters: 7547cb94ee6SHong Zhang + ts - the time-step context 7557cb94ee6SHong Zhang - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 7567cb94ee6SHong Zhang 7577cb94ee6SHong Zhang Level: advanced 7587cb94ee6SHong Zhang 7597cb94ee6SHong Zhang .keywords: Sundials, orthogonalization 7607cb94ee6SHong Zhang 761f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 7627cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), 763f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7647cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 765a43b19c4SJed Brown TSSetExactFinalTime() 7667cb94ee6SHong Zhang 7677cb94ee6SHong Zhang @*/ 7687087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type) 7697cb94ee6SHong Zhang { 7704ac538c5SBarry Smith PetscErrorCode ierr; 7717cb94ee6SHong Zhang 7727cb94ee6SHong Zhang PetscFunctionBegin; 7734ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr); 7747cb94ee6SHong Zhang PetscFunctionReturn(0); 7757cb94ee6SHong Zhang } 7767cb94ee6SHong Zhang 7777cb94ee6SHong Zhang #undef __FUNCT__ 7787cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance" 7797cb94ee6SHong Zhang /*@ 7807cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 7817cb94ee6SHong Zhang Sundials for error control. 7827cb94ee6SHong Zhang 7833f9fe445SBarry Smith Logically Collective on TS 7847cb94ee6SHong Zhang 7857cb94ee6SHong Zhang Input parameters: 7867cb94ee6SHong Zhang + ts - the time-step context 7877cb94ee6SHong Zhang . aabs - the absolute tolerance 7887cb94ee6SHong Zhang - rel - the relative tolerance 7897cb94ee6SHong Zhang 7907cb94ee6SHong Zhang See the Cvode/Sundials users manual for exact details on these parameters. Essentially 7917cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 7927cb94ee6SHong Zhang 7937cb94ee6SHong Zhang Level: intermediate 7947cb94ee6SHong Zhang 7957cb94ee6SHong Zhang .keywords: Sundials, tolerance 7967cb94ee6SHong Zhang 797f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(), 7987cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 799f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 8007cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 801a43b19c4SJed Brown TSSetExactFinalTime() 8027cb94ee6SHong Zhang 8037cb94ee6SHong Zhang @*/ 8047087cfbeSBarry Smith PetscErrorCode TSSundialsSetTolerance(TS ts,double aabs,double rel) 8057cb94ee6SHong Zhang { 8064ac538c5SBarry Smith PetscErrorCode ierr; 8077cb94ee6SHong Zhang 8087cb94ee6SHong Zhang PetscFunctionBegin; 8094ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr); 8107cb94ee6SHong Zhang PetscFunctionReturn(0); 8117cb94ee6SHong Zhang } 8127cb94ee6SHong Zhang 8137cb94ee6SHong Zhang #undef __FUNCT__ 8147cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC" 8157cb94ee6SHong Zhang /*@ 8167cb94ee6SHong Zhang TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 8177cb94ee6SHong Zhang 8187cb94ee6SHong Zhang Input Parameter: 8197cb94ee6SHong Zhang . ts - the time-step context 8207cb94ee6SHong Zhang 8217cb94ee6SHong Zhang Output Parameter: 8227cb94ee6SHong Zhang . pc - the preconditioner context 8237cb94ee6SHong Zhang 8247cb94ee6SHong Zhang Level: advanced 8257cb94ee6SHong Zhang 826f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 8277cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 828f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 8297cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance() 8307cb94ee6SHong Zhang @*/ 8317087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC(TS ts,PC *pc) 8327cb94ee6SHong Zhang { 8334ac538c5SBarry Smith PetscErrorCode ierr; 8347cb94ee6SHong Zhang 8357cb94ee6SHong Zhang PetscFunctionBegin; 8364ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC*),(ts,pc));CHKERRQ(ierr); 8377cb94ee6SHong Zhang PetscFunctionReturn(0); 8387cb94ee6SHong Zhang } 8397cb94ee6SHong Zhang 8407cb94ee6SHong Zhang #undef __FUNCT__ 841f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep" 842f1cd61daSJed Brown /*@ 843f1cd61daSJed Brown TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 844f1cd61daSJed Brown 845f1cd61daSJed Brown Input Parameter: 846f1cd61daSJed Brown + ts - the time-step context 847f1cd61daSJed Brown - mindt - lowest time step if positive, negative to deactivate 848f1cd61daSJed Brown 849fc6b9e64SJed Brown Note: 850fc6b9e64SJed Brown Sundials will error if it is not possible to keep the estimated truncation error below 851fc6b9e64SJed Brown the tolerance set with TSSundialsSetTolerance() without going below this step size. 852fc6b9e64SJed Brown 853f1cd61daSJed Brown Level: beginner 854f1cd61daSJed Brown 855f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 856f1cd61daSJed Brown @*/ 8577087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep(TS ts,PetscReal mindt) 858f1cd61daSJed Brown { 8594ac538c5SBarry Smith PetscErrorCode ierr; 860f1cd61daSJed Brown 861f1cd61daSJed Brown PetscFunctionBegin; 8624ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr); 863f1cd61daSJed Brown PetscFunctionReturn(0); 864f1cd61daSJed Brown } 865f1cd61daSJed Brown 866f1cd61daSJed Brown #undef __FUNCT__ 867f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep" 868f1cd61daSJed Brown /*@ 869f1cd61daSJed Brown TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 870f1cd61daSJed Brown 871f1cd61daSJed Brown Input Parameter: 872f1cd61daSJed Brown + ts - the time-step context 873f1cd61daSJed Brown - maxdt - lowest time step if positive, negative to deactivate 874f1cd61daSJed Brown 875f1cd61daSJed Brown Level: beginner 876f1cd61daSJed Brown 877f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 878f1cd61daSJed Brown @*/ 8797087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt) 880f1cd61daSJed Brown { 8814ac538c5SBarry Smith PetscErrorCode ierr; 882f1cd61daSJed Brown 883f1cd61daSJed Brown PetscFunctionBegin; 8844ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 885f1cd61daSJed Brown PetscFunctionReturn(0); 886f1cd61daSJed Brown } 887f1cd61daSJed Brown 888f1cd61daSJed Brown #undef __FUNCT__ 8892bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps" 8902bfc04deSHong Zhang /*@ 8912bfc04deSHong Zhang TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 8922bfc04deSHong Zhang 8932bfc04deSHong Zhang Input Parameter: 8942bfc04deSHong Zhang + ts - the time-step context 8952bfc04deSHong Zhang - ft - PETSC_TRUE if monitor, else PETSC_FALSE 8962bfc04deSHong Zhang 8972bfc04deSHong Zhang Level: beginner 8982bfc04deSHong Zhang 899f61b2b6aSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 9002bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 901f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 9022bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 9032bfc04deSHong Zhang @*/ 9047087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps(TS ts,PetscBool ft) 9052bfc04deSHong Zhang { 9064ac538c5SBarry Smith PetscErrorCode ierr; 9072bfc04deSHong Zhang 9082bfc04deSHong Zhang PetscFunctionBegin; 9094ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr); 9102bfc04deSHong Zhang PetscFunctionReturn(0); 9112bfc04deSHong Zhang } 9127cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 9137cb94ee6SHong Zhang /*MC 91496f5712cSJed Brown TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 9157cb94ee6SHong Zhang 9167cb94ee6SHong Zhang Options Database: 9177cb94ee6SHong Zhang + -ts_sundials_type <bdf,adams> 9187cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 9197cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 9207cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 9217cb94ee6SHong Zhang . -ts_sundials_linear_tolerance <tol> 922f61b2b6aSHong Zhang . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace 92316016685SHong Zhang - -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps 92416016685SHong Zhang 9257cb94ee6SHong Zhang 9267cb94ee6SHong Zhang Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply 9277cb94ee6SHong Zhang only PETSc PC options 9287cb94ee6SHong Zhang 9297cb94ee6SHong Zhang Level: beginner 9307cb94ee6SHong Zhang 931f61b2b6aSHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(), 932a43b19c4SJed Brown TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime() 9337cb94ee6SHong Zhang 9347cb94ee6SHong Zhang M*/ 9357cb94ee6SHong Zhang #undef __FUNCT__ 9367cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials" 9378cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts) 9387cb94ee6SHong Zhang { 9397cb94ee6SHong Zhang TS_Sundials *cvode; 9407cb94ee6SHong Zhang PetscErrorCode ierr; 94142528757SHong Zhang PC pc; 9427cb94ee6SHong Zhang 9437cb94ee6SHong Zhang PetscFunctionBegin; 944277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Sundials; 94528adb3f7SHong Zhang ts->ops->destroy = TSDestroy_Sundials; 94628adb3f7SHong Zhang ts->ops->view = TSView_Sundials; 947214bc6a2SJed Brown ts->ops->setup = TSSetUp_Sundials; 948b4eba00bSSean Farley ts->ops->step = TSStep_Sundials; 949b4eba00bSSean Farley ts->ops->interpolate = TSInterpolate_Sundials; 950214bc6a2SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Sundials; 9517cb94ee6SHong Zhang 952b00a9115SJed Brown ierr = PetscNewLog(ts,&cvode);CHKERRQ(ierr); 953bbd56ea5SKarl Rupp 9547cb94ee6SHong Zhang ts->data = (void*)cvode; 9556fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF; 9566fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS; 957f61b2b6aSHong Zhang cvode->maxl = 5; 9587cb94ee6SHong Zhang cvode->linear_tol = .05; 959b4eba00bSSean Farley cvode->monitorstep = PETSC_TRUE; 9607cb94ee6SHong Zhang 961ce94432eSBarry Smith ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)ts),&(cvode->comm_sundials));CHKERRQ(ierr); 962f1cd61daSJed Brown 963f1cd61daSJed Brown cvode->mindt = -1.; 964f1cd61daSJed Brown cvode->maxdt = -1.; 965f1cd61daSJed Brown 9667cb94ee6SHong Zhang /* set tolerance for Sundials */ 9677cb94ee6SHong Zhang cvode->reltol = 1e-6; 9682c823083SHong Zhang cvode->abstol = 1e-6; 9697cb94ee6SHong Zhang 97042528757SHong Zhang /* set PCNONE as default pctype */ 97142528757SHong Zhang ierr = TSSundialsGetPC_Sundials(ts,&pc);CHKERRQ(ierr); 97242528757SHong Zhang ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 97342528757SHong Zhang 974bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",TSSundialsSetType_Sundials);CHKERRQ(ierr); 975bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",TSSundialsSetMaxl_Sundials);CHKERRQ(ierr); 976bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr); 977bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr); 978bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",TSSundialsSetTolerance_Sundials);CHKERRQ(ierr); 979bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr); 980bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr); 981bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",TSSundialsGetPC_Sundials);CHKERRQ(ierr); 982bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",TSSundialsGetIterations_Sundials);CHKERRQ(ierr); 983bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr); 9847cb94ee6SHong Zhang PetscFunctionReturn(0); 9857cb94ee6SHong Zhang } 986