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 */ 14bbd56ea5SKarl Rupp PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,booleantype jok,booleantype *jcurPtr, 15bbd56ea5SKarl Rupp realtype _gamma,void *P_data,N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3) 167cb94ee6SHong Zhang { 177cb94ee6SHong Zhang TS ts = (TS) P_data; 187cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 19dcbc6d53SSean Farley PC pc; 200679a0aeSJed Brown Mat J,P; 210679a0aeSJed Brown Vec yy = cvode->w1,yydot = cvode->ydot; 220679a0aeSJed Brown PetscReal gm = (PetscReal)_gamma; 237cb94ee6SHong Zhang PetscScalar *y_data; 247cb94ee6SHong Zhang 257cb94ee6SHong Zhang PetscFunctionBegin; 269566063dSJacob Faibussowitsch PetscCall(TSGetIJacobian(ts,&J,&P,NULL,NULL)); 277cb94ee6SHong Zhang y_data = (PetscScalar*) N_VGetArrayPointer(y); 289566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(yy,y_data)); 299566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(yydot)); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */ 300679a0aeSJed Brown /* compute the shifted Jacobian (1/gm)*I + Jrest */ 319566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,J,P,PETSC_FALSE)); 329566063dSJacob Faibussowitsch PetscCall(VecResetArray(yy)); 339566063dSJacob Faibussowitsch PetscCall(MatScale(P,gm)); /* turn into I-gm*Jrest, J is not used by Sundials */ 347cb94ee6SHong Zhang *jcurPtr = TRUE; 359566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts,&pc)); 369566063dSJacob Faibussowitsch PetscCall(PCSetOperators(pc,J,P)); 377cb94ee6SHong Zhang PetscFunctionReturn(0); 387cb94ee6SHong Zhang } 397cb94ee6SHong Zhang 407cb94ee6SHong Zhang /* 417cb94ee6SHong Zhang TSPSolve_Sundials - routine that we provide to Sundials that applies the preconditioner. 427cb94ee6SHong Zhang */ 434ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z, 444ac7836bSHong Zhang realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp) 457cb94ee6SHong Zhang { 467cb94ee6SHong Zhang TS ts = (TS) P_data; 477cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 48dcbc6d53SSean Farley PC pc; 497cb94ee6SHong Zhang Vec rr = cvode->w1,zz = cvode->w2; 507cb94ee6SHong Zhang PetscScalar *r_data,*z_data; 517cb94ee6SHong Zhang 527cb94ee6SHong Zhang PetscFunctionBegin; 537cb94ee6SHong Zhang /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/ 547cb94ee6SHong Zhang r_data = (PetscScalar*) N_VGetArrayPointer(r); 557cb94ee6SHong Zhang z_data = (PetscScalar*) N_VGetArrayPointer(z); 569566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(rr,r_data)); 579566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(zz,z_data)); 584ac7836bSHong Zhang 597cb94ee6SHong Zhang /* Solve the Px=r and put the result in zz */ 609566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts,&pc)); 619566063dSJacob Faibussowitsch PetscCall(PCApply(pc,rr,zz)); 629566063dSJacob Faibussowitsch PetscCall(VecResetArray(rr)); 639566063dSJacob Faibussowitsch PetscCall(VecResetArray(zz)); 647cb94ee6SHong Zhang PetscFunctionReturn(0); 657cb94ee6SHong Zhang } 667cb94ee6SHong Zhang 677cb94ee6SHong Zhang /* 687cb94ee6SHong Zhang TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side. 697cb94ee6SHong Zhang */ 704ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx) 717cb94ee6SHong Zhang { 727cb94ee6SHong Zhang TS ts = (TS) ctx; 735fcef5e4SJed Brown DM dm; 74942e3340SBarry Smith DMTS tsdm; 755fcef5e4SJed Brown TSIFunction ifunction; 76ce94432eSBarry Smith MPI_Comm comm; 777cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 780679a0aeSJed Brown Vec yy = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot; 797cb94ee6SHong Zhang PetscScalar *y_data,*ydot_data; 807cb94ee6SHong Zhang 817cb94ee6SHong Zhang PetscFunctionBegin; 829566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts,&comm)); 835ff03cf7SDinesh Kaushik /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/ 847cb94ee6SHong Zhang y_data = (PetscScalar*) N_VGetArrayPointer(y); 857cb94ee6SHong Zhang ydot_data = (PetscScalar*) N_VGetArrayPointer(ydot); 869566063dSJacob Faibussowitsch PetscCallAbort(comm,VecPlaceArray(yy,y_data)); 879566063dSJacob Faibussowitsch PetscCallAbort(comm,VecPlaceArray(yyd,ydot_data)); 884ac7836bSHong Zhang 895fcef5e4SJed Brown /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */ 909566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 919566063dSJacob Faibussowitsch PetscCall(DMGetDMTS(dm,&tsdm)); 929566063dSJacob Faibussowitsch PetscCall(DMTSGetIFunction(dm,&ifunction,NULL)); 935fcef5e4SJed Brown if (!ifunction) { 949566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts,t,yy,yyd)); 95bc0cc02bSJed Brown } else { /* If rhsfunction is also set, this computes both parts and shifts them to the right */ 969566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(yydot)); 979566063dSJacob Faibussowitsch PetscCallAbort(comm,TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE)); 989566063dSJacob Faibussowitsch PetscCall(VecScale(yyd,-1.)); 99bc0cc02bSJed Brown } 1009566063dSJacob Faibussowitsch PetscCallAbort(comm,VecResetArray(yy)); 1019566063dSJacob Faibussowitsch PetscCallAbort(comm,VecResetArray(yyd)); 1024ac7836bSHong Zhang PetscFunctionReturn(0); 1037cb94ee6SHong Zhang } 1047cb94ee6SHong Zhang 1057cb94ee6SHong Zhang /* 106b4eba00bSSean Farley TSStep_Sundials - Calls Sundials to integrate the ODE. 1077cb94ee6SHong Zhang */ 108b4eba00bSSean Farley PetscErrorCode TSStep_Sundials(TS ts) 1097cb94ee6SHong Zhang { 1107cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 111b4eba00bSSean Farley PetscInt flag; 112be5899b3SLisandro Dalcin long int nits,lits,nsteps; 1137cb94ee6SHong Zhang realtype t,tout; 1147cb94ee6SHong Zhang PetscScalar *y_data; 1157cb94ee6SHong Zhang void *mem; 1167cb94ee6SHong Zhang 1177cb94ee6SHong Zhang PetscFunctionBegin; 11816016685SHong Zhang mem = cvode->mem; 1197cb94ee6SHong Zhang tout = ts->max_time; 1209566063dSJacob Faibussowitsch PetscCall(VecGetArray(ts->vec_sol,&y_data)); 1217cb94ee6SHong Zhang N_VSetArrayPointer((realtype*)y_data,cvode->y); 1229566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ts->vec_sol,NULL)); 123186e87acSLisandro Dalcin 1249687d888SLisandro Dalcin /* We would like to TSPreStage() and TSPostStage() 1259687d888SLisandro Dalcin * before each stage solve but CVode does not appear to support this. */ 126be5899b3SLisandro Dalcin if (cvode->monitorstep) 127be5899b3SLisandro Dalcin flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP); 128be5899b3SLisandro Dalcin else 129be5899b3SLisandro Dalcin flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL); 1309f94935aSHong Zhang 1319f94935aSHong Zhang if (flag) { /* display error message */ 1329f94935aSHong Zhang switch (flag) { 1339f94935aSHong Zhang case CV_ILL_INPUT: 1349f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT"); 1359f94935aSHong Zhang break; 1369f94935aSHong Zhang case CV_TOO_CLOSE: 1379f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE"); 1389f94935aSHong Zhang break; 1393c7fefeeSJed Brown case CV_TOO_MUCH_WORK: { 1409f94935aSHong Zhang PetscReal tcur; 1419566063dSJacob Faibussowitsch PetscCall(CVodeGetNumSteps(mem,&nsteps)); 1429566063dSJacob Faibussowitsch PetscCall(CVodeGetCurrentTime(mem,&tcur)); 14363a3b9bcSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_WORK. At t=%g, nsteps %ld exceeds maxstep %" PetscInt_FMT ". Increase '-ts_max_steps <>' or modify TSSetMaxSteps()",(double)tcur,nsteps,ts->max_steps); 1443c7fefeeSJed Brown } break; 1459f94935aSHong Zhang case CV_TOO_MUCH_ACC: 1469f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC"); 1479f94935aSHong Zhang break; 1489f94935aSHong Zhang case CV_ERR_FAILURE: 1499f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE"); 1509f94935aSHong Zhang break; 1519f94935aSHong Zhang case CV_CONV_FAILURE: 1529f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE"); 1539f94935aSHong Zhang break; 1549f94935aSHong Zhang case CV_LINIT_FAIL: 1559f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL"); 1569f94935aSHong Zhang break; 1579f94935aSHong Zhang case CV_LSETUP_FAIL: 1589f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL"); 1599f94935aSHong Zhang break; 1609f94935aSHong Zhang case CV_LSOLVE_FAIL: 1619f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL"); 1629f94935aSHong Zhang break; 1639f94935aSHong Zhang case CV_RHSFUNC_FAIL: 1649f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL"); 1659f94935aSHong Zhang break; 1669f94935aSHong Zhang case CV_FIRST_RHSFUNC_ERR: 1679f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR"); 1689f94935aSHong Zhang break; 1699f94935aSHong Zhang case CV_REPTD_RHSFUNC_ERR: 1709f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR"); 1719f94935aSHong Zhang break; 1729f94935aSHong Zhang case CV_UNREC_RHSFUNC_ERR: 1739f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR"); 1749f94935aSHong Zhang break; 1759f94935aSHong Zhang case CV_RTFUNC_FAIL: 1769f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL"); 1779f94935aSHong Zhang break; 1789f94935aSHong Zhang default: 17998921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag); 1809f94935aSHong Zhang } 1819f94935aSHong Zhang } 1829f94935aSHong Zhang 183be5899b3SLisandro Dalcin /* log inner nonlinear and linear iterations */ 1849566063dSJacob Faibussowitsch PetscCall(CVodeGetNumNonlinSolvIters(mem,&nits)); 1859566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumLinIters(mem,&lits)); 186be5899b3SLisandro Dalcin ts->snes_its += nits; ts->ksp_its = lits; 187be5899b3SLisandro Dalcin 1887cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */ 1899566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(cvode->w1,y_data)); 1909566063dSJacob Faibussowitsch PetscCall(VecCopy(cvode->w1,cvode->update)); 1919566063dSJacob Faibussowitsch PetscCall(VecResetArray(cvode->w1)); 1929566063dSJacob Faibussowitsch PetscCall(VecCopy(cvode->update,ts->vec_sol)); 193186e87acSLisandro Dalcin 194186e87acSLisandro Dalcin ts->time_step = t - ts->ptime; 195186e87acSLisandro Dalcin ts->ptime = t; 196186e87acSLisandro Dalcin 1979566063dSJacob Faibussowitsch PetscCall(CVodeGetNumSteps(mem,&nsteps)); 1982808aa04SLisandro Dalcin if (!cvode->monitorstep) ts->steps += nsteps - 1; /* TSStep() increments the step counter by one */ 199b4eba00bSSean Farley PetscFunctionReturn(0); 200b4eba00bSSean Farley } 201b4eba00bSSean Farley 202b4eba00bSSean Farley static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X) 203b4eba00bSSean Farley { 204b4eba00bSSean Farley TS_Sundials *cvode = (TS_Sundials*)ts->data; 205b4eba00bSSean Farley N_Vector y; 206b4eba00bSSean Farley PetscScalar *x_data; 207b4eba00bSSean Farley PetscInt glosize,locsize; 208b4eba00bSSean Farley 209b4eba00bSSean Farley PetscFunctionBegin; 210b4eba00bSSean Farley /* get the vector size */ 2119566063dSJacob Faibussowitsch PetscCall(VecGetSize(X,&glosize)); 2129566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X,&locsize)); 2139566063dSJacob Faibussowitsch PetscCall(VecGetArray(X,&x_data)); 214b4eba00bSSean Farley 2154d78c9e0SClaas Abert /* Initialize N_Vec y with x_data */ 216918687b8SPatrick Sanan if (cvode->use_dense) { 217918687b8SPatrick Sanan PetscMPIInt size; 218918687b8SPatrick Sanan 2199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 2203c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"TSSUNDIALS only supports a dense solve in the serial case"); 221918687b8SPatrick Sanan y = N_VMake_Serial(locsize,(realtype*)x_data); 222918687b8SPatrick Sanan } else { 2234d78c9e0SClaas Abert y = N_VMake_Parallel(cvode->comm_sundials,locsize,glosize,(realtype*)x_data); 224918687b8SPatrick Sanan } 225918687b8SPatrick Sanan 2263c633725SBarry Smith PetscCheck(y,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Interpolated y is not allocated"); 227b4eba00bSSean Farley 2289566063dSJacob Faibussowitsch PetscCall(CVodeGetDky(cvode->mem,t,0,y)); 2299566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X,&x_data)); 2307cb94ee6SHong Zhang PetscFunctionReturn(0); 2317cb94ee6SHong Zhang } 2327cb94ee6SHong Zhang 233277b19d0SLisandro Dalcin PetscErrorCode TSReset_Sundials(TS ts) 234277b19d0SLisandro Dalcin { 235277b19d0SLisandro Dalcin TS_Sundials *cvode = (TS_Sundials*)ts->data; 236277b19d0SLisandro Dalcin 237277b19d0SLisandro Dalcin PetscFunctionBegin; 2389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->update)); 2399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->ydot)); 2409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->w1)); 2419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->w2)); 242bbd56ea5SKarl Rupp if (cvode->mem) CVodeFree(&cvode->mem); 243277b19d0SLisandro Dalcin PetscFunctionReturn(0); 244277b19d0SLisandro Dalcin } 245277b19d0SLisandro Dalcin 2467cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts) 2477cb94ee6SHong Zhang { 2487cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2497cb94ee6SHong Zhang 2507cb94ee6SHong Zhang PetscFunctionBegin; 2519566063dSJacob Faibussowitsch PetscCall(TSReset_Sundials(ts)); 2529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&(cvode->comm_sundials))); 2539566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 2549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",NULL)); 2559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",NULL)); 2569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",NULL)); 2579566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",NULL)); 2589566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",NULL)); 2599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",NULL)); 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",NULL)); 2619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",NULL)); 2629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",NULL)); 2639566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",NULL)); 264*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetUseDense_C",NULL)); 2657cb94ee6SHong Zhang PetscFunctionReturn(0); 2667cb94ee6SHong Zhang } 2677cb94ee6SHong Zhang 268214bc6a2SJed Brown PetscErrorCode TSSetUp_Sundials(TS ts) 2697cb94ee6SHong Zhang { 2707cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 27116016685SHong Zhang PetscInt glosize,locsize,i,flag; 2727cb94ee6SHong Zhang PetscScalar *y_data,*parray; 27316016685SHong Zhang void *mem; 274dcbc6d53SSean Farley PC pc; 27519fd82e9SBarry Smith PCType pctype; 276ace3abfcSBarry Smith PetscBool pcnone; 2777cb94ee6SHong Zhang 2787cb94ee6SHong Zhang PetscFunctionBegin; 27908401ef6SPierre Jolivet PetscCheck(ts->exact_final_time != TS_EXACTFINALTIME_MATCHSTEP,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for exact final time option 'MATCHSTEP' when using Sundials"); 280236c45dcSShri Abhyankar 2817cb94ee6SHong Zhang /* get the vector size */ 2829566063dSJacob Faibussowitsch PetscCall(VecGetSize(ts->vec_sol,&glosize)); 2839566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ts->vec_sol,&locsize)); 2847cb94ee6SHong Zhang 2857cb94ee6SHong Zhang /* allocate the memory for N_Vec y */ 286918687b8SPatrick Sanan if (cvode->use_dense) { 287918687b8SPatrick Sanan PetscMPIInt size; 288918687b8SPatrick Sanan 2899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 2903c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"TSSUNDIALS only supports a dense solve in the serial case"); 291918687b8SPatrick Sanan cvode->y = N_VNew_Serial(locsize); 292918687b8SPatrick Sanan } else { 293a07356b8SSatish Balay cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 294918687b8SPatrick Sanan } 2953c633725SBarry Smith PetscCheck(cvode->y,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cvode->y is not allocated"); 2967cb94ee6SHong Zhang 29728adb3f7SHong Zhang /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 2989566063dSJacob Faibussowitsch PetscCall(VecGetArray(ts->vec_sol,&parray)); 2997cb94ee6SHong Zhang y_data = (PetscScalar*) N_VGetArrayPointer(cvode->y); 3007cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 3019566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ts->vec_sol,NULL)); 30242528757SHong Zhang 3039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&cvode->update)); 3049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&cvode->ydot)); 3059566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->update)); 3069566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->ydot)); 3077cb94ee6SHong Zhang 3087cb94ee6SHong Zhang /* 3097cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are 3107cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided 3117cb94ee6SHong Zhang by Sundials and set using VecPlaceArray(). 3127cb94ee6SHong Zhang */ 3139566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,NULL,&cvode->w1)); 3149566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,NULL,&cvode->w2)); 3159566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w1)); 3169566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w2)); 31716016685SHong Zhang 31816016685SHong Zhang /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 31916016685SHong Zhang mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 3203c633725SBarry Smith PetscCheck(mem,PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails"); 32116016685SHong Zhang cvode->mem = mem; 32216016685SHong Zhang 32316016685SHong Zhang /* Set the pointer to user-defined data */ 32416016685SHong Zhang flag = CVodeSetUserData(mem, ts); 32528b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails"); 32616016685SHong Zhang 327fc6b9e64SJed Brown /* Sundials may choose to use a smaller initial step, but will never use a larger step. */ 328cdbf8f93SLisandro Dalcin flag = CVodeSetInitStep(mem,(realtype)ts->time_step); 32928b400f6SJacob Faibussowitsch PetscCheck(!flag,PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetInitStep() failed"); 330f1cd61daSJed Brown if (cvode->mindt > 0) { 331f1cd61daSJed Brown flag = CVodeSetMinStep(mem,(realtype)cvode->mindt); 3329f94935aSHong Zhang if (flag) { 33308401ef6SPierre Jolivet PetscCheck(flag != CV_MEM_NULL,PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL"); 33408401ef6SPierre Jolivet else PetscCheck(flag != CV_ILL_INPUT,PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size"); 335ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed"); 3369f94935aSHong Zhang } 337f1cd61daSJed Brown } 338f1cd61daSJed Brown if (cvode->maxdt > 0) { 339f1cd61daSJed Brown flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt); 34028b400f6SJacob Faibussowitsch PetscCheck(!flag,PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMaxStep() failed"); 341f1cd61daSJed Brown } 342f1cd61daSJed Brown 34316016685SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the 344a5b23f4aSJose E. Roman * user's right hand side function in u'=f(t,u), the initial time T0, and 34516016685SHong Zhang * the initial dependent variable vector cvode->y */ 34616016685SHong Zhang flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y); 34728b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag); 34816016685SHong Zhang 3499f94935aSHong Zhang /* specifies scalar relative and absolute tolerances */ 35016016685SHong Zhang flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol); 35128b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag); 35216016685SHong Zhang 353c4e80e11SFlorian /* Specify max order of BDF / ADAMS method */ 354c4e80e11SFlorian if (cvode->maxord != PETSC_DEFAULT) { 355c4e80e11SFlorian flag = CVodeSetMaxOrd(mem,cvode->maxord); 35628b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetMaxOrd() fails, flag %d",flag); 357c4e80e11SFlorian } 358c4e80e11SFlorian 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); 36128b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetMaxNumSteps() fails, flag %d",flag); 36216016685SHong Zhang 363918687b8SPatrick Sanan if (cvode->use_dense) { 364918687b8SPatrick Sanan /* call CVDense to use a dense linear solver. */ 365918687b8SPatrick Sanan PetscMPIInt size; 366918687b8SPatrick Sanan 3679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 3683c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"TSSUNDIALS only supports a dense solve in the serial case"); 369918687b8SPatrick Sanan flag = CVDense(mem,locsize); 37028b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVDense() fails, flag %d",flag); 371918687b8SPatrick Sanan } else { 37216016685SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 37316016685SHong Zhang /* setup the ode integrator with the given preconditioner */ 3749566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts,&pc)); 3759566063dSJacob Faibussowitsch PetscCall(PCGetType(pc,&pctype)); 3769566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone)); 37716016685SHong Zhang if (pcnone) { 37816016685SHong Zhang flag = CVSpgmr(mem,PREC_NONE,0); 37928b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 38016016685SHong Zhang } else { 381f61b2b6aSHong Zhang flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl); 38228b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 38316016685SHong Zhang 38416016685SHong Zhang /* Set preconditioner and solve routines Precond and PSolve, 38516016685SHong Zhang and the pointer to the user-defined block data */ 38616016685SHong Zhang flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials); 38728b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag); 38816016685SHong Zhang } 389918687b8SPatrick Sanan } 3907cb94ee6SHong Zhang PetscFunctionReturn(0); 3917cb94ee6SHong Zhang } 3927cb94ee6SHong Zhang 3936fadb2cdSHong Zhang /* type of CVODE linear multistep method */ 394c793f718SLisandro Dalcin const char *const TSSundialsLmmTypes[] = {"","ADAMS","BDF","TSSundialsLmmType","SUNDIALS_",NULL}; 3956fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */ 396c793f718SLisandro Dalcin const char *const TSSundialsGramSchmidtTypes[] = {"","MODIFIED","CLASSICAL","TSSundialsGramSchmidtType","SUNDIALS_",NULL}; 397a04cf4d8SBarry Smith 3984416b707SBarry Smith PetscErrorCode TSSetFromOptions_Sundials(PetscOptionItems *PetscOptionsObject,TS ts) 3997cb94ee6SHong Zhang { 4007cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4017cb94ee6SHong Zhang int indx; 402ace3abfcSBarry Smith PetscBool flag; 4037bda3a07SJed Brown PC pc; 4047cb94ee6SHong Zhang 4057cb94ee6SHong Zhang PetscFunctionBegin; 406d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"SUNDIALS ODE solver options"); 4079566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag)); 4087cb94ee6SHong Zhang if (flag) { 4099566063dSJacob Faibussowitsch PetscCall(TSSundialsSetType(ts,(TSSundialsLmmType)indx)); 4107cb94ee6SHong Zhang } 4119566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag)); 4127cb94ee6SHong Zhang if (flag) { 4139566063dSJacob Faibussowitsch PetscCall(TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx)); 4147cb94ee6SHong Zhang } 4159566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,NULL)); 4169566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,NULL)); 4179566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,NULL)); 4189566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,NULL)); 4199566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,NULL)); 4209566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_sundials_maxord","Max Order for BDF/Adams method","TSSundialsSetMaxOrd",cvode->maxord,&cvode->maxord,NULL)); 4219566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,NULL)); 4229566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internal steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,NULL)); 4239566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_sundials_use_dense","Use dense internal solver in SUNDIALS (serial only)","TSSundialsSetUseDense",cvode->use_dense,&cvode->use_dense,NULL)); 424d0609cedSBarry Smith PetscOptionsHeadEnd(); 4259566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts,&pc)); 4269566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions(pc)); 4277cb94ee6SHong Zhang PetscFunctionReturn(0); 4287cb94ee6SHong Zhang } 4297cb94ee6SHong Zhang 4307cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer) 4317cb94ee6SHong Zhang { 4327cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4337cb94ee6SHong Zhang char *type; 4347cb94ee6SHong Zhang char atype[] = "Adams"; 4357cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 436ace3abfcSBarry Smith PetscBool iascii,isstring; 4372c823083SHong Zhang long int nsteps,its,nfevals,nlinsetups,nfails,itmp; 4382c823083SHong Zhang PetscInt qlast,qcur; 4395d47aee6SHong Zhang PetscReal hinused,hlast,hcur,tcur,tolsfac; 44042528757SHong Zhang PC pc; 4417cb94ee6SHong Zhang 4427cb94ee6SHong Zhang PetscFunctionBegin; 443bbd56ea5SKarl Rupp if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype; 444bbd56ea5SKarl Rupp else type = btype; 4457cb94ee6SHong Zhang 4469566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 4479566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring)); 4487cb94ee6SHong Zhang if (iascii) { 4499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials integrator does not use SNES!\n")); 4509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials integrator type %s\n",type)); 45163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials maxord %" PetscInt_FMT "\n",cvode->maxord)); 4529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",(double)cvode->abstol,(double)cvode->reltol)); 453918687b8SPatrick Sanan if (cvode->use_dense) { 4549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials integrator using a dense linear solve\n")); 455918687b8SPatrick Sanan } else { 4569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",(double)cvode->linear_tol)); 45763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %" PetscInt_FMT "\n",cvode->maxl)); 4587cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 4599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n")); 4607cb94ee6SHong Zhang } else { 4619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n")); 4627cb94ee6SHong Zhang } 463918687b8SPatrick Sanan } 4649566063dSJacob Faibussowitsch if (cvode->mindt > 0) PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",(double)cvode->mindt)); 4659566063dSJacob Faibussowitsch if (cvode->maxdt > 0) PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",(double)cvode->maxdt)); 4662c823083SHong Zhang 4675d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */ 4689566063dSJacob Faibussowitsch PetscCall(CVodeGetTolScaleFactor(cvode->mem,&tolsfac)); 4699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac)); 4709566063dSJacob Faibussowitsch PetscCall(CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,&nlinsetups,&nfails,&qlast,&qcur,&hinused,&hlast,&hcur,&tcur)); 47163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %ld\n",nsteps)); 47263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %ld\n",nfevals)); 47363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %ld\n",nlinsetups)); 47463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %ld\n",nfails)); 4752c823083SHong Zhang 4769566063dSJacob Faibussowitsch PetscCall(CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails)); 47763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %ld\n",its)); 47863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %ld\n",nfails)); 479918687b8SPatrick Sanan if (!cvode->use_dense) { 4809566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumLinIters(cvode->mem, &its)); /* its = no. of calls to TSPrecond_Sundials() */ 48163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %ld\n",its)); 4829566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumConvFails(cvode->mem,&itmp)); 48363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %ld\n",itmp)); 48442528757SHong Zhang 4859566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts,&pc)); 4869566063dSJacob Faibussowitsch PetscCall(PCView(pc,viewer)); 4879566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumPrecEvals(cvode->mem,&itmp)); 48863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %ld\n",itmp)); 4899566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumPrecSolves(cvode->mem,&itmp)); 49063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %ld\n",itmp)); 491918687b8SPatrick Sanan } 4929566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumJtimesEvals(cvode->mem,&itmp)); 49363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %ld\n",itmp)); 4949566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumRhsEvals(cvode->mem,&itmp)); 49563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %ld\n",itmp)); 4967cb94ee6SHong Zhang } else if (isstring) { 4979566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer,"Sundials type %s",type)); 4987cb94ee6SHong Zhang } 4997cb94ee6SHong Zhang PetscFunctionReturn(0); 5007cb94ee6SHong Zhang } 5017cb94ee6SHong Zhang 5027cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 5037087cfbeSBarry Smith PetscErrorCode TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type) 5047cb94ee6SHong Zhang { 5057cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5067cb94ee6SHong Zhang 5077cb94ee6SHong Zhang PetscFunctionBegin; 5087cb94ee6SHong Zhang cvode->cvode_type = type; 5097cb94ee6SHong Zhang PetscFunctionReturn(0); 5107cb94ee6SHong Zhang } 5117cb94ee6SHong Zhang 512f61b2b6aSHong Zhang PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl) 5137cb94ee6SHong Zhang { 5147cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5157cb94ee6SHong Zhang 5167cb94ee6SHong Zhang PetscFunctionBegin; 517f61b2b6aSHong Zhang cvode->maxl = maxl; 5187cb94ee6SHong Zhang PetscFunctionReturn(0); 5197cb94ee6SHong Zhang } 5207cb94ee6SHong Zhang 521590ae059SPatrick Sanan PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts,PetscReal 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 5307087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type) 5317cb94ee6SHong Zhang { 5327cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5337cb94ee6SHong Zhang 5347cb94ee6SHong Zhang PetscFunctionBegin; 5357cb94ee6SHong Zhang cvode->gtype = type; 5367cb94ee6SHong Zhang PetscFunctionReturn(0); 5377cb94ee6SHong Zhang } 5387cb94ee6SHong Zhang 539590ae059SPatrick Sanan PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts,PetscReal aabs,PetscReal rel) 5407cb94ee6SHong Zhang { 5417cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5427cb94ee6SHong Zhang 5437cb94ee6SHong Zhang PetscFunctionBegin; 5447cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 5457cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 5467cb94ee6SHong Zhang PetscFunctionReturn(0); 5477cb94ee6SHong Zhang } 5487cb94ee6SHong Zhang 5497087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt) 550f1cd61daSJed Brown { 551f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 552f1cd61daSJed Brown 553f1cd61daSJed Brown PetscFunctionBegin; 554f1cd61daSJed Brown cvode->mindt = mindt; 555f1cd61daSJed Brown PetscFunctionReturn(0); 556f1cd61daSJed Brown } 557f1cd61daSJed Brown 5587087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt) 559f1cd61daSJed Brown { 560f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 561f1cd61daSJed Brown 562f1cd61daSJed Brown PetscFunctionBegin; 563f1cd61daSJed Brown cvode->maxdt = maxdt; 564f1cd61daSJed Brown PetscFunctionReturn(0); 565f1cd61daSJed Brown } 566918687b8SPatrick Sanan 567918687b8SPatrick Sanan PetscErrorCode TSSundialsSetUseDense_Sundials(TS ts,PetscBool use_dense) 568918687b8SPatrick Sanan { 569918687b8SPatrick Sanan TS_Sundials *cvode = (TS_Sundials*)ts->data; 570918687b8SPatrick Sanan 571918687b8SPatrick Sanan PetscFunctionBegin; 572918687b8SPatrick Sanan cvode->use_dense = use_dense; 573918687b8SPatrick Sanan PetscFunctionReturn(0); 574918687b8SPatrick Sanan } 575918687b8SPatrick Sanan 5767087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC_Sundials(TS ts,PC *pc) 5777cb94ee6SHong Zhang { 578dcbc6d53SSean Farley SNES snes; 579dcbc6d53SSean Farley KSP ksp; 5807cb94ee6SHong Zhang 5817cb94ee6SHong Zhang PetscFunctionBegin; 5829566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts,&snes)); 5839566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes,&ksp)); 5849566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,pc)); 5857cb94ee6SHong Zhang PetscFunctionReturn(0); 5867cb94ee6SHong Zhang } 5877cb94ee6SHong Zhang 5887087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin) 5897cb94ee6SHong Zhang { 5907cb94ee6SHong Zhang PetscFunctionBegin; 5915ef26d82SJed Brown if (nonlin) *nonlin = ts->snes_its; 5925ef26d82SJed Brown if (lin) *lin = ts->ksp_its; 5937cb94ee6SHong Zhang PetscFunctionReturn(0); 5947cb94ee6SHong Zhang } 5957cb94ee6SHong Zhang 5967087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s) 5972bfc04deSHong Zhang { 5982bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5992bfc04deSHong Zhang 6002bfc04deSHong Zhang PetscFunctionBegin; 6012bfc04deSHong Zhang cvode->monitorstep = s; 6022bfc04deSHong Zhang PetscFunctionReturn(0); 6032bfc04deSHong Zhang } 6047cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 6057cb94ee6SHong Zhang 6067cb94ee6SHong Zhang /*@C 6077cb94ee6SHong Zhang TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 6087cb94ee6SHong Zhang 6097cb94ee6SHong Zhang Not Collective 6107cb94ee6SHong Zhang 61197bb3fdcSJose E. Roman Input Parameter: 6127cb94ee6SHong Zhang . ts - the time-step context 6137cb94ee6SHong Zhang 6147cb94ee6SHong Zhang Output Parameters: 6157cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 6167cb94ee6SHong Zhang - lin - number of linear iterations 6177cb94ee6SHong Zhang 6187cb94ee6SHong Zhang Level: advanced 6197cb94ee6SHong Zhang 6207cb94ee6SHong Zhang Notes: 6217cb94ee6SHong Zhang These return the number since the creation of the TS object 6227cb94ee6SHong Zhang 623db781477SPatrick Sanan .seealso: `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 624db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 625db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 626db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsGetPC()`, `TSSetExactFinalTime()` 6277cb94ee6SHong Zhang 6287cb94ee6SHong Zhang @*/ 6297087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations(TS ts,int *nonlin,int *lin) 6307cb94ee6SHong Zhang { 6317cb94ee6SHong Zhang PetscFunctionBegin; 632cac4c232SBarry Smith PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin)); 6337cb94ee6SHong Zhang PetscFunctionReturn(0); 6347cb94ee6SHong Zhang } 6357cb94ee6SHong Zhang 6367cb94ee6SHong Zhang /*@ 6377cb94ee6SHong Zhang TSSundialsSetType - Sets the method that Sundials will use for integration. 6387cb94ee6SHong Zhang 6393f9fe445SBarry Smith Logically Collective on TS 6407cb94ee6SHong Zhang 6418f7d6fe5SPatrick Sanan Input Parameters: 6427cb94ee6SHong Zhang + ts - the time-step context 6437cb94ee6SHong Zhang - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 6447cb94ee6SHong Zhang 6457cb94ee6SHong Zhang Level: intermediate 6467cb94ee6SHong Zhang 647db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetMaxl()`, 648db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 649db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 650db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 651db781477SPatrick Sanan `TSSetExactFinalTime()` 6527cb94ee6SHong Zhang @*/ 6537087cfbeSBarry Smith PetscErrorCode TSSundialsSetType(TS ts,TSSundialsLmmType type) 6547cb94ee6SHong Zhang { 6557cb94ee6SHong Zhang PetscFunctionBegin; 656cac4c232SBarry Smith PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type)); 6577cb94ee6SHong Zhang PetscFunctionReturn(0); 6587cb94ee6SHong Zhang } 6597cb94ee6SHong Zhang 6607cb94ee6SHong Zhang /*@ 661c4e80e11SFlorian TSSundialsSetMaxord - Sets the maximum order for BDF/Adams method used by SUNDIALS. 662c4e80e11SFlorian 663c4e80e11SFlorian Logically Collective on TS 664c4e80e11SFlorian 6658f7d6fe5SPatrick Sanan Input Parameters: 666c4e80e11SFlorian + ts - the time-step context 667c4e80e11SFlorian - maxord - maximum order of BDF / Adams method 668c4e80e11SFlorian 669c4e80e11SFlorian Level: advanced 670c4e80e11SFlorian 671db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, 672db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 673db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 674db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 675db781477SPatrick Sanan `TSSetExactFinalTime()` 676c4e80e11SFlorian 677c4e80e11SFlorian @*/ 678c4e80e11SFlorian PetscErrorCode TSSundialsSetMaxord(TS ts,PetscInt maxord) 679c4e80e11SFlorian { 680c4e80e11SFlorian PetscFunctionBegin; 681c4e80e11SFlorian PetscValidLogicalCollectiveInt(ts,maxord,2); 682cac4c232SBarry Smith PetscTryMethod(ts,"TSSundialsSetMaxOrd_C",(TS,PetscInt),(ts,maxord)); 683c4e80e11SFlorian PetscFunctionReturn(0); 684c4e80e11SFlorian } 685c4e80e11SFlorian 686c4e80e11SFlorian /*@ 687f61b2b6aSHong Zhang TSSundialsSetMaxl - Sets the dimension of the Krylov space used by 6887cb94ee6SHong Zhang GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 689f61b2b6aSHong Zhang this is the maximum number of GMRES steps that will be used. 6907cb94ee6SHong Zhang 6913f9fe445SBarry Smith Logically Collective on TS 6927cb94ee6SHong Zhang 6938f7d6fe5SPatrick Sanan Input Parameters: 6947cb94ee6SHong Zhang + ts - the time-step context 695f61b2b6aSHong Zhang - maxl - number of direction vectors (the dimension of Krylov subspace). 6967cb94ee6SHong Zhang 6977cb94ee6SHong Zhang Level: advanced 6987cb94ee6SHong Zhang 699db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, 700db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 701db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 702db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 703db781477SPatrick Sanan `TSSetExactFinalTime()` 7047cb94ee6SHong Zhang 7057cb94ee6SHong Zhang @*/ 706f61b2b6aSHong Zhang PetscErrorCode TSSundialsSetMaxl(TS ts,PetscInt maxl) 7077cb94ee6SHong Zhang { 7087cb94ee6SHong Zhang PetscFunctionBegin; 709f61b2b6aSHong Zhang PetscValidLogicalCollectiveInt(ts,maxl,2); 710cac4c232SBarry Smith PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl)); 7117cb94ee6SHong Zhang PetscFunctionReturn(0); 7127cb94ee6SHong Zhang } 7137cb94ee6SHong Zhang 7147cb94ee6SHong Zhang /*@ 7157cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 7167cb94ee6SHong Zhang system by SUNDIALS. 7177cb94ee6SHong Zhang 7183f9fe445SBarry Smith Logically Collective on TS 7197cb94ee6SHong Zhang 7208f7d6fe5SPatrick Sanan Input Parameters: 7217cb94ee6SHong Zhang + ts - the time-step context 7227cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 7237cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 7247cb94ee6SHong Zhang 7257cb94ee6SHong Zhang Level: advanced 7267cb94ee6SHong Zhang 727db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 728db781477SPatrick Sanan `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 729db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 730db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 731db781477SPatrick Sanan `TSSetExactFinalTime()` 7327cb94ee6SHong Zhang 7337cb94ee6SHong Zhang @*/ 734590ae059SPatrick Sanan PetscErrorCode TSSundialsSetLinearTolerance(TS ts,PetscReal tol) 7357cb94ee6SHong Zhang { 7367cb94ee6SHong Zhang PetscFunctionBegin; 737c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,tol,2); 738cac4c232SBarry Smith PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,PetscReal),(ts,tol)); 7397cb94ee6SHong Zhang PetscFunctionReturn(0); 7407cb94ee6SHong Zhang } 7417cb94ee6SHong Zhang 7427cb94ee6SHong Zhang /*@ 7437cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 7447cb94ee6SHong Zhang in GMRES method by SUNDIALS linear solver. 7457cb94ee6SHong Zhang 7463f9fe445SBarry Smith Logically Collective on TS 7477cb94ee6SHong Zhang 7488f7d6fe5SPatrick Sanan Input Parameters: 7497cb94ee6SHong Zhang + ts - the time-step context 7507cb94ee6SHong Zhang - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 7517cb94ee6SHong Zhang 7527cb94ee6SHong Zhang Level: advanced 7537cb94ee6SHong Zhang 754db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 755db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, 756db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 757db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 758db781477SPatrick Sanan `TSSetExactFinalTime()` 7597cb94ee6SHong Zhang 7607cb94ee6SHong Zhang @*/ 7617087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type) 7627cb94ee6SHong Zhang { 7637cb94ee6SHong Zhang PetscFunctionBegin; 764cac4c232SBarry Smith PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type)); 7657cb94ee6SHong Zhang PetscFunctionReturn(0); 7667cb94ee6SHong Zhang } 7677cb94ee6SHong Zhang 7687cb94ee6SHong Zhang /*@ 7697cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 7707cb94ee6SHong Zhang Sundials for error control. 7717cb94ee6SHong Zhang 7723f9fe445SBarry Smith Logically Collective on TS 7737cb94ee6SHong Zhang 7748f7d6fe5SPatrick Sanan Input Parameters: 7757cb94ee6SHong Zhang + ts - the time-step context 7767cb94ee6SHong Zhang . aabs - the absolute tolerance 7777cb94ee6SHong Zhang - rel - the relative tolerance 7787cb94ee6SHong Zhang 7797cb94ee6SHong Zhang See the Cvode/Sundials users manual for exact details on these parameters. Essentially 7807cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 7817cb94ee6SHong Zhang 7827cb94ee6SHong Zhang Level: intermediate 7837cb94ee6SHong Zhang 784db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetGMRESMaxl()`, 785db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, 786db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 787db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 788db781477SPatrick Sanan `TSSetExactFinalTime()` 7897cb94ee6SHong Zhang 7907cb94ee6SHong Zhang @*/ 791590ae059SPatrick Sanan PetscErrorCode TSSundialsSetTolerance(TS ts,PetscReal aabs,PetscReal rel) 7927cb94ee6SHong Zhang { 7937cb94ee6SHong Zhang PetscFunctionBegin; 794cac4c232SBarry Smith PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,PetscReal,PetscReal),(ts,aabs,rel)); 7957cb94ee6SHong Zhang PetscFunctionReturn(0); 7967cb94ee6SHong Zhang } 7977cb94ee6SHong Zhang 7987cb94ee6SHong Zhang /*@ 7997cb94ee6SHong Zhang TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 8007cb94ee6SHong Zhang 8017cb94ee6SHong Zhang Input Parameter: 8027cb94ee6SHong Zhang . ts - the time-step context 8037cb94ee6SHong Zhang 8047cb94ee6SHong Zhang Output Parameter: 8057cb94ee6SHong Zhang . pc - the preconditioner context 8067cb94ee6SHong Zhang 8077cb94ee6SHong Zhang Level: advanced 8087cb94ee6SHong Zhang 809db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 810db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 811db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 812db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()` 8137cb94ee6SHong Zhang @*/ 8147087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC(TS ts,PC *pc) 8157cb94ee6SHong Zhang { 8167cb94ee6SHong Zhang PetscFunctionBegin; 817cac4c232SBarry Smith PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC*),(ts,pc)); 8187cb94ee6SHong Zhang PetscFunctionReturn(0); 8197cb94ee6SHong Zhang } 8207cb94ee6SHong Zhang 821f1cd61daSJed Brown /*@ 822f1cd61daSJed Brown TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 823f1cd61daSJed Brown 824d8d19677SJose E. Roman Input Parameters: 825f1cd61daSJed Brown + ts - the time-step context 826f1cd61daSJed Brown - mindt - lowest time step if positive, negative to deactivate 827f1cd61daSJed Brown 828fc6b9e64SJed Brown Note: 829fc6b9e64SJed Brown Sundials will error if it is not possible to keep the estimated truncation error below 830fc6b9e64SJed Brown the tolerance set with TSSundialsSetTolerance() without going below this step size. 831fc6b9e64SJed Brown 832f1cd61daSJed Brown Level: beginner 833f1cd61daSJed Brown 834db781477SPatrick Sanan .seealso: `TSSundialsSetType()`, `TSSundialsSetTolerance()`, 835f1cd61daSJed Brown @*/ 8367087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep(TS ts,PetscReal mindt) 837f1cd61daSJed Brown { 838f1cd61daSJed Brown PetscFunctionBegin; 839cac4c232SBarry Smith PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt)); 840f1cd61daSJed Brown PetscFunctionReturn(0); 841f1cd61daSJed Brown } 842f1cd61daSJed Brown 843f1cd61daSJed Brown /*@ 844f1cd61daSJed Brown TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 845f1cd61daSJed Brown 846d8d19677SJose E. Roman Input Parameters: 847f1cd61daSJed Brown + ts - the time-step context 848f1cd61daSJed Brown - maxdt - lowest time step if positive, negative to deactivate 849f1cd61daSJed Brown 850f1cd61daSJed Brown Level: beginner 851f1cd61daSJed Brown 852db781477SPatrick Sanan .seealso: `TSSundialsSetType()`, `TSSundialsSetTolerance()`, 853f1cd61daSJed Brown @*/ 8547087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt) 855f1cd61daSJed Brown { 856f1cd61daSJed Brown PetscFunctionBegin; 857cac4c232SBarry Smith PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt)); 858f1cd61daSJed Brown PetscFunctionReturn(0); 859f1cd61daSJed Brown } 860f1cd61daSJed Brown 8612bfc04deSHong Zhang /*@ 8622bfc04deSHong Zhang TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 8632bfc04deSHong Zhang 864d8d19677SJose E. Roman Input Parameters: 8652bfc04deSHong Zhang + ts - the time-step context 8662bfc04deSHong Zhang - ft - PETSC_TRUE if monitor, else PETSC_FALSE 8672bfc04deSHong Zhang 8682bfc04deSHong Zhang Level: beginner 8692bfc04deSHong Zhang 870f61b2b6aSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 8712bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 872f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 8732bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 8742bfc04deSHong Zhang @*/ 8757087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps(TS ts,PetscBool ft) 8762bfc04deSHong Zhang { 8772bfc04deSHong Zhang PetscFunctionBegin; 878cac4c232SBarry Smith PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft)); 8792bfc04deSHong Zhang PetscFunctionReturn(0); 8802bfc04deSHong Zhang } 881918687b8SPatrick Sanan 882918687b8SPatrick Sanan /*@ 883918687b8SPatrick Sanan TSSundialsSetUseDense - Set a flag to use a dense linear solver in SUNDIALS (serial only) 884918687b8SPatrick Sanan 885918687b8SPatrick Sanan Logically Collective 886918687b8SPatrick Sanan 887918687b8SPatrick Sanan Input Parameters: 888918687b8SPatrick Sanan + ts - the time-step context 889918687b8SPatrick Sanan - use_dense - PETSC_TRUE to use the dense solver 890918687b8SPatrick Sanan 891918687b8SPatrick Sanan Level: advanced 892918687b8SPatrick Sanan 893db781477SPatrick Sanan .seealso: `TSSUNDIALS` 894918687b8SPatrick Sanan 895918687b8SPatrick Sanan @*/ 896918687b8SPatrick Sanan PetscErrorCode TSSundialsSetUseDense(TS ts,PetscBool use_dense) 897918687b8SPatrick Sanan { 898918687b8SPatrick Sanan PetscFunctionBegin; 899918687b8SPatrick Sanan PetscValidLogicalCollectiveInt(ts,use_dense,2); 900cac4c232SBarry Smith PetscTryMethod(ts,"TSSundialsSetUseDense_C",(TS,PetscBool),(ts,use_dense)); 901918687b8SPatrick Sanan PetscFunctionReturn(0); 902918687b8SPatrick Sanan } 903918687b8SPatrick Sanan 9047cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 9057cb94ee6SHong Zhang /*MC 90696f5712cSJed Brown TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 9077cb94ee6SHong Zhang 9087cb94ee6SHong Zhang Options Database: 909897c9f78SHong Zhang + -ts_sundials_type <bdf,adams> - 9107cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 9117cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 9127cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 913897c9f78SHong Zhang . -ts_sundials_linear_tolerance <tol> - 914f61b2b6aSHong Zhang . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace 915918687b8SPatrick Sanan . -ts_sundials_monitor_steps - Monitor SUNDIALS internal steps 916918687b8SPatrick Sanan - -ts_sundials_use_dense - Use a dense linear solver within CVODE (serial only) 91716016685SHong Zhang 91895452b02SPatrick Sanan Notes: 91995452b02SPatrick Sanan This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply, 920897c9f78SHong Zhang only PETSc PC options. 9217cb94ee6SHong Zhang 9227cb94ee6SHong Zhang Level: beginner 9237cb94ee6SHong Zhang 924db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, `TSSundialsSetLinearTolerance()`, 925db781477SPatrick Sanan `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, `TSSundialsGetIterations()`, `TSSetExactFinalTime()` 9267cb94ee6SHong Zhang 9277cb94ee6SHong Zhang M*/ 9288cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts) 9297cb94ee6SHong Zhang { 9307cb94ee6SHong Zhang TS_Sundials *cvode; 93142528757SHong Zhang PC pc; 9327cb94ee6SHong Zhang 9337cb94ee6SHong Zhang PetscFunctionBegin; 934277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Sundials; 93528adb3f7SHong Zhang ts->ops->destroy = TSDestroy_Sundials; 93628adb3f7SHong Zhang ts->ops->view = TSView_Sundials; 937214bc6a2SJed Brown ts->ops->setup = TSSetUp_Sundials; 938b4eba00bSSean Farley ts->ops->step = TSStep_Sundials; 939b4eba00bSSean Farley ts->ops->interpolate = TSInterpolate_Sundials; 940214bc6a2SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Sundials; 9412ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 9427cb94ee6SHong Zhang 9439566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts,&cvode)); 944bbd56ea5SKarl Rupp 945efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 946efd4aadfSBarry Smith 9477cb94ee6SHong Zhang ts->data = (void*)cvode; 9486fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF; 9496fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS; 950f61b2b6aSHong Zhang cvode->maxl = 5; 951c4e80e11SFlorian cvode->maxord = PETSC_DEFAULT; 9527cb94ee6SHong Zhang cvode->linear_tol = .05; 953b4eba00bSSean Farley cvode->monitorstep = PETSC_TRUE; 954918687b8SPatrick Sanan cvode->use_dense = PETSC_FALSE; 9557cb94ee6SHong Zhang 9569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)ts),&(cvode->comm_sundials))); 957f1cd61daSJed Brown 958f1cd61daSJed Brown cvode->mindt = -1.; 959f1cd61daSJed Brown cvode->maxdt = -1.; 960f1cd61daSJed Brown 9617cb94ee6SHong Zhang /* set tolerance for Sundials */ 9627cb94ee6SHong Zhang cvode->reltol = 1e-6; 9632c823083SHong Zhang cvode->abstol = 1e-6; 9647cb94ee6SHong Zhang 96542528757SHong Zhang /* set PCNONE as default pctype */ 9669566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC_Sundials(ts,&pc)); 9679566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCNONE)); 96842528757SHong Zhang 9699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",TSSundialsSetType_Sundials)); 9709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",TSSundialsSetMaxl_Sundials)); 9719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",TSSundialsSetLinearTolerance_Sundials)); 9729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",TSSundialsSetGramSchmidtType_Sundials)); 9739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",TSSundialsSetTolerance_Sundials)); 9749566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",TSSundialsSetMinTimeStep_Sundials)); 9759566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",TSSundialsSetMaxTimeStep_Sundials)); 9769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",TSSundialsGetPC_Sundials)); 9779566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",TSSundialsGetIterations_Sundials)); 9789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",TSSundialsMonitorInternalSteps_Sundials)); 9799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetUseDense_C",TSSundialsSetUseDense_Sundials)); 9807cb94ee6SHong Zhang PetscFunctionReturn(0); 9817cb94ee6SHong Zhang } 982