17cb94ee6SHong Zhang /* 2bcf0153eSBarry Smith Provides a PETSc interface to version 2.5 of SUNDIALS/CVODE solver (a very old version) 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 */ 143ba16761SJacob Faibussowitsch static PetscErrorCode TSPrecond_Sundials_Petsc(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) 15d71ae5a4SJacob 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)); 32bcf0153eSBarry Smith 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)); 363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 373ba16761SJacob Faibussowitsch } 383ba16761SJacob Faibussowitsch 393ba16761SJacob Faibussowitsch /* Sundial expects an int (*)(args...) but PetscErrorCode is an enum. Instead of switching out 403ba16761SJacob Faibussowitsch all the PetscCalls in TSPrecond_Sundials_Petsc we just wrap it */ 413ba16761SJacob Faibussowitsch static int TSPrecond_Sundials_Private(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) 423ba16761SJacob Faibussowitsch { 433ba16761SJacob Faibussowitsch return (int)TSPrecond_Sundials_Petsc(tn, y, fy, jok, jcurPtr, _gamma, P_data, vtemp1, vtemp2, vtemp3); 447cb94ee6SHong Zhang } 457cb94ee6SHong Zhang 467cb94ee6SHong Zhang /* 47bcf0153eSBarry Smith TSPSolve_Sundials - routine that we provide to SUNDIALS that applies the preconditioner. 487cb94ee6SHong Zhang */ 493ba16761SJacob Faibussowitsch static PetscErrorCode TSPSolve_Sundials_Petsc(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) 50d71ae5a4SJacob Faibussowitsch { 517cb94ee6SHong Zhang TS ts = (TS)P_data; 527cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 53dcbc6d53SSean Farley PC pc; 547cb94ee6SHong Zhang Vec rr = cvode->w1, zz = cvode->w2; 557cb94ee6SHong Zhang PetscScalar *r_data, *z_data; 567cb94ee6SHong Zhang 577cb94ee6SHong Zhang PetscFunctionBegin; 587cb94ee6SHong Zhang /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/ 597cb94ee6SHong Zhang r_data = (PetscScalar *)N_VGetArrayPointer(r); 607cb94ee6SHong Zhang z_data = (PetscScalar *)N_VGetArrayPointer(z); 619566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(rr, r_data)); 629566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(zz, z_data)); 634ac7836bSHong Zhang 647cb94ee6SHong Zhang /* Solve the Px=r and put the result in zz */ 659566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc)); 669566063dSJacob Faibussowitsch PetscCall(PCApply(pc, rr, zz)); 679566063dSJacob Faibussowitsch PetscCall(VecResetArray(rr)); 689566063dSJacob Faibussowitsch PetscCall(VecResetArray(zz)); 693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 703ba16761SJacob Faibussowitsch } 713ba16761SJacob Faibussowitsch 723ba16761SJacob Faibussowitsch /* See TSPrecond_Sundials_Private() */ 733ba16761SJacob Faibussowitsch static int TSPSolve_Sundials_Private(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) 743ba16761SJacob Faibussowitsch { 753ba16761SJacob Faibussowitsch return (int)TSPSolve_Sundials_Petsc(tn, y, fy, r, z, _gamma, delta, lr, P_data, vtemp); 767cb94ee6SHong Zhang } 777cb94ee6SHong Zhang 787cb94ee6SHong Zhang /* 79bcf0153eSBarry Smith TSFunction_Sundials - routine that we provide to SUNDIALS that applies the right hand side. 807cb94ee6SHong Zhang */ 8166976f2fSJacob Faibussowitsch static int TSFunction_Sundials(realtype t, N_Vector y, N_Vector ydot, void *ctx) 82d71ae5a4SJacob Faibussowitsch { 837cb94ee6SHong Zhang TS ts = (TS)ctx; 845fcef5e4SJed Brown DM dm; 85942e3340SBarry Smith DMTS tsdm; 868434afd1SBarry Smith TSIFunctionFn *ifunction; 87ce94432eSBarry Smith MPI_Comm comm; 887cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 890679a0aeSJed Brown Vec yy = cvode->w1, yyd = cvode->w2, yydot = cvode->ydot; 907cb94ee6SHong Zhang PetscScalar *y_data, *ydot_data; 917cb94ee6SHong Zhang 927cb94ee6SHong Zhang PetscFunctionBegin; 939566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 945ff03cf7SDinesh Kaushik /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/ 957cb94ee6SHong Zhang y_data = (PetscScalar *)N_VGetArrayPointer(y); 967cb94ee6SHong Zhang ydot_data = (PetscScalar *)N_VGetArrayPointer(ydot); 979566063dSJacob Faibussowitsch PetscCallAbort(comm, VecPlaceArray(yy, y_data)); 989566063dSJacob Faibussowitsch PetscCallAbort(comm, VecPlaceArray(yyd, ydot_data)); 994ac7836bSHong Zhang 1005fcef5e4SJed Brown /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */ 1019566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 1029566063dSJacob Faibussowitsch PetscCall(DMGetDMTS(dm, &tsdm)); 1039566063dSJacob Faibussowitsch PetscCall(DMTSGetIFunction(dm, &ifunction, NULL)); 1045fcef5e4SJed Brown if (!ifunction) { 1059566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, t, yy, yyd)); 106bc0cc02bSJed Brown } else { /* If rhsfunction is also set, this computes both parts and shifts them to the right */ 1079566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(yydot)); 1089566063dSJacob Faibussowitsch PetscCallAbort(comm, TSComputeIFunction(ts, t, yy, yydot, yyd, PETSC_FALSE)); 1099566063dSJacob Faibussowitsch PetscCall(VecScale(yyd, -1.)); 110bc0cc02bSJed Brown } 1119566063dSJacob Faibussowitsch PetscCallAbort(comm, VecResetArray(yy)); 1129566063dSJacob Faibussowitsch PetscCallAbort(comm, VecResetArray(yyd)); 1133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1147cb94ee6SHong Zhang } 1157cb94ee6SHong Zhang 1167cb94ee6SHong Zhang /* 117bcf0153eSBarry Smith TSStep_Sundials - Calls SUNDIALS to integrate the ODE. 1187cb94ee6SHong Zhang */ 11966976f2fSJacob Faibussowitsch static PetscErrorCode TSStep_Sundials(TS ts) 120d71ae5a4SJacob Faibussowitsch { 1217cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 122b4eba00bSSean Farley PetscInt flag; 123be5899b3SLisandro Dalcin long int nits, lits, nsteps; 1247cb94ee6SHong Zhang realtype t, tout; 1257cb94ee6SHong Zhang PetscScalar *y_data; 1267cb94ee6SHong Zhang void *mem; 1277cb94ee6SHong Zhang 1287cb94ee6SHong Zhang PetscFunctionBegin; 12916016685SHong Zhang mem = cvode->mem; 1307cb94ee6SHong Zhang tout = ts->max_time; 1319566063dSJacob Faibussowitsch PetscCall(VecGetArray(ts->vec_sol, &y_data)); 1327cb94ee6SHong Zhang N_VSetArrayPointer((realtype *)y_data, cvode->y); 1339566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ts->vec_sol, NULL)); 134186e87acSLisandro Dalcin 1359687d888SLisandro Dalcin /* We would like to TSPreStage() and TSPostStage() 1369687d888SLisandro Dalcin * before each stage solve but CVode does not appear to support this. */ 1379371c9d4SSatish Balay if (cvode->monitorstep) flag = CVode(mem, tout, cvode->y, &t, CV_ONE_STEP); 1389371c9d4SSatish Balay else flag = CVode(mem, tout, cvode->y, &t, CV_NORMAL); 1399f94935aSHong Zhang 1409f94935aSHong Zhang if (flag) { /* display error message */ 1419f94935aSHong Zhang switch (flag) { 142d71ae5a4SJacob Faibussowitsch case CV_ILL_INPUT: 143d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_ILL_INPUT"); 144d71ae5a4SJacob Faibussowitsch break; 145d71ae5a4SJacob Faibussowitsch case CV_TOO_CLOSE: 146d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_TOO_CLOSE"); 147d71ae5a4SJacob Faibussowitsch break; 1483c7fefeeSJed Brown case CV_TOO_MUCH_WORK: { 1499f94935aSHong Zhang PetscReal tcur; 1503ba16761SJacob Faibussowitsch PetscCallExternal(CVodeGetNumSteps, mem, &nsteps); 1513ba16761SJacob Faibussowitsch PetscCallExternal(CVodeGetCurrentTime, mem, &tcur); 15263a3b9bcSJacob 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); 1533c7fefeeSJed Brown } break; 154d71ae5a4SJacob Faibussowitsch case CV_TOO_MUCH_ACC: 155d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_TOO_MUCH_ACC"); 156d71ae5a4SJacob Faibussowitsch break; 157d71ae5a4SJacob Faibussowitsch case CV_ERR_FAILURE: 158d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_ERR_FAILURE"); 159d71ae5a4SJacob Faibussowitsch break; 160d71ae5a4SJacob Faibussowitsch case CV_CONV_FAILURE: 161d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_CONV_FAILURE"); 162d71ae5a4SJacob Faibussowitsch break; 163d71ae5a4SJacob Faibussowitsch case CV_LINIT_FAIL: 164d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LINIT_FAIL"); 165d71ae5a4SJacob Faibussowitsch break; 166d71ae5a4SJacob Faibussowitsch case CV_LSETUP_FAIL: 167d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LSETUP_FAIL"); 168d71ae5a4SJacob Faibussowitsch break; 169d71ae5a4SJacob Faibussowitsch case CV_LSOLVE_FAIL: 170d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LSOLVE_FAIL"); 171d71ae5a4SJacob Faibussowitsch break; 172d71ae5a4SJacob Faibussowitsch case CV_RHSFUNC_FAIL: 173d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_RHSFUNC_FAIL"); 174d71ae5a4SJacob Faibussowitsch break; 175d71ae5a4SJacob Faibussowitsch case CV_FIRST_RHSFUNC_ERR: 176d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_FIRST_RHSFUNC_ERR"); 177d71ae5a4SJacob Faibussowitsch break; 178d71ae5a4SJacob Faibussowitsch case CV_REPTD_RHSFUNC_ERR: 179d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_REPTD_RHSFUNC_ERR"); 180d71ae5a4SJacob Faibussowitsch break; 181d71ae5a4SJacob Faibussowitsch case CV_UNREC_RHSFUNC_ERR: 182d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_UNREC_RHSFUNC_ERR"); 183d71ae5a4SJacob Faibussowitsch break; 184d71ae5a4SJacob Faibussowitsch case CV_RTFUNC_FAIL: 185d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_RTFUNC_FAIL"); 186d71ae5a4SJacob Faibussowitsch break; 187d71ae5a4SJacob Faibussowitsch default: 188d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, flag %d", flag); 1899f94935aSHong Zhang } 1909f94935aSHong Zhang } 1919f94935aSHong Zhang 192be5899b3SLisandro Dalcin /* log inner nonlinear and linear iterations */ 1933ba16761SJacob Faibussowitsch PetscCallExternal(CVodeGetNumNonlinSolvIters, mem, &nits); 1943ba16761SJacob Faibussowitsch PetscCallExternal(CVSpilsGetNumLinIters, mem, &lits); 1959371c9d4SSatish Balay ts->snes_its += nits; 1969371c9d4SSatish Balay ts->ksp_its = lits; 197be5899b3SLisandro Dalcin 1987cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */ 1999566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(cvode->w1, y_data)); 2009566063dSJacob Faibussowitsch PetscCall(VecCopy(cvode->w1, cvode->update)); 2019566063dSJacob Faibussowitsch PetscCall(VecResetArray(cvode->w1)); 2029566063dSJacob Faibussowitsch PetscCall(VecCopy(cvode->update, ts->vec_sol)); 203186e87acSLisandro Dalcin 204186e87acSLisandro Dalcin ts->time_step = t - ts->ptime; 205186e87acSLisandro Dalcin ts->ptime = t; 206186e87acSLisandro Dalcin 2073ba16761SJacob Faibussowitsch PetscCallExternal(CVodeGetNumSteps, mem, &nsteps); 2082808aa04SLisandro Dalcin if (!cvode->monitorstep) ts->steps += nsteps - 1; /* TSStep() increments the step counter by one */ 2093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 210b4eba00bSSean Farley } 211b4eba00bSSean Farley 212d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_Sundials(TS ts, PetscReal t, Vec X) 213d71ae5a4SJacob Faibussowitsch { 214b4eba00bSSean Farley TS_Sundials *cvode = (TS_Sundials *)ts->data; 215b4eba00bSSean Farley N_Vector y; 216b4eba00bSSean Farley PetscScalar *x_data; 217b4eba00bSSean Farley PetscInt glosize, locsize; 218b4eba00bSSean Farley 219b4eba00bSSean Farley PetscFunctionBegin; 220b4eba00bSSean Farley /* get the vector size */ 2219566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &glosize)); 2229566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &locsize)); 2239566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x_data)); 224b4eba00bSSean Farley 2254d78c9e0SClaas Abert /* Initialize N_Vec y with x_data */ 226918687b8SPatrick Sanan if (cvode->use_dense) { 227918687b8SPatrick Sanan PetscMPIInt size; 228918687b8SPatrick Sanan 2299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 2303c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "TSSUNDIALS only supports a dense solve in the serial case"); 231918687b8SPatrick Sanan y = N_VMake_Serial(locsize, (realtype *)x_data); 232918687b8SPatrick Sanan } else { 2334d78c9e0SClaas Abert y = N_VMake_Parallel(cvode->comm_sundials, locsize, glosize, (realtype *)x_data); 234918687b8SPatrick Sanan } 235918687b8SPatrick Sanan 2363c633725SBarry Smith PetscCheck(y, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Interpolated y is not allocated"); 237b4eba00bSSean Farley 2383ba16761SJacob Faibussowitsch PetscCallExternal(CVodeGetDky, cvode->mem, t, 0, y); 2399566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x_data)); 2403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2417cb94ee6SHong Zhang } 2427cb94ee6SHong Zhang 24366976f2fSJacob Faibussowitsch static PetscErrorCode TSReset_Sundials(TS ts) 244d71ae5a4SJacob Faibussowitsch { 245277b19d0SLisandro Dalcin TS_Sundials *cvode = (TS_Sundials *)ts->data; 246277b19d0SLisandro Dalcin 247277b19d0SLisandro Dalcin PetscFunctionBegin; 2489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->update)); 2499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->ydot)); 2509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->w1)); 2519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->w2)); 252bbd56ea5SKarl Rupp if (cvode->mem) CVodeFree(&cvode->mem); 2533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 254277b19d0SLisandro Dalcin } 255277b19d0SLisandro Dalcin 25666976f2fSJacob Faibussowitsch static PetscErrorCode TSDestroy_Sundials(TS ts) 257d71ae5a4SJacob Faibussowitsch { 2587cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 2597cb94ee6SHong Zhang 2607cb94ee6SHong Zhang PetscFunctionBegin; 2619566063dSJacob Faibussowitsch PetscCall(TSReset_Sundials(ts)); 262*f4f49eeaSPierre Jolivet PetscCallMPI(MPI_Comm_free(&cvode->comm_sundials)); 2639566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetType_C", NULL)); 2659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxl_C", NULL)); 2669566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetLinearTolerance_C", NULL)); 2679566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetGramSchmidtType_C", NULL)); 2689566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetTolerance_C", NULL)); 2699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMinTimeStep_C", NULL)); 2709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxTimeStep_C", NULL)); 2719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetPC_C", NULL)); 2729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetIterations_C", NULL)); 2739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsMonitorInternalSteps_C", NULL)); 2742e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetUseDense_C", NULL)); 2753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2767cb94ee6SHong Zhang } 2777cb94ee6SHong Zhang 27866976f2fSJacob Faibussowitsch static PetscErrorCode TSSetUp_Sundials(TS ts) 279d71ae5a4SJacob Faibussowitsch { 2807cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 2813ba16761SJacob Faibussowitsch PetscInt glosize, locsize, i; 2827cb94ee6SHong Zhang PetscScalar *y_data, *parray; 283dcbc6d53SSean Farley PC pc; 28419fd82e9SBarry Smith PCType pctype; 285ace3abfcSBarry Smith PetscBool pcnone; 2867cb94ee6SHong Zhang 2877cb94ee6SHong Zhang PetscFunctionBegin; 288bcf0153eSBarry Smith 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"); 289236c45dcSShri Abhyankar 2907cb94ee6SHong Zhang /* get the vector size */ 2919566063dSJacob Faibussowitsch PetscCall(VecGetSize(ts->vec_sol, &glosize)); 2929566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ts->vec_sol, &locsize)); 2937cb94ee6SHong Zhang 2947cb94ee6SHong Zhang /* allocate the memory for N_Vec y */ 295918687b8SPatrick Sanan if (cvode->use_dense) { 296918687b8SPatrick Sanan PetscMPIInt size; 297918687b8SPatrick Sanan 2989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 2993c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "TSSUNDIALS only supports a dense solve in the serial case"); 300918687b8SPatrick Sanan cvode->y = N_VNew_Serial(locsize); 301918687b8SPatrick Sanan } else { 302a07356b8SSatish Balay cvode->y = N_VNew_Parallel(cvode->comm_sundials, locsize, glosize); 303918687b8SPatrick Sanan } 3043c633725SBarry Smith PetscCheck(cvode->y, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "cvode->y is not allocated"); 3057cb94ee6SHong Zhang 30628adb3f7SHong Zhang /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 3079566063dSJacob Faibussowitsch PetscCall(VecGetArray(ts->vec_sol, &parray)); 3087cb94ee6SHong Zhang y_data = (PetscScalar *)N_VGetArrayPointer(cvode->y); 3097cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 3109566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ts->vec_sol, NULL)); 31142528757SHong Zhang 3129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &cvode->update)); 3139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &cvode->ydot)); 3147cb94ee6SHong Zhang 3157cb94ee6SHong Zhang /* 3167cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are 3177cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided 318bcf0153eSBarry Smith by SUNDIALS and set using VecPlaceArray(). 3197cb94ee6SHong Zhang */ 3209566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts), 1, locsize, PETSC_DECIDE, NULL, &cvode->w1)); 3219566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts), 1, locsize, PETSC_DECIDE, NULL, &cvode->w2)); 32216016685SHong Zhang 32316016685SHong Zhang /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 3243ba16761SJacob Faibussowitsch cvode->mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 3253ba16761SJacob Faibussowitsch PetscCheck(cvode->mem, PETSC_COMM_SELF, PETSC_ERR_MEM, "CVodeCreate() fails"); 32616016685SHong Zhang 32716016685SHong Zhang /* Set the pointer to user-defined data */ 3283ba16761SJacob Faibussowitsch PetscCallExternal(CVodeSetUserData, cvode->mem, ts); 32916016685SHong Zhang 330bcf0153eSBarry Smith /* SUNDIALS may choose to use a smaller initial step, but will never use a larger step. */ 3313ba16761SJacob Faibussowitsch PetscCallExternal(CVodeSetInitStep, cvode->mem, (realtype)ts->time_step); 332f1cd61daSJed Brown if (cvode->mindt > 0) { 3333ba16761SJacob Faibussowitsch int flag = CVodeSetMinStep(cvode->mem, (realtype)cvode->mindt); 3349f94935aSHong Zhang if (flag) { 33508401ef6SPierre Jolivet PetscCheck(flag != CV_MEM_NULL, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMinStep() failed, cvode_mem pointer is NULL"); 336f7d195e4SLawrence 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"); 337f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMinStep() failed"); 3389f94935aSHong Zhang } 339f1cd61daSJed Brown } 3403ba16761SJacob Faibussowitsch if (cvode->maxdt > 0) PetscCallExternal(CVodeSetMaxStep, cvode->mem, (realtype)cvode->maxdt); 341f1cd61daSJed Brown 34216016685SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the 343a5b23f4aSJose E. Roman * user's right hand side function in u'=f(t,u), the initial time T0, and 34416016685SHong Zhang * the initial dependent variable vector cvode->y */ 3453ba16761SJacob Faibussowitsch PetscCallExternal(CVodeInit, cvode->mem, TSFunction_Sundials, ts->ptime, cvode->y); 34616016685SHong Zhang 3479f94935aSHong Zhang /* specifies scalar relative and absolute tolerances */ 3483ba16761SJacob Faibussowitsch PetscCallExternal(CVodeSStolerances, cvode->mem, cvode->reltol, cvode->abstol); 34916016685SHong Zhang 350c4e80e11SFlorian /* Specify max order of BDF / ADAMS method */ 3513ba16761SJacob Faibussowitsch if (cvode->maxord != PETSC_DEFAULT) PetscCallExternal(CVodeSetMaxOrd, cvode->mem, cvode->maxord); 352c4e80e11SFlorian 3539f94935aSHong Zhang /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */ 3543ba16761SJacob Faibussowitsch PetscCallExternal(CVodeSetMaxNumSteps, cvode->mem, ts->max_steps); 35516016685SHong Zhang 356918687b8SPatrick Sanan if (cvode->use_dense) { 3573ba16761SJacob Faibussowitsch PetscCallExternal(CVDense, cvode->mem, locsize); 358918687b8SPatrick Sanan } else { 35916016685SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 36016016685SHong Zhang /* setup the ode integrator with the given preconditioner */ 3619566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc)); 3629566063dSJacob Faibussowitsch PetscCall(PCGetType(pc, &pctype)); 3639566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCNONE, &pcnone)); 36416016685SHong Zhang if (pcnone) { 3653ba16761SJacob Faibussowitsch PetscCallExternal(CVSpgmr, cvode->mem, PREC_NONE, 0); 36616016685SHong Zhang } else { 3673ba16761SJacob Faibussowitsch PetscCallExternal(CVSpgmr, cvode->mem, PREC_LEFT, cvode->maxl); 36816016685SHong Zhang 36916016685SHong Zhang /* Set preconditioner and solve routines Precond and PSolve, 37016016685SHong Zhang and the pointer to the user-defined block data */ 3713ba16761SJacob Faibussowitsch PetscCallExternal(CVSpilsSetPreconditioner, cvode->mem, TSPrecond_Sundials_Private, TSPSolve_Sundials_Private); 37216016685SHong Zhang } 373918687b8SPatrick Sanan } 3743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3757cb94ee6SHong Zhang } 3767cb94ee6SHong Zhang 3776fadb2cdSHong Zhang /* type of CVODE linear multistep method */ 378c793f718SLisandro Dalcin const char *const TSSundialsLmmTypes[] = {"", "ADAMS", "BDF", "TSSundialsLmmType", "SUNDIALS_", NULL}; 3796fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */ 380c793f718SLisandro Dalcin const char *const TSSundialsGramSchmidtTypes[] = {"", "MODIFIED", "CLASSICAL", "TSSundialsGramSchmidtType", "SUNDIALS_", NULL}; 381a04cf4d8SBarry Smith 38266976f2fSJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_Sundials(TS ts, PetscOptionItems *PetscOptionsObject) 383d71ae5a4SJacob Faibussowitsch { 3847cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 3857cb94ee6SHong Zhang int indx; 386ace3abfcSBarry Smith PetscBool flag; 3877bda3a07SJed Brown PC pc; 3887cb94ee6SHong Zhang 3897cb94ee6SHong Zhang PetscFunctionBegin; 390d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SUNDIALS ODE solver options"); 3919566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_sundials_type", "Scheme", "TSSundialsSetType", TSSundialsLmmTypes, 3, TSSundialsLmmTypes[cvode->cvode_type], &indx, &flag)); 3921baa6e33SBarry Smith if (flag) PetscCall(TSSundialsSetType(ts, (TSSundialsLmmType)indx)); 3939566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_sundials_gramschmidt_type", "Type of orthogonalization", "TSSundialsSetGramSchmidtType", TSSundialsGramSchmidtTypes, 3, TSSundialsGramSchmidtTypes[cvode->gtype], &indx, &flag)); 3941baa6e33SBarry Smith if (flag) PetscCall(TSSundialsSetGramSchmidtType(ts, (TSSundialsGramSchmidtType)indx)); 3959566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_atol", "Absolute tolerance for convergence", "TSSundialsSetTolerance", cvode->abstol, &cvode->abstol, NULL)); 3969566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_rtol", "Relative tolerance for convergence", "TSSundialsSetTolerance", cvode->reltol, &cvode->reltol, NULL)); 3979566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_mindt", "Minimum step size", "TSSundialsSetMinTimeStep", cvode->mindt, &cvode->mindt, NULL)); 3989566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_maxdt", "Maximum step size", "TSSundialsSetMaxTimeStep", cvode->maxdt, &cvode->maxdt, NULL)); 3999566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_linear_tolerance", "Convergence tolerance for linear solve", "TSSundialsSetLinearTolerance", cvode->linear_tol, &cvode->linear_tol, NULL)); 4009566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_sundials_maxord", "Max Order for BDF/Adams method", "TSSundialsSetMaxOrd", cvode->maxord, &cvode->maxord, NULL)); 4019566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_sundials_maxl", "Max dimension of the Krylov subspace", "TSSundialsSetMaxl", cvode->maxl, &cvode->maxl, NULL)); 4029566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_sundials_monitor_steps", "Monitor SUNDIALS internal steps", "TSSundialsMonitorInternalSteps", cvode->monitorstep, &cvode->monitorstep, NULL)); 4039566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_sundials_use_dense", "Use dense internal solver in SUNDIALS (serial only)", "TSSundialsSetUseDense", cvode->use_dense, &cvode->use_dense, NULL)); 404d0609cedSBarry Smith PetscOptionsHeadEnd(); 4059566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc)); 4069566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions(pc)); 4073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4087cb94ee6SHong Zhang } 4097cb94ee6SHong Zhang 41066976f2fSJacob Faibussowitsch static PetscErrorCode TSView_Sundials(TS ts, PetscViewer viewer) 411d71ae5a4SJacob Faibussowitsch { 4127cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 4137cb94ee6SHong Zhang char *type; 4147cb94ee6SHong Zhang char atype[] = "Adams"; 4157cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 416ace3abfcSBarry Smith PetscBool iascii, isstring; 4172c823083SHong Zhang long int nsteps, its, nfevals, nlinsetups, nfails, itmp; 4182c823083SHong Zhang PetscInt qlast, qcur; 4195d47aee6SHong Zhang PetscReal hinused, hlast, hcur, tcur, tolsfac; 42042528757SHong Zhang PC pc; 4217cb94ee6SHong Zhang 4227cb94ee6SHong Zhang PetscFunctionBegin; 423bbd56ea5SKarl Rupp if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype; 424bbd56ea5SKarl Rupp else type = btype; 4257cb94ee6SHong Zhang 4269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 4279566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 4287cb94ee6SHong Zhang if (iascii) { 429bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS integrator does not use SNES!\n")); 430bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS integrator type %s\n", type)); 431bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS maxord %" PetscInt_FMT "\n", cvode->maxord)); 432bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS abs tol %g rel tol %g\n", (double)cvode->abstol, (double)cvode->reltol)); 433918687b8SPatrick Sanan if (cvode->use_dense) { 434bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS integrator using a dense linear solve\n")); 435918687b8SPatrick Sanan } else { 436bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS linear solver tolerance factor %g\n", (double)cvode->linear_tol)); 437bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS max dimension of Krylov subspace %" PetscInt_FMT "\n", cvode->maxl)); 4387cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 439bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS using modified Gram-Schmidt for orthogonalization in GMRES\n")); 4407cb94ee6SHong Zhang } else { 441bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n")); 4427cb94ee6SHong Zhang } 443918687b8SPatrick Sanan } 444bcf0153eSBarry Smith if (cvode->mindt > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS minimum time step %g\n", (double)cvode->mindt)); 445bcf0153eSBarry Smith if (cvode->maxdt > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS maximum time step %g\n", (double)cvode->maxdt)); 4462c823083SHong Zhang 4475d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */ 4483ba16761SJacob Faibussowitsch PetscCallExternal(CVodeGetTolScaleFactor, cvode->mem, &tolsfac); 449bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS suggested factor for tolerance scaling %g\n", tolsfac)); 4503ba16761SJacob Faibussowitsch PetscCallExternal(CVodeGetIntegratorStats, cvode->mem, &nsteps, &nfevals, &nlinsetups, &nfails, &qlast, &qcur, &hinused, &hlast, &hcur, &tcur); 451bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS cumulative number of internal steps %ld\n", nsteps)); 452bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of calls to rhs function %ld\n", nfevals)); 453bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of calls to linear solver setup function %ld\n", nlinsetups)); 454bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of error test failures %ld\n", nfails)); 4552c823083SHong Zhang 4563ba16761SJacob Faibussowitsch PetscCallExternal(CVodeGetNonlinSolvStats, cvode->mem, &its, &nfails); 457bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of nonlinear solver iterations %ld\n", its)); 458bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of nonlinear convergence failure %ld\n", nfails)); 459918687b8SPatrick Sanan if (!cvode->use_dense) { 4603ba16761SJacob Faibussowitsch PetscCallExternal(CVSpilsGetNumLinIters, cvode->mem, &its); /* its = no. of calls to TSPrecond_Sundials() */ 461bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of linear iterations %ld\n", its)); 4623ba16761SJacob Faibussowitsch PetscCallExternal(CVSpilsGetNumConvFails, cvode->mem, &itmp); 463bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of linear convergence failures %ld\n", itmp)); 46442528757SHong Zhang 4659566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc)); 4669566063dSJacob Faibussowitsch PetscCall(PCView(pc, viewer)); 4673ba16761SJacob Faibussowitsch PetscCallExternal(CVSpilsGetNumPrecEvals, cvode->mem, &itmp); 468bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of preconditioner evaluations %ld\n", itmp)); 4693ba16761SJacob Faibussowitsch PetscCallExternal(CVSpilsGetNumPrecSolves, cvode->mem, &itmp); 470bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of preconditioner solves %ld\n", itmp)); 471918687b8SPatrick Sanan } 4723ba16761SJacob Faibussowitsch PetscCallExternal(CVSpilsGetNumJtimesEvals, cvode->mem, &itmp); 473bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of Jacobian-vector product evaluations %ld\n", itmp)); 4743ba16761SJacob Faibussowitsch PetscCallExternal(CVSpilsGetNumRhsEvals, cvode->mem, &itmp); 475bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of rhs calls for finite diff. Jacobian-vector evals %ld\n", itmp)); 4767cb94ee6SHong Zhang } else if (isstring) { 477bcf0153eSBarry Smith PetscCall(PetscViewerStringSPrintf(viewer, "SUNDIALS type %s", type)); 4787cb94ee6SHong Zhang } 4793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4807cb94ee6SHong Zhang } 4817cb94ee6SHong Zhang 4827cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 48366976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsSetType_Sundials(TS ts, TSSundialsLmmType type) 484d71ae5a4SJacob Faibussowitsch { 4857cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 4867cb94ee6SHong Zhang 4877cb94ee6SHong Zhang PetscFunctionBegin; 4887cb94ee6SHong Zhang cvode->cvode_type = type; 4893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4907cb94ee6SHong Zhang } 4917cb94ee6SHong Zhang 49266976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts, PetscInt maxl) 493d71ae5a4SJacob Faibussowitsch { 4947cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 4957cb94ee6SHong Zhang 4967cb94ee6SHong Zhang PetscFunctionBegin; 497f61b2b6aSHong Zhang cvode->maxl = maxl; 4983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4997cb94ee6SHong Zhang } 5007cb94ee6SHong Zhang 50166976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts, PetscReal tol) 502d71ae5a4SJacob Faibussowitsch { 5037cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 5047cb94ee6SHong Zhang 5057cb94ee6SHong Zhang PetscFunctionBegin; 5067cb94ee6SHong Zhang cvode->linear_tol = tol; 5073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5087cb94ee6SHong Zhang } 5097cb94ee6SHong Zhang 51066976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts, TSSundialsGramSchmidtType type) 511d71ae5a4SJacob Faibussowitsch { 5127cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 5137cb94ee6SHong Zhang 5147cb94ee6SHong Zhang PetscFunctionBegin; 5157cb94ee6SHong Zhang cvode->gtype = type; 5163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5177cb94ee6SHong Zhang } 5187cb94ee6SHong Zhang 51966976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts, PetscReal aabs, PetscReal rel) 520d71ae5a4SJacob Faibussowitsch { 5217cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 5227cb94ee6SHong Zhang 5237cb94ee6SHong Zhang PetscFunctionBegin; 5247cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 5257cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 5263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5277cb94ee6SHong Zhang } 5287cb94ee6SHong Zhang 52966976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts, PetscReal mindt) 530d71ae5a4SJacob Faibussowitsch { 531f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials *)ts->data; 532f1cd61daSJed Brown 533f1cd61daSJed Brown PetscFunctionBegin; 534f1cd61daSJed Brown cvode->mindt = mindt; 5353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 536f1cd61daSJed Brown } 537f1cd61daSJed Brown 53866976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts, PetscReal maxdt) 539d71ae5a4SJacob Faibussowitsch { 540f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials *)ts->data; 541f1cd61daSJed Brown 542f1cd61daSJed Brown PetscFunctionBegin; 543f1cd61daSJed Brown cvode->maxdt = maxdt; 5443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 545f1cd61daSJed Brown } 546918687b8SPatrick Sanan 54766976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsSetUseDense_Sundials(TS ts, PetscBool use_dense) 548d71ae5a4SJacob Faibussowitsch { 549918687b8SPatrick Sanan TS_Sundials *cvode = (TS_Sundials *)ts->data; 550918687b8SPatrick Sanan 551918687b8SPatrick Sanan PetscFunctionBegin; 552918687b8SPatrick Sanan cvode->use_dense = use_dense; 5533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 554918687b8SPatrick Sanan } 555918687b8SPatrick Sanan 55666976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsGetPC_Sundials(TS ts, PC *pc) 557d71ae5a4SJacob Faibussowitsch { 558dcbc6d53SSean Farley SNES snes; 559dcbc6d53SSean Farley KSP ksp; 5607cb94ee6SHong Zhang 5617cb94ee6SHong Zhang PetscFunctionBegin; 5629566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 5639566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 5649566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, pc)); 5653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5667cb94ee6SHong Zhang } 5677cb94ee6SHong Zhang 56866976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsGetIterations_Sundials(TS ts, int *nonlin, int *lin) 569d71ae5a4SJacob Faibussowitsch { 5707cb94ee6SHong Zhang PetscFunctionBegin; 5715ef26d82SJed Brown if (nonlin) *nonlin = ts->snes_its; 5725ef26d82SJed Brown if (lin) *lin = ts->ksp_its; 5733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5747cb94ee6SHong Zhang } 5757cb94ee6SHong Zhang 57666976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts, PetscBool s) 577d71ae5a4SJacob Faibussowitsch { 5782bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 5792bfc04deSHong Zhang 5802bfc04deSHong Zhang PetscFunctionBegin; 5812bfc04deSHong Zhang cvode->monitorstep = s; 5823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5832bfc04deSHong Zhang } 5847cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 5857cb94ee6SHong Zhang 5867cb94ee6SHong Zhang /*@C 587bcf0153eSBarry Smith TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by `TSSUNDIALS`. 5887cb94ee6SHong Zhang 5897cb94ee6SHong Zhang Not Collective 5907cb94ee6SHong Zhang 59197bb3fdcSJose E. Roman Input Parameter: 5927cb94ee6SHong Zhang . ts - the time-step context 5937cb94ee6SHong Zhang 5947cb94ee6SHong Zhang Output Parameters: 5957cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 5967cb94ee6SHong Zhang - lin - number of linear iterations 5977cb94ee6SHong Zhang 5987cb94ee6SHong Zhang Level: advanced 5997cb94ee6SHong Zhang 600bcf0153eSBarry Smith Note: 601bcf0153eSBarry Smith These return the number since the creation of the `TS` object 6027cb94ee6SHong Zhang 6031cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 604db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 60542747ad1SJacob Faibussowitsch `TSSundialsGetPC()`, `TSSetExactFinalTime()` 6067cb94ee6SHong Zhang @*/ 607d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsGetIterations(TS ts, int *nonlin, int *lin) 608d71ae5a4SJacob Faibussowitsch { 6097cb94ee6SHong Zhang PetscFunctionBegin; 610cac4c232SBarry Smith PetscUseMethod(ts, "TSSundialsGetIterations_C", (TS, int *, int *), (ts, nonlin, lin)); 6113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6127cb94ee6SHong Zhang } 6137cb94ee6SHong Zhang 6147cb94ee6SHong Zhang /*@ 615bcf0153eSBarry Smith TSSundialsSetType - Sets the method that `TSSUNDIALS` will use for integration. 6167cb94ee6SHong Zhang 617c3339decSBarry Smith Logically Collective 6187cb94ee6SHong Zhang 6198f7d6fe5SPatrick Sanan Input Parameters: 6207cb94ee6SHong Zhang + ts - the time-step context 621bcf0153eSBarry Smith - type - one of `SUNDIALS_ADAMS` or `SUNDIALS_BDF` 6227cb94ee6SHong Zhang 6237cb94ee6SHong Zhang Level: intermediate 6247cb94ee6SHong Zhang 6251cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsGetIterations()`, `TSSundialsSetMaxl()`, 62642747ad1SJacob Faibussowitsch `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, 62742747ad1SJacob Faibussowitsch `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, `TSSetExactFinalTime()` 6287cb94ee6SHong Zhang @*/ 629d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetType(TS ts, TSSundialsLmmType type) 630d71ae5a4SJacob Faibussowitsch { 6317cb94ee6SHong Zhang PetscFunctionBegin; 632cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetType_C", (TS, TSSundialsLmmType), (ts, type)); 6333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6347cb94ee6SHong Zhang } 6357cb94ee6SHong Zhang 6367cb94ee6SHong Zhang /*@ 637bcf0153eSBarry Smith TSSundialsSetMaxord - Sets the maximum order for BDF/Adams method used by `TSSUNDIALS`. 638c4e80e11SFlorian 639c3339decSBarry Smith Logically Collective 640c4e80e11SFlorian 6418f7d6fe5SPatrick Sanan Input Parameters: 642c4e80e11SFlorian + ts - the time-step context 643c4e80e11SFlorian - maxord - maximum order of BDF / Adams method 644c4e80e11SFlorian 645c4e80e11SFlorian Level: advanced 646c4e80e11SFlorian 6471cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, 64842747ad1SJacob Faibussowitsch `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, 64942747ad1SJacob Faibussowitsch `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, `TSSetExactFinalTime()` 650c4e80e11SFlorian @*/ 651d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetMaxord(TS ts, PetscInt maxord) 652d71ae5a4SJacob Faibussowitsch { 653c4e80e11SFlorian PetscFunctionBegin; 654c4e80e11SFlorian PetscValidLogicalCollectiveInt(ts, maxord, 2); 655cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetMaxOrd_C", (TS, PetscInt), (ts, maxord)); 6563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 657c4e80e11SFlorian } 658c4e80e11SFlorian 659c4e80e11SFlorian /*@ 660f61b2b6aSHong Zhang TSSundialsSetMaxl - Sets the dimension of the Krylov space used by 661bcf0153eSBarry Smith GMRES in the linear solver in `TSSUNDIALS`. `TSSUNDIALS` DOES NOT use restarted GMRES so 662f61b2b6aSHong Zhang this is the maximum number of GMRES steps that will be used. 6637cb94ee6SHong Zhang 664c3339decSBarry Smith Logically Collective 6657cb94ee6SHong Zhang 6668f7d6fe5SPatrick Sanan Input Parameters: 6677cb94ee6SHong Zhang + ts - the time-step context 668f61b2b6aSHong Zhang - maxl - number of direction vectors (the dimension of Krylov subspace). 6697cb94ee6SHong Zhang 6707cb94ee6SHong Zhang Level: advanced 6717cb94ee6SHong Zhang 6721cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, 673db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 674b43aa488SJacob Faibussowitsch `TSSundialsGetPC()`, `TSSetExactFinalTime()` 6757cb94ee6SHong Zhang @*/ 676d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetMaxl(TS ts, PetscInt maxl) 677d71ae5a4SJacob Faibussowitsch { 6787cb94ee6SHong Zhang PetscFunctionBegin; 679f61b2b6aSHong Zhang PetscValidLogicalCollectiveInt(ts, maxl, 2); 680cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetMaxl_C", (TS, PetscInt), (ts, maxl)); 6813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6827cb94ee6SHong Zhang } 6837cb94ee6SHong Zhang 6847cb94ee6SHong Zhang /*@ 6857cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 686bcf0153eSBarry Smith system by `TSSUNDIALS`. 6877cb94ee6SHong Zhang 688c3339decSBarry Smith Logically Collective 6897cb94ee6SHong Zhang 6908f7d6fe5SPatrick Sanan Input Parameters: 6917cb94ee6SHong Zhang + ts - the time-step context 6927cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 6937cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 6947cb94ee6SHong Zhang 6957cb94ee6SHong Zhang Level: advanced 6967cb94ee6SHong Zhang 6971cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 698db781477SPatrick Sanan `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 69942747ad1SJacob Faibussowitsch `TSSundialsGetPC()`, 700db781477SPatrick Sanan `TSSetExactFinalTime()` 7017cb94ee6SHong Zhang @*/ 702d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetLinearTolerance(TS ts, PetscReal tol) 703d71ae5a4SJacob Faibussowitsch { 7047cb94ee6SHong Zhang PetscFunctionBegin; 705c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts, tol, 2); 706cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetLinearTolerance_C", (TS, PetscReal), (ts, tol)); 7073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7087cb94ee6SHong Zhang } 7097cb94ee6SHong Zhang 7107cb94ee6SHong Zhang /*@ 7117cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 712bcf0153eSBarry Smith in GMRES method by `TSSUNDIALS` linear solver. 7137cb94ee6SHong Zhang 714c3339decSBarry Smith Logically Collective 7157cb94ee6SHong Zhang 7168f7d6fe5SPatrick Sanan Input Parameters: 7177cb94ee6SHong Zhang + ts - the time-step context 718bcf0153eSBarry Smith - type - either `SUNDIALS_MODIFIED_GS` or `SUNDIALS_CLASSICAL_GS` 7197cb94ee6SHong Zhang 7207cb94ee6SHong Zhang Level: advanced 7217cb94ee6SHong Zhang 7221cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 723db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, 724b43aa488SJacob Faibussowitsch `TSSundialsGetPC()`, 725db781477SPatrick Sanan `TSSetExactFinalTime()` 7267cb94ee6SHong Zhang @*/ 727d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetGramSchmidtType(TS ts, TSSundialsGramSchmidtType type) 728d71ae5a4SJacob Faibussowitsch { 7297cb94ee6SHong Zhang PetscFunctionBegin; 730cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetGramSchmidtType_C", (TS, TSSundialsGramSchmidtType), (ts, type)); 7313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7327cb94ee6SHong Zhang } 7337cb94ee6SHong Zhang 7347cb94ee6SHong Zhang /*@ 7357cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 736bcf0153eSBarry Smith `TSSUNDIALS` for error control. 7377cb94ee6SHong Zhang 738c3339decSBarry Smith Logically Collective 7397cb94ee6SHong Zhang 7408f7d6fe5SPatrick Sanan Input Parameters: 7417cb94ee6SHong Zhang + ts - the time-step context 7427cb94ee6SHong Zhang . aabs - the absolute tolerance 7437cb94ee6SHong Zhang - rel - the relative tolerance 7447cb94ee6SHong Zhang 745bcf0153eSBarry Smith See the CVODE/SUNDIALS users manual for exact details on these parameters. Essentially 7467cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 7477cb94ee6SHong Zhang 7487cb94ee6SHong Zhang Level: intermediate 7497cb94ee6SHong Zhang 7501cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetGMRESMaxl()`, 751db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, 75242747ad1SJacob Faibussowitsch `TSSundialsGetPC()`, 753db781477SPatrick Sanan `TSSetExactFinalTime()` 7547cb94ee6SHong Zhang @*/ 755d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetTolerance(TS ts, PetscReal aabs, PetscReal rel) 756d71ae5a4SJacob Faibussowitsch { 7577cb94ee6SHong Zhang PetscFunctionBegin; 758cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetTolerance_C", (TS, PetscReal, PetscReal), (ts, aabs, rel)); 7593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7607cb94ee6SHong Zhang } 7617cb94ee6SHong Zhang 7627cb94ee6SHong Zhang /*@ 763bcf0153eSBarry Smith TSSundialsGetPC - Extract the PC context from a time-step context for `TSSUNDIALS`. 7647cb94ee6SHong Zhang 7657cb94ee6SHong Zhang Input Parameter: 7667cb94ee6SHong Zhang . ts - the time-step context 7677cb94ee6SHong Zhang 7687cb94ee6SHong Zhang Output Parameter: 7697cb94ee6SHong Zhang . pc - the preconditioner context 7707cb94ee6SHong Zhang 7717cb94ee6SHong Zhang Level: advanced 7727cb94ee6SHong Zhang 7731cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 774b43aa488SJacob Faibussowitsch `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()` 7757cb94ee6SHong Zhang @*/ 776d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsGetPC(TS ts, PC *pc) 777d71ae5a4SJacob Faibussowitsch { 7787cb94ee6SHong Zhang PetscFunctionBegin; 779cac4c232SBarry Smith PetscUseMethod(ts, "TSSundialsGetPC_C", (TS, PC *), (ts, pc)); 7803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7817cb94ee6SHong Zhang } 7827cb94ee6SHong Zhang 783f1cd61daSJed Brown /*@ 784f1cd61daSJed Brown TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 785f1cd61daSJed Brown 786d8d19677SJose E. Roman Input Parameters: 787f1cd61daSJed Brown + ts - the time-step context 788f1cd61daSJed Brown - mindt - lowest time step if positive, negative to deactivate 789f1cd61daSJed Brown 790fc6b9e64SJed Brown Note: 791bcf0153eSBarry Smith `TSSUNDIALS` will error if it is not possible to keep the estimated truncation error below 792bcf0153eSBarry Smith the tolerance set with `TSSundialsSetTolerance()` without going below this step size. 793fc6b9e64SJed Brown 794f1cd61daSJed Brown Level: beginner 795f1cd61daSJed Brown 7961cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsSetType()`, `TSSundialsSetTolerance()`, 797f1cd61daSJed Brown @*/ 798d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetMinTimeStep(TS ts, PetscReal mindt) 799d71ae5a4SJacob Faibussowitsch { 800f1cd61daSJed Brown PetscFunctionBegin; 801cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetMinTimeStep_C", (TS, PetscReal), (ts, mindt)); 8023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 803f1cd61daSJed Brown } 804f1cd61daSJed Brown 805f1cd61daSJed Brown /*@ 806f1cd61daSJed Brown TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 807f1cd61daSJed Brown 808d8d19677SJose E. Roman Input Parameters: 809f1cd61daSJed Brown + ts - the time-step context 810f1cd61daSJed Brown - maxdt - lowest time step if positive, negative to deactivate 811f1cd61daSJed Brown 812f1cd61daSJed Brown Level: beginner 813f1cd61daSJed Brown 8141cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsSetType()`, `TSSundialsSetTolerance()`, 815f1cd61daSJed Brown @*/ 816d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetMaxTimeStep(TS ts, PetscReal maxdt) 817d71ae5a4SJacob Faibussowitsch { 818f1cd61daSJed Brown PetscFunctionBegin; 819cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt)); 8203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 821f1cd61daSJed Brown } 822f1cd61daSJed Brown 8232bfc04deSHong Zhang /*@ 824bcf0153eSBarry Smith TSSundialsMonitorInternalSteps - Monitor `TSSUNDIALS` internal steps (Defaults to false). 8252bfc04deSHong Zhang 826d8d19677SJose E. Roman Input Parameters: 8272bfc04deSHong Zhang + ts - the time-step context 828bcf0153eSBarry Smith - ft - `PETSC_TRUE` if monitor, else `PETSC_FALSE` 8292bfc04deSHong Zhang 8302bfc04deSHong Zhang Level: beginner 8312bfc04deSHong Zhang 8321cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 833bcf0153eSBarry Smith `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 834b43aa488SJacob Faibussowitsch `TSSundialsGetPC()` 8352bfc04deSHong Zhang @*/ 836d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsMonitorInternalSteps(TS ts, PetscBool ft) 837d71ae5a4SJacob Faibussowitsch { 8382bfc04deSHong Zhang PetscFunctionBegin; 839cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsMonitorInternalSteps_C", (TS, PetscBool), (ts, ft)); 8403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8412bfc04deSHong Zhang } 842918687b8SPatrick Sanan 843918687b8SPatrick Sanan /*@ 844bcf0153eSBarry Smith TSSundialsSetUseDense - Set a flag to use a dense linear solver in `TSSUNDIALS` (serial only) 845918687b8SPatrick Sanan 846918687b8SPatrick Sanan Logically Collective 847918687b8SPatrick Sanan 848918687b8SPatrick Sanan Input Parameters: 849918687b8SPatrick Sanan + ts - the time-step context 850bcf0153eSBarry Smith - use_dense - `PETSC_TRUE` to use the dense solver 851918687b8SPatrick Sanan 852918687b8SPatrick Sanan Level: advanced 853918687b8SPatrick Sanan 8541cc06b55SBarry Smith .seealso: [](ch_ts), `TSSUNDIALS` 855918687b8SPatrick Sanan @*/ 856d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetUseDense(TS ts, PetscBool use_dense) 857d71ae5a4SJacob Faibussowitsch { 858918687b8SPatrick Sanan PetscFunctionBegin; 859918687b8SPatrick Sanan PetscValidLogicalCollectiveInt(ts, use_dense, 2); 860cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetUseDense_C", (TS, PetscBool), (ts, use_dense)); 8613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 862918687b8SPatrick Sanan } 863918687b8SPatrick Sanan 8647cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 8657cb94ee6SHong Zhang /*MC 866bcf0153eSBarry Smith TSSUNDIALS - ODE solver using a very old version of the LLNL CVODE/SUNDIALS package, version 2.5 (now called SUNDIALS). Requires ./configure --download-sundials 8677cb94ee6SHong Zhang 868bcf0153eSBarry Smith Options Database Keys: 869897c9f78SHong Zhang + -ts_sundials_type <bdf,adams> - 8707cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 8717cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 8727cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 873897c9f78SHong Zhang . -ts_sundials_linear_tolerance <tol> - 874f61b2b6aSHong Zhang . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace 875918687b8SPatrick Sanan . -ts_sundials_monitor_steps - Monitor SUNDIALS internal steps 876918687b8SPatrick Sanan - -ts_sundials_use_dense - Use a dense linear solver within CVODE (serial only) 87716016685SHong Zhang 8787cb94ee6SHong Zhang Level: beginner 8797cb94ee6SHong Zhang 880bcf0153eSBarry Smith Note: 881bcf0153eSBarry Smith This uses its own nonlinear solver and Krylov method so PETSc `SNES` and `KSP` options do not apply, 882bcf0153eSBarry Smith only PETSc `PC` options. 8837cb94ee6SHong Zhang 8841cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, `TSSundialsSetLinearTolerance()`, `TSType`, 885bcf0153eSBarry Smith `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, `TSSundialsGetIterations()`, `TSSetExactFinalTime()` 8867cb94ee6SHong Zhang M*/ 887d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts) 888d71ae5a4SJacob Faibussowitsch { 8897cb94ee6SHong Zhang TS_Sundials *cvode; 89042528757SHong Zhang PC pc; 8917cb94ee6SHong Zhang 8927cb94ee6SHong Zhang PetscFunctionBegin; 893277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Sundials; 89428adb3f7SHong Zhang ts->ops->destroy = TSDestroy_Sundials; 89528adb3f7SHong Zhang ts->ops->view = TSView_Sundials; 896214bc6a2SJed Brown ts->ops->setup = TSSetUp_Sundials; 897b4eba00bSSean Farley ts->ops->step = TSStep_Sundials; 898b4eba00bSSean Farley ts->ops->interpolate = TSInterpolate_Sundials; 899214bc6a2SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Sundials; 9002ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 9017cb94ee6SHong Zhang 9024dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&cvode)); 903bbd56ea5SKarl Rupp 904efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 905efd4aadfSBarry Smith 9067cb94ee6SHong Zhang ts->data = (void *)cvode; 9076fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF; 9086fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS; 909f61b2b6aSHong Zhang cvode->maxl = 5; 910c4e80e11SFlorian cvode->maxord = PETSC_DEFAULT; 9117cb94ee6SHong Zhang cvode->linear_tol = .05; 912b4eba00bSSean Farley cvode->monitorstep = PETSC_TRUE; 913918687b8SPatrick Sanan cvode->use_dense = PETSC_FALSE; 9147cb94ee6SHong Zhang 915*f4f49eeaSPierre Jolivet PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)ts), &cvode->comm_sundials)); 916f1cd61daSJed Brown 917f1cd61daSJed Brown cvode->mindt = -1.; 918f1cd61daSJed Brown cvode->maxdt = -1.; 919f1cd61daSJed Brown 920bcf0153eSBarry Smith /* set tolerance for SUNDIALS */ 9217cb94ee6SHong Zhang cvode->reltol = 1e-6; 9222c823083SHong Zhang cvode->abstol = 1e-6; 9237cb94ee6SHong Zhang 92442528757SHong Zhang /* set PCNONE as default pctype */ 9259566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC_Sundials(ts, &pc)); 9269566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 92742528757SHong Zhang 9289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetType_C", TSSundialsSetType_Sundials)); 9299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxl_C", TSSundialsSetMaxl_Sundials)); 9309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetLinearTolerance_C", TSSundialsSetLinearTolerance_Sundials)); 9319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetGramSchmidtType_C", TSSundialsSetGramSchmidtType_Sundials)); 9329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetTolerance_C", TSSundialsSetTolerance_Sundials)); 9339566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMinTimeStep_C", TSSundialsSetMinTimeStep_Sundials)); 9349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxTimeStep_C", TSSundialsSetMaxTimeStep_Sundials)); 9359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetPC_C", TSSundialsGetPC_Sundials)); 9369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetIterations_C", TSSundialsGetIterations_Sundials)); 9379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsMonitorInternalSteps_C", TSSundialsMonitorInternalSteps_Sundials)); 9389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetUseDense_C", TSSundialsSetUseDense_Sundials)); 9393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9407cb94ee6SHong Zhang } 941