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 */ 14*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPrecond_Sundials(realtype tn, N_Vector y, N_Vector fy, booleantype jok, booleantype *jcurPtr, realtype _gamma, void *P_data, N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3) 15*d71ae5a4SJacob Faibussowitsch { 167cb94ee6SHong Zhang TS ts = (TS)P_data; 177cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 18dcbc6d53SSean Farley PC pc; 190679a0aeSJed Brown Mat J, P; 200679a0aeSJed Brown Vec yy = cvode->w1, yydot = cvode->ydot; 210679a0aeSJed Brown PetscReal gm = (PetscReal)_gamma; 227cb94ee6SHong Zhang PetscScalar *y_data; 237cb94ee6SHong Zhang 247cb94ee6SHong Zhang PetscFunctionBegin; 259566063dSJacob Faibussowitsch PetscCall(TSGetIJacobian(ts, &J, &P, NULL, NULL)); 267cb94ee6SHong Zhang y_data = (PetscScalar *)N_VGetArrayPointer(y); 279566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(yy, y_data)); 289566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(yydot)); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */ 290679a0aeSJed Brown /* compute the shifted Jacobian (1/gm)*I + Jrest */ 309566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, ts->ptime, yy, yydot, 1 / gm, J, P, PETSC_FALSE)); 319566063dSJacob Faibussowitsch PetscCall(VecResetArray(yy)); 329566063dSJacob Faibussowitsch PetscCall(MatScale(P, gm)); /* turn into I-gm*Jrest, J is not used by Sundials */ 337cb94ee6SHong Zhang *jcurPtr = TRUE; 349566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc)); 359566063dSJacob Faibussowitsch PetscCall(PCSetOperators(pc, J, P)); 367cb94ee6SHong Zhang PetscFunctionReturn(0); 377cb94ee6SHong Zhang } 387cb94ee6SHong Zhang 397cb94ee6SHong Zhang /* 407cb94ee6SHong Zhang TSPSolve_Sundials - routine that we provide to Sundials that applies the preconditioner. 417cb94ee6SHong Zhang */ 42*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPSolve_Sundials(realtype tn, N_Vector y, N_Vector fy, N_Vector r, N_Vector z, realtype _gamma, realtype delta, int lr, void *P_data, N_Vector vtemp) 43*d71ae5a4SJacob Faibussowitsch { 447cb94ee6SHong Zhang TS ts = (TS)P_data; 457cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 46dcbc6d53SSean Farley PC pc; 477cb94ee6SHong Zhang Vec rr = cvode->w1, zz = cvode->w2; 487cb94ee6SHong Zhang PetscScalar *r_data, *z_data; 497cb94ee6SHong Zhang 507cb94ee6SHong Zhang PetscFunctionBegin; 517cb94ee6SHong Zhang /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/ 527cb94ee6SHong Zhang r_data = (PetscScalar *)N_VGetArrayPointer(r); 537cb94ee6SHong Zhang z_data = (PetscScalar *)N_VGetArrayPointer(z); 549566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(rr, r_data)); 559566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(zz, z_data)); 564ac7836bSHong Zhang 577cb94ee6SHong Zhang /* Solve the Px=r and put the result in zz */ 589566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc)); 599566063dSJacob Faibussowitsch PetscCall(PCApply(pc, rr, zz)); 609566063dSJacob Faibussowitsch PetscCall(VecResetArray(rr)); 619566063dSJacob Faibussowitsch PetscCall(VecResetArray(zz)); 627cb94ee6SHong Zhang PetscFunctionReturn(0); 637cb94ee6SHong Zhang } 647cb94ee6SHong Zhang 657cb94ee6SHong Zhang /* 667cb94ee6SHong Zhang TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side. 677cb94ee6SHong Zhang */ 68*d71ae5a4SJacob Faibussowitsch int TSFunction_Sundials(realtype t, N_Vector y, N_Vector ydot, void *ctx) 69*d71ae5a4SJacob Faibussowitsch { 707cb94ee6SHong Zhang TS ts = (TS)ctx; 715fcef5e4SJed Brown DM dm; 72942e3340SBarry Smith DMTS tsdm; 735fcef5e4SJed Brown TSIFunction ifunction; 74ce94432eSBarry Smith MPI_Comm comm; 757cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 760679a0aeSJed Brown Vec yy = cvode->w1, yyd = cvode->w2, yydot = cvode->ydot; 777cb94ee6SHong Zhang PetscScalar *y_data, *ydot_data; 787cb94ee6SHong Zhang 797cb94ee6SHong Zhang PetscFunctionBegin; 809566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 815ff03cf7SDinesh Kaushik /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/ 827cb94ee6SHong Zhang y_data = (PetscScalar *)N_VGetArrayPointer(y); 837cb94ee6SHong Zhang ydot_data = (PetscScalar *)N_VGetArrayPointer(ydot); 849566063dSJacob Faibussowitsch PetscCallAbort(comm, VecPlaceArray(yy, y_data)); 859566063dSJacob Faibussowitsch PetscCallAbort(comm, VecPlaceArray(yyd, ydot_data)); 864ac7836bSHong Zhang 875fcef5e4SJed Brown /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */ 889566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 899566063dSJacob Faibussowitsch PetscCall(DMGetDMTS(dm, &tsdm)); 909566063dSJacob Faibussowitsch PetscCall(DMTSGetIFunction(dm, &ifunction, NULL)); 915fcef5e4SJed Brown if (!ifunction) { 929566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, t, yy, yyd)); 93bc0cc02bSJed Brown } else { /* If rhsfunction is also set, this computes both parts and shifts them to the right */ 949566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(yydot)); 959566063dSJacob Faibussowitsch PetscCallAbort(comm, TSComputeIFunction(ts, t, yy, yydot, yyd, PETSC_FALSE)); 969566063dSJacob Faibussowitsch PetscCall(VecScale(yyd, -1.)); 97bc0cc02bSJed Brown } 989566063dSJacob Faibussowitsch PetscCallAbort(comm, VecResetArray(yy)); 999566063dSJacob Faibussowitsch PetscCallAbort(comm, VecResetArray(yyd)); 1004ac7836bSHong Zhang PetscFunctionReturn(0); 1017cb94ee6SHong Zhang } 1027cb94ee6SHong Zhang 1037cb94ee6SHong Zhang /* 104b4eba00bSSean Farley TSStep_Sundials - Calls Sundials to integrate the ODE. 1057cb94ee6SHong Zhang */ 106*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSStep_Sundials(TS ts) 107*d71ae5a4SJacob Faibussowitsch { 1087cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 109b4eba00bSSean Farley PetscInt flag; 110be5899b3SLisandro Dalcin long int nits, lits, nsteps; 1117cb94ee6SHong Zhang realtype t, tout; 1127cb94ee6SHong Zhang PetscScalar *y_data; 1137cb94ee6SHong Zhang void *mem; 1147cb94ee6SHong Zhang 1157cb94ee6SHong Zhang PetscFunctionBegin; 11616016685SHong Zhang mem = cvode->mem; 1177cb94ee6SHong Zhang tout = ts->max_time; 1189566063dSJacob Faibussowitsch PetscCall(VecGetArray(ts->vec_sol, &y_data)); 1197cb94ee6SHong Zhang N_VSetArrayPointer((realtype *)y_data, cvode->y); 1209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ts->vec_sol, NULL)); 121186e87acSLisandro Dalcin 1229687d888SLisandro Dalcin /* We would like to TSPreStage() and TSPostStage() 1239687d888SLisandro Dalcin * before each stage solve but CVode does not appear to support this. */ 1249371c9d4SSatish Balay if (cvode->monitorstep) flag = CVode(mem, tout, cvode->y, &t, CV_ONE_STEP); 1259371c9d4SSatish Balay else flag = CVode(mem, tout, cvode->y, &t, CV_NORMAL); 1269f94935aSHong Zhang 1279f94935aSHong Zhang if (flag) { /* display error message */ 1289f94935aSHong Zhang switch (flag) { 129*d71ae5a4SJacob Faibussowitsch case CV_ILL_INPUT: 130*d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_ILL_INPUT"); 131*d71ae5a4SJacob Faibussowitsch break; 132*d71ae5a4SJacob Faibussowitsch case CV_TOO_CLOSE: 133*d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_TOO_CLOSE"); 134*d71ae5a4SJacob Faibussowitsch break; 1353c7fefeeSJed Brown case CV_TOO_MUCH_WORK: { 1369f94935aSHong Zhang PetscReal tcur; 1379566063dSJacob Faibussowitsch PetscCall(CVodeGetNumSteps(mem, &nsteps)); 1389566063dSJacob Faibussowitsch PetscCall(CVodeGetCurrentTime(mem, &tcur)); 13963a3b9bcSJacob 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); 1403c7fefeeSJed Brown } break; 141*d71ae5a4SJacob Faibussowitsch case CV_TOO_MUCH_ACC: 142*d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_TOO_MUCH_ACC"); 143*d71ae5a4SJacob Faibussowitsch break; 144*d71ae5a4SJacob Faibussowitsch case CV_ERR_FAILURE: 145*d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_ERR_FAILURE"); 146*d71ae5a4SJacob Faibussowitsch break; 147*d71ae5a4SJacob Faibussowitsch case CV_CONV_FAILURE: 148*d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_CONV_FAILURE"); 149*d71ae5a4SJacob Faibussowitsch break; 150*d71ae5a4SJacob Faibussowitsch case CV_LINIT_FAIL: 151*d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LINIT_FAIL"); 152*d71ae5a4SJacob Faibussowitsch break; 153*d71ae5a4SJacob Faibussowitsch case CV_LSETUP_FAIL: 154*d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LSETUP_FAIL"); 155*d71ae5a4SJacob Faibussowitsch break; 156*d71ae5a4SJacob Faibussowitsch case CV_LSOLVE_FAIL: 157*d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LSOLVE_FAIL"); 158*d71ae5a4SJacob Faibussowitsch break; 159*d71ae5a4SJacob Faibussowitsch case CV_RHSFUNC_FAIL: 160*d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_RHSFUNC_FAIL"); 161*d71ae5a4SJacob Faibussowitsch break; 162*d71ae5a4SJacob Faibussowitsch case CV_FIRST_RHSFUNC_ERR: 163*d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_FIRST_RHSFUNC_ERR"); 164*d71ae5a4SJacob Faibussowitsch break; 165*d71ae5a4SJacob Faibussowitsch case CV_REPTD_RHSFUNC_ERR: 166*d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_REPTD_RHSFUNC_ERR"); 167*d71ae5a4SJacob Faibussowitsch break; 168*d71ae5a4SJacob Faibussowitsch case CV_UNREC_RHSFUNC_ERR: 169*d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_UNREC_RHSFUNC_ERR"); 170*d71ae5a4SJacob Faibussowitsch break; 171*d71ae5a4SJacob Faibussowitsch case CV_RTFUNC_FAIL: 172*d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_RTFUNC_FAIL"); 173*d71ae5a4SJacob Faibussowitsch break; 174*d71ae5a4SJacob Faibussowitsch default: 175*d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, flag %d", flag); 1769f94935aSHong Zhang } 1779f94935aSHong Zhang } 1789f94935aSHong Zhang 179be5899b3SLisandro Dalcin /* log inner nonlinear and linear iterations */ 1809566063dSJacob Faibussowitsch PetscCall(CVodeGetNumNonlinSolvIters(mem, &nits)); 1819566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumLinIters(mem, &lits)); 1829371c9d4SSatish Balay ts->snes_its += nits; 1839371c9d4SSatish Balay ts->ksp_its = lits; 184be5899b3SLisandro Dalcin 1857cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */ 1869566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(cvode->w1, y_data)); 1879566063dSJacob Faibussowitsch PetscCall(VecCopy(cvode->w1, cvode->update)); 1889566063dSJacob Faibussowitsch PetscCall(VecResetArray(cvode->w1)); 1899566063dSJacob Faibussowitsch PetscCall(VecCopy(cvode->update, ts->vec_sol)); 190186e87acSLisandro Dalcin 191186e87acSLisandro Dalcin ts->time_step = t - ts->ptime; 192186e87acSLisandro Dalcin ts->ptime = t; 193186e87acSLisandro Dalcin 1949566063dSJacob Faibussowitsch PetscCall(CVodeGetNumSteps(mem, &nsteps)); 1952808aa04SLisandro Dalcin if (!cvode->monitorstep) ts->steps += nsteps - 1; /* TSStep() increments the step counter by one */ 196b4eba00bSSean Farley PetscFunctionReturn(0); 197b4eba00bSSean Farley } 198b4eba00bSSean Farley 199*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_Sundials(TS ts, PetscReal t, Vec X) 200*d71ae5a4SJacob Faibussowitsch { 201b4eba00bSSean Farley TS_Sundials *cvode = (TS_Sundials *)ts->data; 202b4eba00bSSean Farley N_Vector y; 203b4eba00bSSean Farley PetscScalar *x_data; 204b4eba00bSSean Farley PetscInt glosize, locsize; 205b4eba00bSSean Farley 206b4eba00bSSean Farley PetscFunctionBegin; 207b4eba00bSSean Farley /* get the vector size */ 2089566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &glosize)); 2099566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &locsize)); 2109566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x_data)); 211b4eba00bSSean Farley 2124d78c9e0SClaas Abert /* Initialize N_Vec y with x_data */ 213918687b8SPatrick Sanan if (cvode->use_dense) { 214918687b8SPatrick Sanan PetscMPIInt size; 215918687b8SPatrick Sanan 2169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 2173c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "TSSUNDIALS only supports a dense solve in the serial case"); 218918687b8SPatrick Sanan y = N_VMake_Serial(locsize, (realtype *)x_data); 219918687b8SPatrick Sanan } else { 2204d78c9e0SClaas Abert y = N_VMake_Parallel(cvode->comm_sundials, locsize, glosize, (realtype *)x_data); 221918687b8SPatrick Sanan } 222918687b8SPatrick Sanan 2233c633725SBarry Smith PetscCheck(y, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Interpolated y is not allocated"); 224b4eba00bSSean Farley 2259566063dSJacob Faibussowitsch PetscCall(CVodeGetDky(cvode->mem, t, 0, y)); 2269566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x_data)); 2277cb94ee6SHong Zhang PetscFunctionReturn(0); 2287cb94ee6SHong Zhang } 2297cb94ee6SHong Zhang 230*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSReset_Sundials(TS ts) 231*d71ae5a4SJacob Faibussowitsch { 232277b19d0SLisandro Dalcin TS_Sundials *cvode = (TS_Sundials *)ts->data; 233277b19d0SLisandro Dalcin 234277b19d0SLisandro Dalcin PetscFunctionBegin; 2359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->update)); 2369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->ydot)); 2379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->w1)); 2389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->w2)); 239bbd56ea5SKarl Rupp if (cvode->mem) CVodeFree(&cvode->mem); 240277b19d0SLisandro Dalcin PetscFunctionReturn(0); 241277b19d0SLisandro Dalcin } 242277b19d0SLisandro Dalcin 243*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDestroy_Sundials(TS ts) 244*d71ae5a4SJacob Faibussowitsch { 2457cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 2467cb94ee6SHong Zhang 2477cb94ee6SHong Zhang PetscFunctionBegin; 2489566063dSJacob Faibussowitsch PetscCall(TSReset_Sundials(ts)); 2499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&(cvode->comm_sundials))); 2509566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 2519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetType_C", NULL)); 2529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxl_C", NULL)); 2539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetLinearTolerance_C", NULL)); 2549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetGramSchmidtType_C", NULL)); 2559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetTolerance_C", NULL)); 2569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMinTimeStep_C", NULL)); 2579566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxTimeStep_C", NULL)); 2589566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetPC_C", NULL)); 2599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetIterations_C", NULL)); 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsMonitorInternalSteps_C", NULL)); 2612e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetUseDense_C", NULL)); 2627cb94ee6SHong Zhang PetscFunctionReturn(0); 2637cb94ee6SHong Zhang } 2647cb94ee6SHong Zhang 265*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetUp_Sundials(TS ts) 266*d71ae5a4SJacob Faibussowitsch { 2677cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 26816016685SHong Zhang PetscInt glosize, locsize, i, flag; 2697cb94ee6SHong Zhang PetscScalar *y_data, *parray; 27016016685SHong Zhang void *mem; 271dcbc6d53SSean Farley PC pc; 27219fd82e9SBarry Smith PCType pctype; 273ace3abfcSBarry Smith PetscBool pcnone; 2747cb94ee6SHong Zhang 2757cb94ee6SHong Zhang PetscFunctionBegin; 27608401ef6SPierre 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"); 277236c45dcSShri Abhyankar 2787cb94ee6SHong Zhang /* get the vector size */ 2799566063dSJacob Faibussowitsch PetscCall(VecGetSize(ts->vec_sol, &glosize)); 2809566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ts->vec_sol, &locsize)); 2817cb94ee6SHong Zhang 2827cb94ee6SHong Zhang /* allocate the memory for N_Vec y */ 283918687b8SPatrick Sanan if (cvode->use_dense) { 284918687b8SPatrick Sanan PetscMPIInt size; 285918687b8SPatrick Sanan 2869566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 2873c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "TSSUNDIALS only supports a dense solve in the serial case"); 288918687b8SPatrick Sanan cvode->y = N_VNew_Serial(locsize); 289918687b8SPatrick Sanan } else { 290a07356b8SSatish Balay cvode->y = N_VNew_Parallel(cvode->comm_sundials, locsize, glosize); 291918687b8SPatrick Sanan } 2923c633725SBarry Smith PetscCheck(cvode->y, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "cvode->y is not allocated"); 2937cb94ee6SHong Zhang 29428adb3f7SHong Zhang /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 2959566063dSJacob Faibussowitsch PetscCall(VecGetArray(ts->vec_sol, &parray)); 2967cb94ee6SHong Zhang y_data = (PetscScalar *)N_VGetArrayPointer(cvode->y); 2977cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 2989566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ts->vec_sol, NULL)); 29942528757SHong Zhang 3009566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &cvode->update)); 3019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &cvode->ydot)); 3027cb94ee6SHong Zhang 3037cb94ee6SHong Zhang /* 3047cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are 3057cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided 3067cb94ee6SHong Zhang by Sundials and set using VecPlaceArray(). 3077cb94ee6SHong Zhang */ 3089566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts), 1, locsize, PETSC_DECIDE, NULL, &cvode->w1)); 3099566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts), 1, locsize, PETSC_DECIDE, NULL, &cvode->w2)); 31016016685SHong Zhang 31116016685SHong Zhang /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 31216016685SHong Zhang mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 3133c633725SBarry Smith PetscCheck(mem, PETSC_COMM_SELF, PETSC_ERR_MEM, "CVodeCreate() fails"); 31416016685SHong Zhang cvode->mem = mem; 31516016685SHong Zhang 31616016685SHong Zhang /* Set the pointer to user-defined data */ 31716016685SHong Zhang flag = CVodeSetUserData(mem, ts); 31828b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeSetUserData() fails"); 31916016685SHong Zhang 320fc6b9e64SJed Brown /* Sundials may choose to use a smaller initial step, but will never use a larger step. */ 321cdbf8f93SLisandro Dalcin flag = CVodeSetInitStep(mem, (realtype)ts->time_step); 32228b400f6SJacob Faibussowitsch PetscCheck(!flag, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetInitStep() failed"); 323f1cd61daSJed Brown if (cvode->mindt > 0) { 324f1cd61daSJed Brown flag = CVodeSetMinStep(mem, (realtype)cvode->mindt); 3259f94935aSHong Zhang if (flag) { 32608401ef6SPierre Jolivet PetscCheck(flag != CV_MEM_NULL, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMinStep() failed, cvode_mem pointer is NULL"); 327f7d195e4SLawrence Mitchell PetscCheck(flag != CV_ILL_INPUT, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size"); 328f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMinStep() failed"); 3299f94935aSHong Zhang } 330f1cd61daSJed Brown } 331f1cd61daSJed Brown if (cvode->maxdt > 0) { 332f1cd61daSJed Brown flag = CVodeSetMaxStep(mem, (realtype)cvode->maxdt); 33328b400f6SJacob Faibussowitsch PetscCheck(!flag, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMaxStep() failed"); 334f1cd61daSJed Brown } 335f1cd61daSJed Brown 33616016685SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the 337a5b23f4aSJose E. Roman * user's right hand side function in u'=f(t,u), the initial time T0, and 33816016685SHong Zhang * the initial dependent variable vector cvode->y */ 33916016685SHong Zhang flag = CVodeInit(mem, TSFunction_Sundials, ts->ptime, cvode->y); 34028b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeInit() fails, flag %d", flag); 34116016685SHong Zhang 3429f94935aSHong Zhang /* specifies scalar relative and absolute tolerances */ 34316016685SHong Zhang flag = CVodeSStolerances(mem, cvode->reltol, cvode->abstol); 34428b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeSStolerances() fails, flag %d", flag); 34516016685SHong Zhang 346c4e80e11SFlorian /* Specify max order of BDF / ADAMS method */ 347c4e80e11SFlorian if (cvode->maxord != PETSC_DEFAULT) { 348c4e80e11SFlorian flag = CVodeSetMaxOrd(mem, cvode->maxord); 34928b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeSetMaxOrd() fails, flag %d", flag); 350c4e80e11SFlorian } 351c4e80e11SFlorian 3529f94935aSHong Zhang /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */ 3539f94935aSHong Zhang flag = CVodeSetMaxNumSteps(mem, ts->max_steps); 35428b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeSetMaxNumSteps() fails, flag %d", flag); 35516016685SHong Zhang 356918687b8SPatrick Sanan if (cvode->use_dense) { 357918687b8SPatrick Sanan /* call CVDense to use a dense linear solver. */ 358918687b8SPatrick Sanan PetscMPIInt size; 359918687b8SPatrick Sanan 3609566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 3613c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "TSSUNDIALS only supports a dense solve in the serial case"); 362918687b8SPatrick Sanan flag = CVDense(mem, locsize); 36328b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVDense() fails, flag %d", flag); 364918687b8SPatrick Sanan } else { 36516016685SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 36616016685SHong Zhang /* setup the ode integrator with the given preconditioner */ 3679566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc)); 3689566063dSJacob Faibussowitsch PetscCall(PCGetType(pc, &pctype)); 3699566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCNONE, &pcnone)); 37016016685SHong Zhang if (pcnone) { 37116016685SHong Zhang flag = CVSpgmr(mem, PREC_NONE, 0); 37228b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVSpgmr() fails, flag %d", flag); 37316016685SHong Zhang } else { 374f61b2b6aSHong Zhang flag = CVSpgmr(mem, PREC_LEFT, cvode->maxl); 37528b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVSpgmr() fails, flag %d", flag); 37616016685SHong Zhang 37716016685SHong Zhang /* Set preconditioner and solve routines Precond and PSolve, 37816016685SHong Zhang and the pointer to the user-defined block data */ 37916016685SHong Zhang flag = CVSpilsSetPreconditioner(mem, TSPrecond_Sundials, TSPSolve_Sundials); 38028b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVSpilsSetPreconditioner() fails, flag %d", flag); 38116016685SHong Zhang } 382918687b8SPatrick Sanan } 3837cb94ee6SHong Zhang PetscFunctionReturn(0); 3847cb94ee6SHong Zhang } 3857cb94ee6SHong Zhang 3866fadb2cdSHong Zhang /* type of CVODE linear multistep method */ 387c793f718SLisandro Dalcin const char *const TSSundialsLmmTypes[] = {"", "ADAMS", "BDF", "TSSundialsLmmType", "SUNDIALS_", NULL}; 3886fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */ 389c793f718SLisandro Dalcin const char *const TSSundialsGramSchmidtTypes[] = {"", "MODIFIED", "CLASSICAL", "TSSundialsGramSchmidtType", "SUNDIALS_", NULL}; 390a04cf4d8SBarry Smith 391*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetFromOptions_Sundials(TS ts, PetscOptionItems *PetscOptionsObject) 392*d71ae5a4SJacob Faibussowitsch { 3937cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 3947cb94ee6SHong Zhang int indx; 395ace3abfcSBarry Smith PetscBool flag; 3967bda3a07SJed Brown PC pc; 3977cb94ee6SHong Zhang 3987cb94ee6SHong Zhang PetscFunctionBegin; 399d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SUNDIALS ODE solver options"); 4009566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_sundials_type", "Scheme", "TSSundialsSetType", TSSundialsLmmTypes, 3, TSSundialsLmmTypes[cvode->cvode_type], &indx, &flag)); 4011baa6e33SBarry Smith if (flag) PetscCall(TSSundialsSetType(ts, (TSSundialsLmmType)indx)); 4029566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_sundials_gramschmidt_type", "Type of orthogonalization", "TSSundialsSetGramSchmidtType", TSSundialsGramSchmidtTypes, 3, TSSundialsGramSchmidtTypes[cvode->gtype], &indx, &flag)); 4031baa6e33SBarry Smith if (flag) PetscCall(TSSundialsSetGramSchmidtType(ts, (TSSundialsGramSchmidtType)indx)); 4049566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_atol", "Absolute tolerance for convergence", "TSSundialsSetTolerance", cvode->abstol, &cvode->abstol, NULL)); 4059566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_rtol", "Relative tolerance for convergence", "TSSundialsSetTolerance", cvode->reltol, &cvode->reltol, NULL)); 4069566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_mindt", "Minimum step size", "TSSundialsSetMinTimeStep", cvode->mindt, &cvode->mindt, NULL)); 4079566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_maxdt", "Maximum step size", "TSSundialsSetMaxTimeStep", cvode->maxdt, &cvode->maxdt, NULL)); 4089566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_linear_tolerance", "Convergence tolerance for linear solve", "TSSundialsSetLinearTolerance", cvode->linear_tol, &cvode->linear_tol, NULL)); 4099566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_sundials_maxord", "Max Order for BDF/Adams method", "TSSundialsSetMaxOrd", cvode->maxord, &cvode->maxord, NULL)); 4109566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_sundials_maxl", "Max dimension of the Krylov subspace", "TSSundialsSetMaxl", cvode->maxl, &cvode->maxl, NULL)); 4119566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_sundials_monitor_steps", "Monitor SUNDIALS internal steps", "TSSundialsMonitorInternalSteps", cvode->monitorstep, &cvode->monitorstep, NULL)); 4129566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_sundials_use_dense", "Use dense internal solver in SUNDIALS (serial only)", "TSSundialsSetUseDense", cvode->use_dense, &cvode->use_dense, NULL)); 413d0609cedSBarry Smith PetscOptionsHeadEnd(); 4149566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc)); 4159566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions(pc)); 4167cb94ee6SHong Zhang PetscFunctionReturn(0); 4177cb94ee6SHong Zhang } 4187cb94ee6SHong Zhang 419*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSView_Sundials(TS ts, PetscViewer viewer) 420*d71ae5a4SJacob Faibussowitsch { 4217cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 4227cb94ee6SHong Zhang char *type; 4237cb94ee6SHong Zhang char atype[] = "Adams"; 4247cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 425ace3abfcSBarry Smith PetscBool iascii, isstring; 4262c823083SHong Zhang long int nsteps, its, nfevals, nlinsetups, nfails, itmp; 4272c823083SHong Zhang PetscInt qlast, qcur; 4285d47aee6SHong Zhang PetscReal hinused, hlast, hcur, tcur, tolsfac; 42942528757SHong Zhang PC pc; 4307cb94ee6SHong Zhang 4317cb94ee6SHong Zhang PetscFunctionBegin; 432bbd56ea5SKarl Rupp if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype; 433bbd56ea5SKarl Rupp else type = btype; 4347cb94ee6SHong Zhang 4359566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 4369566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 4377cb94ee6SHong Zhang if (iascii) { 4389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials integrator does not use SNES!\n")); 4399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials integrator type %s\n", type)); 44063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials maxord %" PetscInt_FMT "\n", cvode->maxord)); 4419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials abs tol %g rel tol %g\n", (double)cvode->abstol, (double)cvode->reltol)); 442918687b8SPatrick Sanan if (cvode->use_dense) { 4439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials integrator using a dense linear solve\n")); 444918687b8SPatrick Sanan } else { 4459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials linear solver tolerance factor %g\n", (double)cvode->linear_tol)); 44663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials max dimension of Krylov subspace %" PetscInt_FMT "\n", cvode->maxl)); 4477cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 4489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n")); 4497cb94ee6SHong Zhang } else { 4509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n")); 4517cb94ee6SHong Zhang } 452918687b8SPatrick Sanan } 4539566063dSJacob Faibussowitsch if (cvode->mindt > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials minimum time step %g\n", (double)cvode->mindt)); 4549566063dSJacob Faibussowitsch if (cvode->maxdt > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials maximum time step %g\n", (double)cvode->maxdt)); 4552c823083SHong Zhang 4565d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */ 4579566063dSJacob Faibussowitsch PetscCall(CVodeGetTolScaleFactor(cvode->mem, &tolsfac)); 4589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials suggested factor for tolerance scaling %g\n", tolsfac)); 4599566063dSJacob Faibussowitsch PetscCall(CVodeGetIntegratorStats(cvode->mem, &nsteps, &nfevals, &nlinsetups, &nfails, &qlast, &qcur, &hinused, &hlast, &hcur, &tcur)); 46063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials cumulative number of internal steps %ld\n", nsteps)); 46163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of calls to rhs function %ld\n", nfevals)); 46263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of calls to linear solver setup function %ld\n", nlinsetups)); 46363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of error test failures %ld\n", nfails)); 4642c823083SHong Zhang 4659566063dSJacob Faibussowitsch PetscCall(CVodeGetNonlinSolvStats(cvode->mem, &its, &nfails)); 46663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of nonlinear solver iterations %ld\n", its)); 46763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of nonlinear convergence failure %ld\n", nfails)); 468918687b8SPatrick Sanan if (!cvode->use_dense) { 4699566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumLinIters(cvode->mem, &its)); /* its = no. of calls to TSPrecond_Sundials() */ 47063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of linear iterations %ld\n", its)); 4719566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumConvFails(cvode->mem, &itmp)); 47263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of linear convergence failures %ld\n", itmp)); 47342528757SHong Zhang 4749566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc)); 4759566063dSJacob Faibussowitsch PetscCall(PCView(pc, viewer)); 4769566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumPrecEvals(cvode->mem, &itmp)); 47763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of preconditioner evaluations %ld\n", itmp)); 4789566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumPrecSolves(cvode->mem, &itmp)); 47963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of preconditioner solves %ld\n", itmp)); 480918687b8SPatrick Sanan } 4819566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumJtimesEvals(cvode->mem, &itmp)); 48263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of Jacobian-vector product evaluations %ld\n", itmp)); 4839566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumRhsEvals(cvode->mem, &itmp)); 48463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of rhs calls for finite diff. Jacobian-vector evals %ld\n", itmp)); 4857cb94ee6SHong Zhang } else if (isstring) { 4869566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer, "Sundials type %s", type)); 4877cb94ee6SHong Zhang } 4887cb94ee6SHong Zhang PetscFunctionReturn(0); 4897cb94ee6SHong Zhang } 4907cb94ee6SHong Zhang 4917cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 492*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetType_Sundials(TS ts, TSSundialsLmmType type) 493*d71ae5a4SJacob Faibussowitsch { 4947cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 4957cb94ee6SHong Zhang 4967cb94ee6SHong Zhang PetscFunctionBegin; 4977cb94ee6SHong Zhang cvode->cvode_type = type; 4987cb94ee6SHong Zhang PetscFunctionReturn(0); 4997cb94ee6SHong Zhang } 5007cb94ee6SHong Zhang 501*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts, PetscInt maxl) 502*d71ae5a4SJacob Faibussowitsch { 5037cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 5047cb94ee6SHong Zhang 5057cb94ee6SHong Zhang PetscFunctionBegin; 506f61b2b6aSHong Zhang cvode->maxl = maxl; 5077cb94ee6SHong Zhang PetscFunctionReturn(0); 5087cb94ee6SHong Zhang } 5097cb94ee6SHong Zhang 510*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts, PetscReal tol) 511*d71ae5a4SJacob Faibussowitsch { 5127cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 5137cb94ee6SHong Zhang 5147cb94ee6SHong Zhang PetscFunctionBegin; 5157cb94ee6SHong Zhang cvode->linear_tol = tol; 5167cb94ee6SHong Zhang PetscFunctionReturn(0); 5177cb94ee6SHong Zhang } 5187cb94ee6SHong Zhang 519*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts, TSSundialsGramSchmidtType type) 520*d71ae5a4SJacob Faibussowitsch { 5217cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 5227cb94ee6SHong Zhang 5237cb94ee6SHong Zhang PetscFunctionBegin; 5247cb94ee6SHong Zhang cvode->gtype = type; 5257cb94ee6SHong Zhang PetscFunctionReturn(0); 5267cb94ee6SHong Zhang } 5277cb94ee6SHong Zhang 528*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts, PetscReal aabs, PetscReal rel) 529*d71ae5a4SJacob Faibussowitsch { 5307cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 5317cb94ee6SHong Zhang 5327cb94ee6SHong Zhang PetscFunctionBegin; 5337cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 5347cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 5357cb94ee6SHong Zhang PetscFunctionReturn(0); 5367cb94ee6SHong Zhang } 5377cb94ee6SHong Zhang 538*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts, PetscReal mindt) 539*d71ae5a4SJacob Faibussowitsch { 540f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials *)ts->data; 541f1cd61daSJed Brown 542f1cd61daSJed Brown PetscFunctionBegin; 543f1cd61daSJed Brown cvode->mindt = mindt; 544f1cd61daSJed Brown PetscFunctionReturn(0); 545f1cd61daSJed Brown } 546f1cd61daSJed Brown 547*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts, PetscReal maxdt) 548*d71ae5a4SJacob Faibussowitsch { 549f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials *)ts->data; 550f1cd61daSJed Brown 551f1cd61daSJed Brown PetscFunctionBegin; 552f1cd61daSJed Brown cvode->maxdt = maxdt; 553f1cd61daSJed Brown PetscFunctionReturn(0); 554f1cd61daSJed Brown } 555918687b8SPatrick Sanan 556*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetUseDense_Sundials(TS ts, PetscBool use_dense) 557*d71ae5a4SJacob Faibussowitsch { 558918687b8SPatrick Sanan TS_Sundials *cvode = (TS_Sundials *)ts->data; 559918687b8SPatrick Sanan 560918687b8SPatrick Sanan PetscFunctionBegin; 561918687b8SPatrick Sanan cvode->use_dense = use_dense; 562918687b8SPatrick Sanan PetscFunctionReturn(0); 563918687b8SPatrick Sanan } 564918687b8SPatrick Sanan 565*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsGetPC_Sundials(TS ts, PC *pc) 566*d71ae5a4SJacob Faibussowitsch { 567dcbc6d53SSean Farley SNES snes; 568dcbc6d53SSean Farley KSP ksp; 5697cb94ee6SHong Zhang 5707cb94ee6SHong Zhang PetscFunctionBegin; 5719566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 5729566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 5739566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, pc)); 5747cb94ee6SHong Zhang PetscFunctionReturn(0); 5757cb94ee6SHong Zhang } 5767cb94ee6SHong Zhang 577*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsGetIterations_Sundials(TS ts, int *nonlin, int *lin) 578*d71ae5a4SJacob Faibussowitsch { 5797cb94ee6SHong Zhang PetscFunctionBegin; 5805ef26d82SJed Brown if (nonlin) *nonlin = ts->snes_its; 5815ef26d82SJed Brown if (lin) *lin = ts->ksp_its; 5827cb94ee6SHong Zhang PetscFunctionReturn(0); 5837cb94ee6SHong Zhang } 5847cb94ee6SHong Zhang 585*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts, PetscBool s) 586*d71ae5a4SJacob Faibussowitsch { 5872bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 5882bfc04deSHong Zhang 5892bfc04deSHong Zhang PetscFunctionBegin; 5902bfc04deSHong Zhang cvode->monitorstep = s; 5912bfc04deSHong Zhang PetscFunctionReturn(0); 5922bfc04deSHong Zhang } 5937cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 5947cb94ee6SHong Zhang 5957cb94ee6SHong Zhang /*@C 5967cb94ee6SHong Zhang TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 5977cb94ee6SHong Zhang 5987cb94ee6SHong Zhang Not Collective 5997cb94ee6SHong Zhang 60097bb3fdcSJose E. Roman Input Parameter: 6017cb94ee6SHong Zhang . ts - the time-step context 6027cb94ee6SHong Zhang 6037cb94ee6SHong Zhang Output Parameters: 6047cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 6057cb94ee6SHong Zhang - lin - number of linear iterations 6067cb94ee6SHong Zhang 6077cb94ee6SHong Zhang Level: advanced 6087cb94ee6SHong Zhang 6097cb94ee6SHong Zhang Notes: 6107cb94ee6SHong Zhang These return the number since the creation of the TS object 6117cb94ee6SHong Zhang 612db781477SPatrick Sanan .seealso: `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 613db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 614db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 615db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsGetPC()`, `TSSetExactFinalTime()` 6167cb94ee6SHong Zhang 6177cb94ee6SHong Zhang @*/ 618*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsGetIterations(TS ts, int *nonlin, int *lin) 619*d71ae5a4SJacob Faibussowitsch { 6207cb94ee6SHong Zhang PetscFunctionBegin; 621cac4c232SBarry Smith PetscUseMethod(ts, "TSSundialsGetIterations_C", (TS, int *, int *), (ts, nonlin, lin)); 6227cb94ee6SHong Zhang PetscFunctionReturn(0); 6237cb94ee6SHong Zhang } 6247cb94ee6SHong Zhang 6257cb94ee6SHong Zhang /*@ 6267cb94ee6SHong Zhang TSSundialsSetType - Sets the method that Sundials will use for integration. 6277cb94ee6SHong Zhang 6283f9fe445SBarry Smith Logically Collective on TS 6297cb94ee6SHong Zhang 6308f7d6fe5SPatrick Sanan Input Parameters: 6317cb94ee6SHong Zhang + ts - the time-step context 6327cb94ee6SHong Zhang - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 6337cb94ee6SHong Zhang 6347cb94ee6SHong Zhang Level: intermediate 6357cb94ee6SHong Zhang 636db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetMaxl()`, 637db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 638db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 639db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 640db781477SPatrick Sanan `TSSetExactFinalTime()` 6417cb94ee6SHong Zhang @*/ 642*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetType(TS ts, TSSundialsLmmType type) 643*d71ae5a4SJacob Faibussowitsch { 6447cb94ee6SHong Zhang PetscFunctionBegin; 645cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetType_C", (TS, TSSundialsLmmType), (ts, type)); 6467cb94ee6SHong Zhang PetscFunctionReturn(0); 6477cb94ee6SHong Zhang } 6487cb94ee6SHong Zhang 6497cb94ee6SHong Zhang /*@ 650c4e80e11SFlorian TSSundialsSetMaxord - Sets the maximum order for BDF/Adams method used by SUNDIALS. 651c4e80e11SFlorian 652c4e80e11SFlorian Logically Collective on TS 653c4e80e11SFlorian 6548f7d6fe5SPatrick Sanan Input Parameters: 655c4e80e11SFlorian + ts - the time-step context 656c4e80e11SFlorian - maxord - maximum order of BDF / Adams method 657c4e80e11SFlorian 658c4e80e11SFlorian Level: advanced 659c4e80e11SFlorian 660db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, 661db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 662db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 663db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 664db781477SPatrick Sanan `TSSetExactFinalTime()` 665c4e80e11SFlorian 666c4e80e11SFlorian @*/ 667*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetMaxord(TS ts, PetscInt maxord) 668*d71ae5a4SJacob Faibussowitsch { 669c4e80e11SFlorian PetscFunctionBegin; 670c4e80e11SFlorian PetscValidLogicalCollectiveInt(ts, maxord, 2); 671cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetMaxOrd_C", (TS, PetscInt), (ts, maxord)); 672c4e80e11SFlorian PetscFunctionReturn(0); 673c4e80e11SFlorian } 674c4e80e11SFlorian 675c4e80e11SFlorian /*@ 676f61b2b6aSHong Zhang TSSundialsSetMaxl - Sets the dimension of the Krylov space used by 6777cb94ee6SHong Zhang GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 678f61b2b6aSHong Zhang this is the maximum number of GMRES steps that will be used. 6797cb94ee6SHong Zhang 6803f9fe445SBarry Smith Logically Collective on TS 6817cb94ee6SHong Zhang 6828f7d6fe5SPatrick Sanan Input Parameters: 6837cb94ee6SHong Zhang + ts - the time-step context 684f61b2b6aSHong Zhang - maxl - number of direction vectors (the dimension of Krylov subspace). 6857cb94ee6SHong Zhang 6867cb94ee6SHong Zhang Level: advanced 6877cb94ee6SHong Zhang 688db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, 689db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 690db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 691db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 692db781477SPatrick Sanan `TSSetExactFinalTime()` 6937cb94ee6SHong Zhang 6947cb94ee6SHong Zhang @*/ 695*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetMaxl(TS ts, PetscInt maxl) 696*d71ae5a4SJacob Faibussowitsch { 6977cb94ee6SHong Zhang PetscFunctionBegin; 698f61b2b6aSHong Zhang PetscValidLogicalCollectiveInt(ts, maxl, 2); 699cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetMaxl_C", (TS, PetscInt), (ts, maxl)); 7007cb94ee6SHong Zhang PetscFunctionReturn(0); 7017cb94ee6SHong Zhang } 7027cb94ee6SHong Zhang 7037cb94ee6SHong Zhang /*@ 7047cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 7057cb94ee6SHong Zhang system by SUNDIALS. 7067cb94ee6SHong Zhang 7073f9fe445SBarry Smith Logically Collective on TS 7087cb94ee6SHong Zhang 7098f7d6fe5SPatrick Sanan Input Parameters: 7107cb94ee6SHong Zhang + ts - the time-step context 7117cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 7127cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 7137cb94ee6SHong Zhang 7147cb94ee6SHong Zhang Level: advanced 7157cb94ee6SHong Zhang 716db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 717db781477SPatrick Sanan `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 718db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 719db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 720db781477SPatrick Sanan `TSSetExactFinalTime()` 7217cb94ee6SHong Zhang 7227cb94ee6SHong Zhang @*/ 723*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetLinearTolerance(TS ts, PetscReal tol) 724*d71ae5a4SJacob Faibussowitsch { 7257cb94ee6SHong Zhang PetscFunctionBegin; 726c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts, tol, 2); 727cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetLinearTolerance_C", (TS, PetscReal), (ts, tol)); 7287cb94ee6SHong Zhang PetscFunctionReturn(0); 7297cb94ee6SHong Zhang } 7307cb94ee6SHong Zhang 7317cb94ee6SHong Zhang /*@ 7327cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 7337cb94ee6SHong Zhang in GMRES method by SUNDIALS linear solver. 7347cb94ee6SHong Zhang 7353f9fe445SBarry Smith Logically Collective on TS 7367cb94ee6SHong Zhang 7378f7d6fe5SPatrick Sanan Input Parameters: 7387cb94ee6SHong Zhang + ts - the time-step context 7397cb94ee6SHong Zhang - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 7407cb94ee6SHong Zhang 7417cb94ee6SHong Zhang Level: advanced 7427cb94ee6SHong Zhang 743db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 744db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, 745db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 746db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 747db781477SPatrick Sanan `TSSetExactFinalTime()` 7487cb94ee6SHong Zhang 7497cb94ee6SHong Zhang @*/ 750*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetGramSchmidtType(TS ts, TSSundialsGramSchmidtType type) 751*d71ae5a4SJacob Faibussowitsch { 7527cb94ee6SHong Zhang PetscFunctionBegin; 753cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetGramSchmidtType_C", (TS, TSSundialsGramSchmidtType), (ts, type)); 7547cb94ee6SHong Zhang PetscFunctionReturn(0); 7557cb94ee6SHong Zhang } 7567cb94ee6SHong Zhang 7577cb94ee6SHong Zhang /*@ 7587cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 7597cb94ee6SHong Zhang Sundials for error control. 7607cb94ee6SHong Zhang 7613f9fe445SBarry Smith Logically Collective on TS 7627cb94ee6SHong Zhang 7638f7d6fe5SPatrick Sanan Input Parameters: 7647cb94ee6SHong Zhang + ts - the time-step context 7657cb94ee6SHong Zhang . aabs - the absolute tolerance 7667cb94ee6SHong Zhang - rel - the relative tolerance 7677cb94ee6SHong Zhang 7687cb94ee6SHong Zhang See the Cvode/Sundials users manual for exact details on these parameters. Essentially 7697cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 7707cb94ee6SHong Zhang 7717cb94ee6SHong Zhang Level: intermediate 7727cb94ee6SHong Zhang 773db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetGMRESMaxl()`, 774db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, 775db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 776db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 777db781477SPatrick Sanan `TSSetExactFinalTime()` 7787cb94ee6SHong Zhang 7797cb94ee6SHong Zhang @*/ 780*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetTolerance(TS ts, PetscReal aabs, PetscReal rel) 781*d71ae5a4SJacob Faibussowitsch { 7827cb94ee6SHong Zhang PetscFunctionBegin; 783cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetTolerance_C", (TS, PetscReal, PetscReal), (ts, aabs, rel)); 7847cb94ee6SHong Zhang PetscFunctionReturn(0); 7857cb94ee6SHong Zhang } 7867cb94ee6SHong Zhang 7877cb94ee6SHong Zhang /*@ 7887cb94ee6SHong Zhang TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 7897cb94ee6SHong Zhang 7907cb94ee6SHong Zhang Input Parameter: 7917cb94ee6SHong Zhang . ts - the time-step context 7927cb94ee6SHong Zhang 7937cb94ee6SHong Zhang Output Parameter: 7947cb94ee6SHong Zhang . pc - the preconditioner context 7957cb94ee6SHong Zhang 7967cb94ee6SHong Zhang Level: advanced 7977cb94ee6SHong Zhang 798db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 799db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 800db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 801db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()` 8027cb94ee6SHong Zhang @*/ 803*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsGetPC(TS ts, PC *pc) 804*d71ae5a4SJacob Faibussowitsch { 8057cb94ee6SHong Zhang PetscFunctionBegin; 806cac4c232SBarry Smith PetscUseMethod(ts, "TSSundialsGetPC_C", (TS, PC *), (ts, pc)); 8077cb94ee6SHong Zhang PetscFunctionReturn(0); 8087cb94ee6SHong Zhang } 8097cb94ee6SHong Zhang 810f1cd61daSJed Brown /*@ 811f1cd61daSJed Brown TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 812f1cd61daSJed Brown 813d8d19677SJose E. Roman Input Parameters: 814f1cd61daSJed Brown + ts - the time-step context 815f1cd61daSJed Brown - mindt - lowest time step if positive, negative to deactivate 816f1cd61daSJed Brown 817fc6b9e64SJed Brown Note: 818fc6b9e64SJed Brown Sundials will error if it is not possible to keep the estimated truncation error below 819fc6b9e64SJed Brown the tolerance set with TSSundialsSetTolerance() without going below this step size. 820fc6b9e64SJed Brown 821f1cd61daSJed Brown Level: beginner 822f1cd61daSJed Brown 823db781477SPatrick Sanan .seealso: `TSSundialsSetType()`, `TSSundialsSetTolerance()`, 824f1cd61daSJed Brown @*/ 825*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetMinTimeStep(TS ts, PetscReal mindt) 826*d71ae5a4SJacob Faibussowitsch { 827f1cd61daSJed Brown PetscFunctionBegin; 828cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetMinTimeStep_C", (TS, PetscReal), (ts, mindt)); 829f1cd61daSJed Brown PetscFunctionReturn(0); 830f1cd61daSJed Brown } 831f1cd61daSJed Brown 832f1cd61daSJed Brown /*@ 833f1cd61daSJed Brown TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 834f1cd61daSJed Brown 835d8d19677SJose E. Roman Input Parameters: 836f1cd61daSJed Brown + ts - the time-step context 837f1cd61daSJed Brown - maxdt - lowest time step if positive, negative to deactivate 838f1cd61daSJed Brown 839f1cd61daSJed Brown Level: beginner 840f1cd61daSJed Brown 841db781477SPatrick Sanan .seealso: `TSSundialsSetType()`, `TSSundialsSetTolerance()`, 842f1cd61daSJed Brown @*/ 843*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetMaxTimeStep(TS ts, PetscReal maxdt) 844*d71ae5a4SJacob Faibussowitsch { 845f1cd61daSJed Brown PetscFunctionBegin; 846cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt)); 847f1cd61daSJed Brown PetscFunctionReturn(0); 848f1cd61daSJed Brown } 849f1cd61daSJed Brown 8502bfc04deSHong Zhang /*@ 8512bfc04deSHong Zhang TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 8522bfc04deSHong Zhang 853d8d19677SJose E. Roman Input Parameters: 8542bfc04deSHong Zhang + ts - the time-step context 8552bfc04deSHong Zhang - ft - PETSC_TRUE if monitor, else PETSC_FALSE 8562bfc04deSHong Zhang 8572bfc04deSHong Zhang Level: beginner 8582bfc04deSHong Zhang 859f61b2b6aSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 8602bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 861f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 8622bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 8632bfc04deSHong Zhang @*/ 864*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsMonitorInternalSteps(TS ts, PetscBool ft) 865*d71ae5a4SJacob Faibussowitsch { 8662bfc04deSHong Zhang PetscFunctionBegin; 867cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsMonitorInternalSteps_C", (TS, PetscBool), (ts, ft)); 8682bfc04deSHong Zhang PetscFunctionReturn(0); 8692bfc04deSHong Zhang } 870918687b8SPatrick Sanan 871918687b8SPatrick Sanan /*@ 872918687b8SPatrick Sanan TSSundialsSetUseDense - Set a flag to use a dense linear solver in SUNDIALS (serial only) 873918687b8SPatrick Sanan 874918687b8SPatrick Sanan Logically Collective 875918687b8SPatrick Sanan 876918687b8SPatrick Sanan Input Parameters: 877918687b8SPatrick Sanan + ts - the time-step context 878918687b8SPatrick Sanan - use_dense - PETSC_TRUE to use the dense solver 879918687b8SPatrick Sanan 880918687b8SPatrick Sanan Level: advanced 881918687b8SPatrick Sanan 882db781477SPatrick Sanan .seealso: `TSSUNDIALS` 883918687b8SPatrick Sanan 884918687b8SPatrick Sanan @*/ 885*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetUseDense(TS ts, PetscBool use_dense) 886*d71ae5a4SJacob Faibussowitsch { 887918687b8SPatrick Sanan PetscFunctionBegin; 888918687b8SPatrick Sanan PetscValidLogicalCollectiveInt(ts, use_dense, 2); 889cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetUseDense_C", (TS, PetscBool), (ts, use_dense)); 890918687b8SPatrick Sanan PetscFunctionReturn(0); 891918687b8SPatrick Sanan } 892918687b8SPatrick Sanan 8937cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 8947cb94ee6SHong Zhang /*MC 89596f5712cSJed Brown TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 8967cb94ee6SHong Zhang 8977cb94ee6SHong Zhang Options Database: 898897c9f78SHong Zhang + -ts_sundials_type <bdf,adams> - 8997cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 9007cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 9017cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 902897c9f78SHong Zhang . -ts_sundials_linear_tolerance <tol> - 903f61b2b6aSHong Zhang . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace 904918687b8SPatrick Sanan . -ts_sundials_monitor_steps - Monitor SUNDIALS internal steps 905918687b8SPatrick Sanan - -ts_sundials_use_dense - Use a dense linear solver within CVODE (serial only) 90616016685SHong Zhang 90795452b02SPatrick Sanan Notes: 90895452b02SPatrick Sanan This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply, 909897c9f78SHong Zhang only PETSc PC options. 9107cb94ee6SHong Zhang 9117cb94ee6SHong Zhang Level: beginner 9127cb94ee6SHong Zhang 913db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, `TSSundialsSetLinearTolerance()`, 914db781477SPatrick Sanan `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, `TSSundialsGetIterations()`, `TSSetExactFinalTime()` 9157cb94ee6SHong Zhang 9167cb94ee6SHong Zhang M*/ 917*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts) 918*d71ae5a4SJacob Faibussowitsch { 9197cb94ee6SHong Zhang TS_Sundials *cvode; 92042528757SHong Zhang PC pc; 9217cb94ee6SHong Zhang 9227cb94ee6SHong Zhang PetscFunctionBegin; 923277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Sundials; 92428adb3f7SHong Zhang ts->ops->destroy = TSDestroy_Sundials; 92528adb3f7SHong Zhang ts->ops->view = TSView_Sundials; 926214bc6a2SJed Brown ts->ops->setup = TSSetUp_Sundials; 927b4eba00bSSean Farley ts->ops->step = TSStep_Sundials; 928b4eba00bSSean Farley ts->ops->interpolate = TSInterpolate_Sundials; 929214bc6a2SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Sundials; 9302ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 9317cb94ee6SHong Zhang 9324dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&cvode)); 933bbd56ea5SKarl Rupp 934efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 935efd4aadfSBarry Smith 9367cb94ee6SHong Zhang ts->data = (void *)cvode; 9376fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF; 9386fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS; 939f61b2b6aSHong Zhang cvode->maxl = 5; 940c4e80e11SFlorian cvode->maxord = PETSC_DEFAULT; 9417cb94ee6SHong Zhang cvode->linear_tol = .05; 942b4eba00bSSean Farley cvode->monitorstep = PETSC_TRUE; 943918687b8SPatrick Sanan cvode->use_dense = PETSC_FALSE; 9447cb94ee6SHong Zhang 9459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)ts), &(cvode->comm_sundials))); 946f1cd61daSJed Brown 947f1cd61daSJed Brown cvode->mindt = -1.; 948f1cd61daSJed Brown cvode->maxdt = -1.; 949f1cd61daSJed Brown 9507cb94ee6SHong Zhang /* set tolerance for Sundials */ 9517cb94ee6SHong Zhang cvode->reltol = 1e-6; 9522c823083SHong Zhang cvode->abstol = 1e-6; 9537cb94ee6SHong Zhang 95442528757SHong Zhang /* set PCNONE as default pctype */ 9559566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC_Sundials(ts, &pc)); 9569566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 95742528757SHong Zhang 9589566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetType_C", TSSundialsSetType_Sundials)); 9599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxl_C", TSSundialsSetMaxl_Sundials)); 9609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetLinearTolerance_C", TSSundialsSetLinearTolerance_Sundials)); 9619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetGramSchmidtType_C", TSSundialsSetGramSchmidtType_Sundials)); 9629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetTolerance_C", TSSundialsSetTolerance_Sundials)); 9639566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMinTimeStep_C", TSSundialsSetMinTimeStep_Sundials)); 9649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxTimeStep_C", TSSundialsSetMaxTimeStep_Sundials)); 9659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetPC_C", TSSundialsGetPC_Sundials)); 9669566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetIterations_C", TSSundialsGetIterations_Sundials)); 9679566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsMonitorInternalSteps_C", TSSundialsMonitorInternalSteps_Sundials)); 9689566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetUseDense_C", TSSundialsSetUseDense_Sundials)); 9697cb94ee6SHong Zhang PetscFunctionReturn(0); 9707cb94ee6SHong Zhang } 971