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 */ 149371c9d4SSatish Balay 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) { 157cb94ee6SHong Zhang TS ts = (TS)P_data; 167cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 17dcbc6d53SSean Farley PC pc; 180679a0aeSJed Brown Mat J, P; 190679a0aeSJed Brown Vec yy = cvode->w1, yydot = cvode->ydot; 200679a0aeSJed Brown PetscReal gm = (PetscReal)_gamma; 217cb94ee6SHong Zhang PetscScalar *y_data; 227cb94ee6SHong Zhang 237cb94ee6SHong Zhang PetscFunctionBegin; 249566063dSJacob Faibussowitsch PetscCall(TSGetIJacobian(ts, &J, &P, NULL, NULL)); 257cb94ee6SHong Zhang y_data = (PetscScalar *)N_VGetArrayPointer(y); 269566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(yy, y_data)); 279566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(yydot)); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */ 280679a0aeSJed Brown /* compute the shifted Jacobian (1/gm)*I + Jrest */ 299566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, ts->ptime, yy, yydot, 1 / gm, J, P, PETSC_FALSE)); 309566063dSJacob Faibussowitsch PetscCall(VecResetArray(yy)); 319566063dSJacob Faibussowitsch PetscCall(MatScale(P, gm)); /* turn into I-gm*Jrest, J is not used by Sundials */ 327cb94ee6SHong Zhang *jcurPtr = TRUE; 339566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc)); 349566063dSJacob Faibussowitsch PetscCall(PCSetOperators(pc, J, P)); 357cb94ee6SHong Zhang PetscFunctionReturn(0); 367cb94ee6SHong Zhang } 377cb94ee6SHong Zhang 387cb94ee6SHong Zhang /* 397cb94ee6SHong Zhang TSPSolve_Sundials - routine that we provide to Sundials that applies the preconditioner. 407cb94ee6SHong Zhang */ 419371c9d4SSatish Balay 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) { 427cb94ee6SHong Zhang TS ts = (TS)P_data; 437cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 44dcbc6d53SSean Farley PC pc; 457cb94ee6SHong Zhang Vec rr = cvode->w1, zz = cvode->w2; 467cb94ee6SHong Zhang PetscScalar *r_data, *z_data; 477cb94ee6SHong Zhang 487cb94ee6SHong Zhang PetscFunctionBegin; 497cb94ee6SHong Zhang /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/ 507cb94ee6SHong Zhang r_data = (PetscScalar *)N_VGetArrayPointer(r); 517cb94ee6SHong Zhang z_data = (PetscScalar *)N_VGetArrayPointer(z); 529566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(rr, r_data)); 539566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(zz, z_data)); 544ac7836bSHong Zhang 557cb94ee6SHong Zhang /* Solve the Px=r and put the result in zz */ 569566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc)); 579566063dSJacob Faibussowitsch PetscCall(PCApply(pc, rr, zz)); 589566063dSJacob Faibussowitsch PetscCall(VecResetArray(rr)); 599566063dSJacob Faibussowitsch PetscCall(VecResetArray(zz)); 607cb94ee6SHong Zhang PetscFunctionReturn(0); 617cb94ee6SHong Zhang } 627cb94ee6SHong Zhang 637cb94ee6SHong Zhang /* 647cb94ee6SHong Zhang TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side. 657cb94ee6SHong Zhang */ 669371c9d4SSatish Balay int TSFunction_Sundials(realtype t, N_Vector y, N_Vector ydot, void *ctx) { 677cb94ee6SHong Zhang TS ts = (TS)ctx; 685fcef5e4SJed Brown DM dm; 69942e3340SBarry Smith DMTS tsdm; 705fcef5e4SJed Brown TSIFunction ifunction; 71ce94432eSBarry Smith MPI_Comm comm; 727cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 730679a0aeSJed Brown Vec yy = cvode->w1, yyd = cvode->w2, yydot = cvode->ydot; 747cb94ee6SHong Zhang PetscScalar *y_data, *ydot_data; 757cb94ee6SHong Zhang 767cb94ee6SHong Zhang PetscFunctionBegin; 779566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 785ff03cf7SDinesh Kaushik /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/ 797cb94ee6SHong Zhang y_data = (PetscScalar *)N_VGetArrayPointer(y); 807cb94ee6SHong Zhang ydot_data = (PetscScalar *)N_VGetArrayPointer(ydot); 819566063dSJacob Faibussowitsch PetscCallAbort(comm, VecPlaceArray(yy, y_data)); 829566063dSJacob Faibussowitsch PetscCallAbort(comm, VecPlaceArray(yyd, ydot_data)); 834ac7836bSHong Zhang 845fcef5e4SJed Brown /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */ 859566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 869566063dSJacob Faibussowitsch PetscCall(DMGetDMTS(dm, &tsdm)); 879566063dSJacob Faibussowitsch PetscCall(DMTSGetIFunction(dm, &ifunction, NULL)); 885fcef5e4SJed Brown if (!ifunction) { 899566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, t, yy, yyd)); 90bc0cc02bSJed Brown } else { /* If rhsfunction is also set, this computes both parts and shifts them to the right */ 919566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(yydot)); 929566063dSJacob Faibussowitsch PetscCallAbort(comm, TSComputeIFunction(ts, t, yy, yydot, yyd, PETSC_FALSE)); 939566063dSJacob Faibussowitsch PetscCall(VecScale(yyd, -1.)); 94bc0cc02bSJed Brown } 959566063dSJacob Faibussowitsch PetscCallAbort(comm, VecResetArray(yy)); 969566063dSJacob Faibussowitsch PetscCallAbort(comm, VecResetArray(yyd)); 974ac7836bSHong Zhang PetscFunctionReturn(0); 987cb94ee6SHong Zhang } 997cb94ee6SHong Zhang 1007cb94ee6SHong Zhang /* 101b4eba00bSSean Farley TSStep_Sundials - Calls Sundials to integrate the ODE. 1027cb94ee6SHong Zhang */ 1039371c9d4SSatish Balay PetscErrorCode TSStep_Sundials(TS ts) { 1047cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 105b4eba00bSSean Farley PetscInt flag; 106be5899b3SLisandro Dalcin long int nits, lits, nsteps; 1077cb94ee6SHong Zhang realtype t, tout; 1087cb94ee6SHong Zhang PetscScalar *y_data; 1097cb94ee6SHong Zhang void *mem; 1107cb94ee6SHong Zhang 1117cb94ee6SHong Zhang PetscFunctionBegin; 11216016685SHong Zhang mem = cvode->mem; 1137cb94ee6SHong Zhang tout = ts->max_time; 1149566063dSJacob Faibussowitsch PetscCall(VecGetArray(ts->vec_sol, &y_data)); 1157cb94ee6SHong Zhang N_VSetArrayPointer((realtype *)y_data, cvode->y); 1169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ts->vec_sol, NULL)); 117186e87acSLisandro Dalcin 1189687d888SLisandro Dalcin /* We would like to TSPreStage() and TSPostStage() 1199687d888SLisandro Dalcin * before each stage solve but CVode does not appear to support this. */ 1209371c9d4SSatish Balay if (cvode->monitorstep) flag = CVode(mem, tout, cvode->y, &t, CV_ONE_STEP); 1219371c9d4SSatish Balay else flag = CVode(mem, tout, cvode->y, &t, CV_NORMAL); 1229f94935aSHong Zhang 1239f94935aSHong Zhang if (flag) { /* display error message */ 1249f94935aSHong Zhang switch (flag) { 1259371c9d4SSatish Balay case CV_ILL_INPUT: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_ILL_INPUT"); break; 1269371c9d4SSatish Balay case CV_TOO_CLOSE: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_TOO_CLOSE"); break; 1273c7fefeeSJed Brown case CV_TOO_MUCH_WORK: { 1289f94935aSHong Zhang PetscReal tcur; 1299566063dSJacob Faibussowitsch PetscCall(CVodeGetNumSteps(mem, &nsteps)); 1309566063dSJacob Faibussowitsch PetscCall(CVodeGetCurrentTime(mem, &tcur)); 13163a3b9bcSJacob 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); 1323c7fefeeSJed Brown } break; 1339371c9d4SSatish Balay case CV_TOO_MUCH_ACC: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_TOO_MUCH_ACC"); break; 1349371c9d4SSatish Balay case CV_ERR_FAILURE: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_ERR_FAILURE"); break; 1359371c9d4SSatish Balay case CV_CONV_FAILURE: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_CONV_FAILURE"); break; 1369371c9d4SSatish Balay case CV_LINIT_FAIL: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LINIT_FAIL"); break; 1379371c9d4SSatish Balay case CV_LSETUP_FAIL: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LSETUP_FAIL"); break; 1389371c9d4SSatish Balay case CV_LSOLVE_FAIL: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LSOLVE_FAIL"); break; 1399371c9d4SSatish Balay case CV_RHSFUNC_FAIL: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_RHSFUNC_FAIL"); break; 1409371c9d4SSatish Balay case CV_FIRST_RHSFUNC_ERR: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_FIRST_RHSFUNC_ERR"); break; 1419371c9d4SSatish Balay case CV_REPTD_RHSFUNC_ERR: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_REPTD_RHSFUNC_ERR"); break; 1429371c9d4SSatish Balay case CV_UNREC_RHSFUNC_ERR: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_UNREC_RHSFUNC_ERR"); break; 1439371c9d4SSatish Balay case CV_RTFUNC_FAIL: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_RTFUNC_FAIL"); break; 1449371c9d4SSatish Balay default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, flag %d", flag); 1459f94935aSHong Zhang } 1469f94935aSHong Zhang } 1479f94935aSHong Zhang 148be5899b3SLisandro Dalcin /* log inner nonlinear and linear iterations */ 1499566063dSJacob Faibussowitsch PetscCall(CVodeGetNumNonlinSolvIters(mem, &nits)); 1509566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumLinIters(mem, &lits)); 1519371c9d4SSatish Balay ts->snes_its += nits; 1529371c9d4SSatish Balay ts->ksp_its = lits; 153be5899b3SLisandro Dalcin 1547cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */ 1559566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(cvode->w1, y_data)); 1569566063dSJacob Faibussowitsch PetscCall(VecCopy(cvode->w1, cvode->update)); 1579566063dSJacob Faibussowitsch PetscCall(VecResetArray(cvode->w1)); 1589566063dSJacob Faibussowitsch PetscCall(VecCopy(cvode->update, ts->vec_sol)); 159186e87acSLisandro Dalcin 160186e87acSLisandro Dalcin ts->time_step = t - ts->ptime; 161186e87acSLisandro Dalcin ts->ptime = t; 162186e87acSLisandro Dalcin 1639566063dSJacob Faibussowitsch PetscCall(CVodeGetNumSteps(mem, &nsteps)); 1642808aa04SLisandro Dalcin if (!cvode->monitorstep) ts->steps += nsteps - 1; /* TSStep() increments the step counter by one */ 165b4eba00bSSean Farley PetscFunctionReturn(0); 166b4eba00bSSean Farley } 167b4eba00bSSean Farley 1689371c9d4SSatish Balay static PetscErrorCode TSInterpolate_Sundials(TS ts, PetscReal t, Vec X) { 169b4eba00bSSean Farley TS_Sundials *cvode = (TS_Sundials *)ts->data; 170b4eba00bSSean Farley N_Vector y; 171b4eba00bSSean Farley PetscScalar *x_data; 172b4eba00bSSean Farley PetscInt glosize, locsize; 173b4eba00bSSean Farley 174b4eba00bSSean Farley PetscFunctionBegin; 175b4eba00bSSean Farley /* get the vector size */ 1769566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &glosize)); 1779566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &locsize)); 1789566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x_data)); 179b4eba00bSSean Farley 1804d78c9e0SClaas Abert /* Initialize N_Vec y with x_data */ 181918687b8SPatrick Sanan if (cvode->use_dense) { 182918687b8SPatrick Sanan PetscMPIInt size; 183918687b8SPatrick Sanan 1849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1853c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "TSSUNDIALS only supports a dense solve in the serial case"); 186918687b8SPatrick Sanan y = N_VMake_Serial(locsize, (realtype *)x_data); 187918687b8SPatrick Sanan } else { 1884d78c9e0SClaas Abert y = N_VMake_Parallel(cvode->comm_sundials, locsize, glosize, (realtype *)x_data); 189918687b8SPatrick Sanan } 190918687b8SPatrick Sanan 1913c633725SBarry Smith PetscCheck(y, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Interpolated y is not allocated"); 192b4eba00bSSean Farley 1939566063dSJacob Faibussowitsch PetscCall(CVodeGetDky(cvode->mem, t, 0, y)); 1949566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x_data)); 1957cb94ee6SHong Zhang PetscFunctionReturn(0); 1967cb94ee6SHong Zhang } 1977cb94ee6SHong Zhang 1989371c9d4SSatish Balay PetscErrorCode TSReset_Sundials(TS ts) { 199277b19d0SLisandro Dalcin TS_Sundials *cvode = (TS_Sundials *)ts->data; 200277b19d0SLisandro Dalcin 201277b19d0SLisandro Dalcin PetscFunctionBegin; 2029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->update)); 2039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->ydot)); 2049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->w1)); 2059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->w2)); 206bbd56ea5SKarl Rupp if (cvode->mem) CVodeFree(&cvode->mem); 207277b19d0SLisandro Dalcin PetscFunctionReturn(0); 208277b19d0SLisandro Dalcin } 209277b19d0SLisandro Dalcin 2109371c9d4SSatish Balay PetscErrorCode TSDestroy_Sundials(TS ts) { 2117cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 2127cb94ee6SHong Zhang 2137cb94ee6SHong Zhang PetscFunctionBegin; 2149566063dSJacob Faibussowitsch PetscCall(TSReset_Sundials(ts)); 2159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&(cvode->comm_sundials))); 2169566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 2179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetType_C", NULL)); 2189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxl_C", NULL)); 2199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetLinearTolerance_C", NULL)); 2209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetGramSchmidtType_C", NULL)); 2219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetTolerance_C", NULL)); 2229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMinTimeStep_C", NULL)); 2239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxTimeStep_C", NULL)); 2249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetPC_C", NULL)); 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetIterations_C", NULL)); 2269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsMonitorInternalSteps_C", NULL)); 2272e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetUseDense_C", NULL)); 2287cb94ee6SHong Zhang PetscFunctionReturn(0); 2297cb94ee6SHong Zhang } 2307cb94ee6SHong Zhang 2319371c9d4SSatish Balay PetscErrorCode TSSetUp_Sundials(TS ts) { 2327cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 23316016685SHong Zhang PetscInt glosize, locsize, i, flag; 2347cb94ee6SHong Zhang PetscScalar *y_data, *parray; 23516016685SHong Zhang void *mem; 236dcbc6d53SSean Farley PC pc; 23719fd82e9SBarry Smith PCType pctype; 238ace3abfcSBarry Smith PetscBool pcnone; 2397cb94ee6SHong Zhang 2407cb94ee6SHong Zhang PetscFunctionBegin; 24108401ef6SPierre 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"); 242236c45dcSShri Abhyankar 2437cb94ee6SHong Zhang /* get the vector size */ 2449566063dSJacob Faibussowitsch PetscCall(VecGetSize(ts->vec_sol, &glosize)); 2459566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ts->vec_sol, &locsize)); 2467cb94ee6SHong Zhang 2477cb94ee6SHong Zhang /* allocate the memory for N_Vec y */ 248918687b8SPatrick Sanan if (cvode->use_dense) { 249918687b8SPatrick Sanan PetscMPIInt size; 250918687b8SPatrick Sanan 2519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 2523c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "TSSUNDIALS only supports a dense solve in the serial case"); 253918687b8SPatrick Sanan cvode->y = N_VNew_Serial(locsize); 254918687b8SPatrick Sanan } else { 255a07356b8SSatish Balay cvode->y = N_VNew_Parallel(cvode->comm_sundials, locsize, glosize); 256918687b8SPatrick Sanan } 2573c633725SBarry Smith PetscCheck(cvode->y, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "cvode->y is not allocated"); 2587cb94ee6SHong Zhang 25928adb3f7SHong Zhang /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 2609566063dSJacob Faibussowitsch PetscCall(VecGetArray(ts->vec_sol, &parray)); 2617cb94ee6SHong Zhang y_data = (PetscScalar *)N_VGetArrayPointer(cvode->y); 2627cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 2639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ts->vec_sol, NULL)); 26442528757SHong Zhang 2659566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &cvode->update)); 2669566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &cvode->ydot)); 2677cb94ee6SHong Zhang 2687cb94ee6SHong Zhang /* 2697cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are 2707cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided 2717cb94ee6SHong Zhang by Sundials and set using VecPlaceArray(). 2727cb94ee6SHong Zhang */ 2739566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts), 1, locsize, PETSC_DECIDE, NULL, &cvode->w1)); 2749566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts), 1, locsize, PETSC_DECIDE, NULL, &cvode->w2)); 27516016685SHong Zhang 27616016685SHong Zhang /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 27716016685SHong Zhang mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 2783c633725SBarry Smith PetscCheck(mem, PETSC_COMM_SELF, PETSC_ERR_MEM, "CVodeCreate() fails"); 27916016685SHong Zhang cvode->mem = mem; 28016016685SHong Zhang 28116016685SHong Zhang /* Set the pointer to user-defined data */ 28216016685SHong Zhang flag = CVodeSetUserData(mem, ts); 28328b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeSetUserData() fails"); 28416016685SHong Zhang 285fc6b9e64SJed Brown /* Sundials may choose to use a smaller initial step, but will never use a larger step. */ 286cdbf8f93SLisandro Dalcin flag = CVodeSetInitStep(mem, (realtype)ts->time_step); 28728b400f6SJacob Faibussowitsch PetscCheck(!flag, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetInitStep() failed"); 288f1cd61daSJed Brown if (cvode->mindt > 0) { 289f1cd61daSJed Brown flag = CVodeSetMinStep(mem, (realtype)cvode->mindt); 2909f94935aSHong Zhang if (flag) { 29108401ef6SPierre Jolivet PetscCheck(flag != CV_MEM_NULL, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMinStep() failed, cvode_mem pointer is NULL"); 292f7d195e4SLawrence 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"); 293f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMinStep() failed"); 2949f94935aSHong Zhang } 295f1cd61daSJed Brown } 296f1cd61daSJed Brown if (cvode->maxdt > 0) { 297f1cd61daSJed Brown flag = CVodeSetMaxStep(mem, (realtype)cvode->maxdt); 29828b400f6SJacob Faibussowitsch PetscCheck(!flag, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMaxStep() failed"); 299f1cd61daSJed Brown } 300f1cd61daSJed Brown 30116016685SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the 302a5b23f4aSJose E. Roman * user's right hand side function in u'=f(t,u), the initial time T0, and 30316016685SHong Zhang * the initial dependent variable vector cvode->y */ 30416016685SHong Zhang flag = CVodeInit(mem, TSFunction_Sundials, ts->ptime, cvode->y); 30528b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeInit() fails, flag %d", flag); 30616016685SHong Zhang 3079f94935aSHong Zhang /* specifies scalar relative and absolute tolerances */ 30816016685SHong Zhang flag = CVodeSStolerances(mem, cvode->reltol, cvode->abstol); 30928b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeSStolerances() fails, flag %d", flag); 31016016685SHong Zhang 311c4e80e11SFlorian /* Specify max order of BDF / ADAMS method */ 312c4e80e11SFlorian if (cvode->maxord != PETSC_DEFAULT) { 313c4e80e11SFlorian flag = CVodeSetMaxOrd(mem, cvode->maxord); 31428b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeSetMaxOrd() fails, flag %d", flag); 315c4e80e11SFlorian } 316c4e80e11SFlorian 3179f94935aSHong Zhang /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */ 3189f94935aSHong Zhang flag = CVodeSetMaxNumSteps(mem, ts->max_steps); 31928b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeSetMaxNumSteps() fails, flag %d", flag); 32016016685SHong Zhang 321918687b8SPatrick Sanan if (cvode->use_dense) { 322918687b8SPatrick Sanan /* call CVDense to use a dense linear solver. */ 323918687b8SPatrick Sanan PetscMPIInt size; 324918687b8SPatrick Sanan 3259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 3263c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "TSSUNDIALS only supports a dense solve in the serial case"); 327918687b8SPatrick Sanan flag = CVDense(mem, locsize); 32828b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVDense() fails, flag %d", flag); 329918687b8SPatrick Sanan } else { 33016016685SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 33116016685SHong Zhang /* setup the ode integrator with the given preconditioner */ 3329566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc)); 3339566063dSJacob Faibussowitsch PetscCall(PCGetType(pc, &pctype)); 3349566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCNONE, &pcnone)); 33516016685SHong Zhang if (pcnone) { 33616016685SHong Zhang flag = CVSpgmr(mem, PREC_NONE, 0); 33728b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVSpgmr() fails, flag %d", flag); 33816016685SHong Zhang } else { 339f61b2b6aSHong Zhang flag = CVSpgmr(mem, PREC_LEFT, cvode->maxl); 34028b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVSpgmr() fails, flag %d", flag); 34116016685SHong Zhang 34216016685SHong Zhang /* Set preconditioner and solve routines Precond and PSolve, 34316016685SHong Zhang and the pointer to the user-defined block data */ 34416016685SHong Zhang flag = CVSpilsSetPreconditioner(mem, TSPrecond_Sundials, TSPSolve_Sundials); 34528b400f6SJacob Faibussowitsch PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVSpilsSetPreconditioner() fails, flag %d", flag); 34616016685SHong Zhang } 347918687b8SPatrick Sanan } 3487cb94ee6SHong Zhang PetscFunctionReturn(0); 3497cb94ee6SHong Zhang } 3507cb94ee6SHong Zhang 3516fadb2cdSHong Zhang /* type of CVODE linear multistep method */ 352c793f718SLisandro Dalcin const char *const TSSundialsLmmTypes[] = {"", "ADAMS", "BDF", "TSSundialsLmmType", "SUNDIALS_", NULL}; 3536fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */ 354c793f718SLisandro Dalcin const char *const TSSundialsGramSchmidtTypes[] = {"", "MODIFIED", "CLASSICAL", "TSSundialsGramSchmidtType", "SUNDIALS_", NULL}; 355a04cf4d8SBarry Smith 3569371c9d4SSatish Balay PetscErrorCode TSSetFromOptions_Sundials(TS ts, PetscOptionItems *PetscOptionsObject) { 3577cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 3587cb94ee6SHong Zhang int indx; 359ace3abfcSBarry Smith PetscBool flag; 3607bda3a07SJed Brown PC pc; 3617cb94ee6SHong Zhang 3627cb94ee6SHong Zhang PetscFunctionBegin; 363d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SUNDIALS ODE solver options"); 3649566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_sundials_type", "Scheme", "TSSundialsSetType", TSSundialsLmmTypes, 3, TSSundialsLmmTypes[cvode->cvode_type], &indx, &flag)); 3651baa6e33SBarry Smith if (flag) PetscCall(TSSundialsSetType(ts, (TSSundialsLmmType)indx)); 3669566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_sundials_gramschmidt_type", "Type of orthogonalization", "TSSundialsSetGramSchmidtType", TSSundialsGramSchmidtTypes, 3, TSSundialsGramSchmidtTypes[cvode->gtype], &indx, &flag)); 3671baa6e33SBarry Smith if (flag) PetscCall(TSSundialsSetGramSchmidtType(ts, (TSSundialsGramSchmidtType)indx)); 3689566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_atol", "Absolute tolerance for convergence", "TSSundialsSetTolerance", cvode->abstol, &cvode->abstol, NULL)); 3699566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_rtol", "Relative tolerance for convergence", "TSSundialsSetTolerance", cvode->reltol, &cvode->reltol, NULL)); 3709566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_mindt", "Minimum step size", "TSSundialsSetMinTimeStep", cvode->mindt, &cvode->mindt, NULL)); 3719566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_maxdt", "Maximum step size", "TSSundialsSetMaxTimeStep", cvode->maxdt, &cvode->maxdt, NULL)); 3729566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_linear_tolerance", "Convergence tolerance for linear solve", "TSSundialsSetLinearTolerance", cvode->linear_tol, &cvode->linear_tol, NULL)); 3739566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_sundials_maxord", "Max Order for BDF/Adams method", "TSSundialsSetMaxOrd", cvode->maxord, &cvode->maxord, NULL)); 3749566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_sundials_maxl", "Max dimension of the Krylov subspace", "TSSundialsSetMaxl", cvode->maxl, &cvode->maxl, NULL)); 3759566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_sundials_monitor_steps", "Monitor SUNDIALS internal steps", "TSSundialsMonitorInternalSteps", cvode->monitorstep, &cvode->monitorstep, NULL)); 3769566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_sundials_use_dense", "Use dense internal solver in SUNDIALS (serial only)", "TSSundialsSetUseDense", cvode->use_dense, &cvode->use_dense, NULL)); 377d0609cedSBarry Smith PetscOptionsHeadEnd(); 3789566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc)); 3799566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions(pc)); 3807cb94ee6SHong Zhang PetscFunctionReturn(0); 3817cb94ee6SHong Zhang } 3827cb94ee6SHong Zhang 3839371c9d4SSatish Balay PetscErrorCode TSView_Sundials(TS ts, PetscViewer viewer) { 3847cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 3857cb94ee6SHong Zhang char *type; 3867cb94ee6SHong Zhang char atype[] = "Adams"; 3877cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 388ace3abfcSBarry Smith PetscBool iascii, isstring; 3892c823083SHong Zhang long int nsteps, its, nfevals, nlinsetups, nfails, itmp; 3902c823083SHong Zhang PetscInt qlast, qcur; 3915d47aee6SHong Zhang PetscReal hinused, hlast, hcur, tcur, tolsfac; 39242528757SHong Zhang PC pc; 3937cb94ee6SHong Zhang 3947cb94ee6SHong Zhang PetscFunctionBegin; 395bbd56ea5SKarl Rupp if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype; 396bbd56ea5SKarl Rupp else type = btype; 3977cb94ee6SHong Zhang 3989566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 3999566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 4007cb94ee6SHong Zhang if (iascii) { 4019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials integrator does not use SNES!\n")); 4029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials integrator type %s\n", type)); 40363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials maxord %" PetscInt_FMT "\n", cvode->maxord)); 4049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials abs tol %g rel tol %g\n", (double)cvode->abstol, (double)cvode->reltol)); 405918687b8SPatrick Sanan if (cvode->use_dense) { 4069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials integrator using a dense linear solve\n")); 407918687b8SPatrick Sanan } else { 4089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials linear solver tolerance factor %g\n", (double)cvode->linear_tol)); 40963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials max dimension of Krylov subspace %" PetscInt_FMT "\n", cvode->maxl)); 4107cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 4119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n")); 4127cb94ee6SHong Zhang } else { 4139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n")); 4147cb94ee6SHong Zhang } 415918687b8SPatrick Sanan } 4169566063dSJacob Faibussowitsch if (cvode->mindt > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials minimum time step %g\n", (double)cvode->mindt)); 4179566063dSJacob Faibussowitsch if (cvode->maxdt > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials maximum time step %g\n", (double)cvode->maxdt)); 4182c823083SHong Zhang 4195d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */ 4209566063dSJacob Faibussowitsch PetscCall(CVodeGetTolScaleFactor(cvode->mem, &tolsfac)); 4219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials suggested factor for tolerance scaling %g\n", tolsfac)); 4229566063dSJacob Faibussowitsch PetscCall(CVodeGetIntegratorStats(cvode->mem, &nsteps, &nfevals, &nlinsetups, &nfails, &qlast, &qcur, &hinused, &hlast, &hcur, &tcur)); 42363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials cumulative number of internal steps %ld\n", nsteps)); 42463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of calls to rhs function %ld\n", nfevals)); 42563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of calls to linear solver setup function %ld\n", nlinsetups)); 42663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of error test failures %ld\n", nfails)); 4272c823083SHong Zhang 4289566063dSJacob Faibussowitsch PetscCall(CVodeGetNonlinSolvStats(cvode->mem, &its, &nfails)); 42963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of nonlinear solver iterations %ld\n", its)); 43063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of nonlinear convergence failure %ld\n", nfails)); 431918687b8SPatrick Sanan if (!cvode->use_dense) { 4329566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumLinIters(cvode->mem, &its)); /* its = no. of calls to TSPrecond_Sundials() */ 43363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of linear iterations %ld\n", its)); 4349566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumConvFails(cvode->mem, &itmp)); 43563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of linear convergence failures %ld\n", itmp)); 43642528757SHong Zhang 4379566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc)); 4389566063dSJacob Faibussowitsch PetscCall(PCView(pc, viewer)); 4399566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumPrecEvals(cvode->mem, &itmp)); 44063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of preconditioner evaluations %ld\n", itmp)); 4419566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumPrecSolves(cvode->mem, &itmp)); 44263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of preconditioner solves %ld\n", itmp)); 443918687b8SPatrick Sanan } 4449566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumJtimesEvals(cvode->mem, &itmp)); 44563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of Jacobian-vector product evaluations %ld\n", itmp)); 4469566063dSJacob Faibussowitsch PetscCall(CVSpilsGetNumRhsEvals(cvode->mem, &itmp)); 44763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of rhs calls for finite diff. Jacobian-vector evals %ld\n", itmp)); 4487cb94ee6SHong Zhang } else if (isstring) { 4499566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer, "Sundials type %s", type)); 4507cb94ee6SHong Zhang } 4517cb94ee6SHong Zhang PetscFunctionReturn(0); 4527cb94ee6SHong Zhang } 4537cb94ee6SHong Zhang 4547cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 4559371c9d4SSatish Balay PetscErrorCode TSSundialsSetType_Sundials(TS ts, TSSundialsLmmType type) { 4567cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 4577cb94ee6SHong Zhang 4587cb94ee6SHong Zhang PetscFunctionBegin; 4597cb94ee6SHong Zhang cvode->cvode_type = type; 4607cb94ee6SHong Zhang PetscFunctionReturn(0); 4617cb94ee6SHong Zhang } 4627cb94ee6SHong Zhang 4639371c9d4SSatish Balay PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts, PetscInt maxl) { 4647cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 4657cb94ee6SHong Zhang 4667cb94ee6SHong Zhang PetscFunctionBegin; 467f61b2b6aSHong Zhang cvode->maxl = maxl; 4687cb94ee6SHong Zhang PetscFunctionReturn(0); 4697cb94ee6SHong Zhang } 4707cb94ee6SHong Zhang 4719371c9d4SSatish Balay PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts, PetscReal tol) { 4727cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 4737cb94ee6SHong Zhang 4747cb94ee6SHong Zhang PetscFunctionBegin; 4757cb94ee6SHong Zhang cvode->linear_tol = tol; 4767cb94ee6SHong Zhang PetscFunctionReturn(0); 4777cb94ee6SHong Zhang } 4787cb94ee6SHong Zhang 4799371c9d4SSatish Balay PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts, TSSundialsGramSchmidtType type) { 4807cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 4817cb94ee6SHong Zhang 4827cb94ee6SHong Zhang PetscFunctionBegin; 4837cb94ee6SHong Zhang cvode->gtype = type; 4847cb94ee6SHong Zhang PetscFunctionReturn(0); 4857cb94ee6SHong Zhang } 4867cb94ee6SHong Zhang 4879371c9d4SSatish Balay PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts, PetscReal aabs, PetscReal rel) { 4887cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 4897cb94ee6SHong Zhang 4907cb94ee6SHong Zhang PetscFunctionBegin; 4917cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 4927cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 4937cb94ee6SHong Zhang PetscFunctionReturn(0); 4947cb94ee6SHong Zhang } 4957cb94ee6SHong Zhang 4969371c9d4SSatish Balay PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts, PetscReal mindt) { 497f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials *)ts->data; 498f1cd61daSJed Brown 499f1cd61daSJed Brown PetscFunctionBegin; 500f1cd61daSJed Brown cvode->mindt = mindt; 501f1cd61daSJed Brown PetscFunctionReturn(0); 502f1cd61daSJed Brown } 503f1cd61daSJed Brown 5049371c9d4SSatish Balay PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts, PetscReal maxdt) { 505f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials *)ts->data; 506f1cd61daSJed Brown 507f1cd61daSJed Brown PetscFunctionBegin; 508f1cd61daSJed Brown cvode->maxdt = maxdt; 509f1cd61daSJed Brown PetscFunctionReturn(0); 510f1cd61daSJed Brown } 511918687b8SPatrick Sanan 5129371c9d4SSatish Balay PetscErrorCode TSSundialsSetUseDense_Sundials(TS ts, PetscBool use_dense) { 513918687b8SPatrick Sanan TS_Sundials *cvode = (TS_Sundials *)ts->data; 514918687b8SPatrick Sanan 515918687b8SPatrick Sanan PetscFunctionBegin; 516918687b8SPatrick Sanan cvode->use_dense = use_dense; 517918687b8SPatrick Sanan PetscFunctionReturn(0); 518918687b8SPatrick Sanan } 519918687b8SPatrick Sanan 5209371c9d4SSatish Balay PetscErrorCode TSSundialsGetPC_Sundials(TS ts, PC *pc) { 521dcbc6d53SSean Farley SNES snes; 522dcbc6d53SSean Farley KSP ksp; 5237cb94ee6SHong Zhang 5247cb94ee6SHong Zhang PetscFunctionBegin; 5259566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 5269566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 5279566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, pc)); 5287cb94ee6SHong Zhang PetscFunctionReturn(0); 5297cb94ee6SHong Zhang } 5307cb94ee6SHong Zhang 5319371c9d4SSatish Balay PetscErrorCode TSSundialsGetIterations_Sundials(TS ts, int *nonlin, int *lin) { 5327cb94ee6SHong Zhang PetscFunctionBegin; 5335ef26d82SJed Brown if (nonlin) *nonlin = ts->snes_its; 5345ef26d82SJed Brown if (lin) *lin = ts->ksp_its; 5357cb94ee6SHong Zhang PetscFunctionReturn(0); 5367cb94ee6SHong Zhang } 5377cb94ee6SHong Zhang 5389371c9d4SSatish Balay PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts, PetscBool s) { 5392bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data; 5402bfc04deSHong Zhang 5412bfc04deSHong Zhang PetscFunctionBegin; 5422bfc04deSHong Zhang cvode->monitorstep = s; 5432bfc04deSHong Zhang PetscFunctionReturn(0); 5442bfc04deSHong Zhang } 5457cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 5467cb94ee6SHong Zhang 5477cb94ee6SHong Zhang /*@C 5487cb94ee6SHong Zhang TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 5497cb94ee6SHong Zhang 5507cb94ee6SHong Zhang Not Collective 5517cb94ee6SHong Zhang 55297bb3fdcSJose E. Roman Input Parameter: 5537cb94ee6SHong Zhang . ts - the time-step context 5547cb94ee6SHong Zhang 5557cb94ee6SHong Zhang Output Parameters: 5567cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 5577cb94ee6SHong Zhang - lin - number of linear iterations 5587cb94ee6SHong Zhang 5597cb94ee6SHong Zhang Level: advanced 5607cb94ee6SHong Zhang 5617cb94ee6SHong Zhang Notes: 5627cb94ee6SHong Zhang These return the number since the creation of the TS object 5637cb94ee6SHong Zhang 564db781477SPatrick Sanan .seealso: `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 565db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 566db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 567db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsGetPC()`, `TSSetExactFinalTime()` 5687cb94ee6SHong Zhang 5697cb94ee6SHong Zhang @*/ 5709371c9d4SSatish Balay PetscErrorCode TSSundialsGetIterations(TS ts, int *nonlin, int *lin) { 5717cb94ee6SHong Zhang PetscFunctionBegin; 572cac4c232SBarry Smith PetscUseMethod(ts, "TSSundialsGetIterations_C", (TS, int *, int *), (ts, nonlin, lin)); 5737cb94ee6SHong Zhang PetscFunctionReturn(0); 5747cb94ee6SHong Zhang } 5757cb94ee6SHong Zhang 5767cb94ee6SHong Zhang /*@ 5777cb94ee6SHong Zhang TSSundialsSetType - Sets the method that Sundials will use for integration. 5787cb94ee6SHong Zhang 5793f9fe445SBarry Smith Logically Collective on TS 5807cb94ee6SHong Zhang 5818f7d6fe5SPatrick Sanan Input Parameters: 5827cb94ee6SHong Zhang + ts - the time-step context 5837cb94ee6SHong Zhang - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 5847cb94ee6SHong Zhang 5857cb94ee6SHong Zhang Level: intermediate 5867cb94ee6SHong Zhang 587db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetMaxl()`, 588db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 589db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 590db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 591db781477SPatrick Sanan `TSSetExactFinalTime()` 5927cb94ee6SHong Zhang @*/ 5939371c9d4SSatish Balay PetscErrorCode TSSundialsSetType(TS ts, TSSundialsLmmType type) { 5947cb94ee6SHong Zhang PetscFunctionBegin; 595cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetType_C", (TS, TSSundialsLmmType), (ts, type)); 5967cb94ee6SHong Zhang PetscFunctionReturn(0); 5977cb94ee6SHong Zhang } 5987cb94ee6SHong Zhang 5997cb94ee6SHong Zhang /*@ 600c4e80e11SFlorian TSSundialsSetMaxord - Sets the maximum order for BDF/Adams method used by SUNDIALS. 601c4e80e11SFlorian 602c4e80e11SFlorian Logically Collective on TS 603c4e80e11SFlorian 6048f7d6fe5SPatrick Sanan Input Parameters: 605c4e80e11SFlorian + ts - the time-step context 606c4e80e11SFlorian - maxord - maximum order of BDF / Adams method 607c4e80e11SFlorian 608c4e80e11SFlorian Level: advanced 609c4e80e11SFlorian 610db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, 611db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 612db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 613db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 614db781477SPatrick Sanan `TSSetExactFinalTime()` 615c4e80e11SFlorian 616c4e80e11SFlorian @*/ 6179371c9d4SSatish Balay PetscErrorCode TSSundialsSetMaxord(TS ts, PetscInt maxord) { 618c4e80e11SFlorian PetscFunctionBegin; 619c4e80e11SFlorian PetscValidLogicalCollectiveInt(ts, maxord, 2); 620cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetMaxOrd_C", (TS, PetscInt), (ts, maxord)); 621c4e80e11SFlorian PetscFunctionReturn(0); 622c4e80e11SFlorian } 623c4e80e11SFlorian 624c4e80e11SFlorian /*@ 625f61b2b6aSHong Zhang TSSundialsSetMaxl - Sets the dimension of the Krylov space used by 6267cb94ee6SHong Zhang GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 627f61b2b6aSHong Zhang this is the maximum number of GMRES steps that will be used. 6287cb94ee6SHong Zhang 6293f9fe445SBarry Smith Logically Collective on TS 6307cb94ee6SHong Zhang 6318f7d6fe5SPatrick Sanan Input Parameters: 6327cb94ee6SHong Zhang + ts - the time-step context 633f61b2b6aSHong Zhang - maxl - number of direction vectors (the dimension of Krylov subspace). 6347cb94ee6SHong Zhang 6357cb94ee6SHong Zhang Level: advanced 6367cb94ee6SHong Zhang 637db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, 638db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 639db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 640db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 641db781477SPatrick Sanan `TSSetExactFinalTime()` 6427cb94ee6SHong Zhang 6437cb94ee6SHong Zhang @*/ 6449371c9d4SSatish Balay PetscErrorCode TSSundialsSetMaxl(TS ts, PetscInt maxl) { 6457cb94ee6SHong Zhang PetscFunctionBegin; 646f61b2b6aSHong Zhang PetscValidLogicalCollectiveInt(ts, maxl, 2); 647cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetMaxl_C", (TS, PetscInt), (ts, maxl)); 6487cb94ee6SHong Zhang PetscFunctionReturn(0); 6497cb94ee6SHong Zhang } 6507cb94ee6SHong Zhang 6517cb94ee6SHong Zhang /*@ 6527cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 6537cb94ee6SHong Zhang system by SUNDIALS. 6547cb94ee6SHong Zhang 6553f9fe445SBarry Smith Logically Collective on TS 6567cb94ee6SHong Zhang 6578f7d6fe5SPatrick Sanan Input Parameters: 6587cb94ee6SHong Zhang + ts - the time-step context 6597cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 6607cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 6617cb94ee6SHong Zhang 6627cb94ee6SHong Zhang Level: advanced 6637cb94ee6SHong Zhang 664db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 665db781477SPatrick Sanan `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 666db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 667db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 668db781477SPatrick Sanan `TSSetExactFinalTime()` 6697cb94ee6SHong Zhang 6707cb94ee6SHong Zhang @*/ 6719371c9d4SSatish Balay PetscErrorCode TSSundialsSetLinearTolerance(TS ts, PetscReal tol) { 6727cb94ee6SHong Zhang PetscFunctionBegin; 673c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts, tol, 2); 674cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetLinearTolerance_C", (TS, PetscReal), (ts, tol)); 6757cb94ee6SHong Zhang PetscFunctionReturn(0); 6767cb94ee6SHong Zhang } 6777cb94ee6SHong Zhang 6787cb94ee6SHong Zhang /*@ 6797cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 6807cb94ee6SHong Zhang in GMRES method by SUNDIALS linear solver. 6817cb94ee6SHong Zhang 6823f9fe445SBarry Smith Logically Collective on TS 6837cb94ee6SHong Zhang 6848f7d6fe5SPatrick Sanan Input Parameters: 6857cb94ee6SHong Zhang + ts - the time-step context 6867cb94ee6SHong Zhang - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 6877cb94ee6SHong Zhang 6887cb94ee6SHong Zhang Level: advanced 6897cb94ee6SHong Zhang 690db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 691db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, 692db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 693db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 694db781477SPatrick Sanan `TSSetExactFinalTime()` 6957cb94ee6SHong Zhang 6967cb94ee6SHong Zhang @*/ 6979371c9d4SSatish Balay PetscErrorCode TSSundialsSetGramSchmidtType(TS ts, TSSundialsGramSchmidtType type) { 6987cb94ee6SHong Zhang PetscFunctionBegin; 699cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetGramSchmidtType_C", (TS, TSSundialsGramSchmidtType), (ts, type)); 7007cb94ee6SHong Zhang PetscFunctionReturn(0); 7017cb94ee6SHong Zhang } 7027cb94ee6SHong Zhang 7037cb94ee6SHong Zhang /*@ 7047cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 7057cb94ee6SHong Zhang Sundials for error control. 7067cb94ee6SHong Zhang 7073f9fe445SBarry Smith Logically Collective on TS 7087cb94ee6SHong Zhang 7098f7d6fe5SPatrick Sanan Input Parameters: 7107cb94ee6SHong Zhang + ts - the time-step context 7117cb94ee6SHong Zhang . aabs - the absolute tolerance 7127cb94ee6SHong Zhang - rel - the relative tolerance 7137cb94ee6SHong Zhang 7147cb94ee6SHong Zhang See the Cvode/Sundials users manual for exact details on these parameters. Essentially 7157cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 7167cb94ee6SHong Zhang 7177cb94ee6SHong Zhang Level: intermediate 7187cb94ee6SHong Zhang 719db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetGMRESMaxl()`, 720db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, 721db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 722db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, 723db781477SPatrick Sanan `TSSetExactFinalTime()` 7247cb94ee6SHong Zhang 7257cb94ee6SHong Zhang @*/ 7269371c9d4SSatish Balay PetscErrorCode TSSundialsSetTolerance(TS ts, PetscReal aabs, PetscReal rel) { 7277cb94ee6SHong Zhang PetscFunctionBegin; 728cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetTolerance_C", (TS, PetscReal, PetscReal), (ts, aabs, rel)); 7297cb94ee6SHong Zhang PetscFunctionReturn(0); 7307cb94ee6SHong Zhang } 7317cb94ee6SHong Zhang 7327cb94ee6SHong Zhang /*@ 7337cb94ee6SHong Zhang TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 7347cb94ee6SHong Zhang 7357cb94ee6SHong Zhang Input Parameter: 7367cb94ee6SHong Zhang . ts - the time-step context 7377cb94ee6SHong Zhang 7387cb94ee6SHong Zhang Output Parameter: 7397cb94ee6SHong Zhang . pc - the preconditioner context 7407cb94ee6SHong Zhang 7417cb94ee6SHong Zhang Level: advanced 7427cb94ee6SHong Zhang 743db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, 744db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, 745db781477SPatrick Sanan `TSSundialsGetIterations()`, `TSSundialsSetType()`, 746db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()` 7477cb94ee6SHong Zhang @*/ 7489371c9d4SSatish Balay PetscErrorCode TSSundialsGetPC(TS ts, PC *pc) { 7497cb94ee6SHong Zhang PetscFunctionBegin; 750cac4c232SBarry Smith PetscUseMethod(ts, "TSSundialsGetPC_C", (TS, PC *), (ts, pc)); 7517cb94ee6SHong Zhang PetscFunctionReturn(0); 7527cb94ee6SHong Zhang } 7537cb94ee6SHong Zhang 754f1cd61daSJed Brown /*@ 755f1cd61daSJed Brown TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 756f1cd61daSJed Brown 757d8d19677SJose E. Roman Input Parameters: 758f1cd61daSJed Brown + ts - the time-step context 759f1cd61daSJed Brown - mindt - lowest time step if positive, negative to deactivate 760f1cd61daSJed Brown 761fc6b9e64SJed Brown Note: 762fc6b9e64SJed Brown Sundials will error if it is not possible to keep the estimated truncation error below 763fc6b9e64SJed Brown the tolerance set with TSSundialsSetTolerance() without going below this step size. 764fc6b9e64SJed Brown 765f1cd61daSJed Brown Level: beginner 766f1cd61daSJed Brown 767db781477SPatrick Sanan .seealso: `TSSundialsSetType()`, `TSSundialsSetTolerance()`, 768f1cd61daSJed Brown @*/ 7699371c9d4SSatish Balay PetscErrorCode TSSundialsSetMinTimeStep(TS ts, PetscReal mindt) { 770f1cd61daSJed Brown PetscFunctionBegin; 771cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetMinTimeStep_C", (TS, PetscReal), (ts, mindt)); 772f1cd61daSJed Brown PetscFunctionReturn(0); 773f1cd61daSJed Brown } 774f1cd61daSJed Brown 775f1cd61daSJed Brown /*@ 776f1cd61daSJed Brown TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 777f1cd61daSJed Brown 778d8d19677SJose E. Roman Input Parameters: 779f1cd61daSJed Brown + ts - the time-step context 780f1cd61daSJed Brown - maxdt - lowest time step if positive, negative to deactivate 781f1cd61daSJed Brown 782f1cd61daSJed Brown Level: beginner 783f1cd61daSJed Brown 784db781477SPatrick Sanan .seealso: `TSSundialsSetType()`, `TSSundialsSetTolerance()`, 785f1cd61daSJed Brown @*/ 7869371c9d4SSatish Balay PetscErrorCode TSSundialsSetMaxTimeStep(TS ts, PetscReal maxdt) { 787f1cd61daSJed Brown PetscFunctionBegin; 788cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt)); 789f1cd61daSJed Brown PetscFunctionReturn(0); 790f1cd61daSJed Brown } 791f1cd61daSJed Brown 7922bfc04deSHong Zhang /*@ 7932bfc04deSHong Zhang TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 7942bfc04deSHong Zhang 795d8d19677SJose E. Roman Input Parameters: 7962bfc04deSHong Zhang + ts - the time-step context 7972bfc04deSHong Zhang - ft - PETSC_TRUE if monitor, else PETSC_FALSE 7982bfc04deSHong Zhang 7992bfc04deSHong Zhang Level: beginner 8002bfc04deSHong Zhang 801f61b2b6aSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 8022bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 803f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 8042bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 8052bfc04deSHong Zhang @*/ 8069371c9d4SSatish Balay PetscErrorCode TSSundialsMonitorInternalSteps(TS ts, PetscBool ft) { 8072bfc04deSHong Zhang PetscFunctionBegin; 808cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsMonitorInternalSteps_C", (TS, PetscBool), (ts, ft)); 8092bfc04deSHong Zhang PetscFunctionReturn(0); 8102bfc04deSHong Zhang } 811918687b8SPatrick Sanan 812918687b8SPatrick Sanan /*@ 813918687b8SPatrick Sanan TSSundialsSetUseDense - Set a flag to use a dense linear solver in SUNDIALS (serial only) 814918687b8SPatrick Sanan 815918687b8SPatrick Sanan Logically Collective 816918687b8SPatrick Sanan 817918687b8SPatrick Sanan Input Parameters: 818918687b8SPatrick Sanan + ts - the time-step context 819918687b8SPatrick Sanan - use_dense - PETSC_TRUE to use the dense solver 820918687b8SPatrick Sanan 821918687b8SPatrick Sanan Level: advanced 822918687b8SPatrick Sanan 823db781477SPatrick Sanan .seealso: `TSSUNDIALS` 824918687b8SPatrick Sanan 825918687b8SPatrick Sanan @*/ 8269371c9d4SSatish Balay PetscErrorCode TSSundialsSetUseDense(TS ts, PetscBool use_dense) { 827918687b8SPatrick Sanan PetscFunctionBegin; 828918687b8SPatrick Sanan PetscValidLogicalCollectiveInt(ts, use_dense, 2); 829cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetUseDense_C", (TS, PetscBool), (ts, use_dense)); 830918687b8SPatrick Sanan PetscFunctionReturn(0); 831918687b8SPatrick Sanan } 832918687b8SPatrick Sanan 8337cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 8347cb94ee6SHong Zhang /*MC 83596f5712cSJed Brown TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 8367cb94ee6SHong Zhang 8377cb94ee6SHong Zhang Options Database: 838897c9f78SHong Zhang + -ts_sundials_type <bdf,adams> - 8397cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 8407cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 8417cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 842897c9f78SHong Zhang . -ts_sundials_linear_tolerance <tol> - 843f61b2b6aSHong Zhang . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace 844918687b8SPatrick Sanan . -ts_sundials_monitor_steps - Monitor SUNDIALS internal steps 845918687b8SPatrick Sanan - -ts_sundials_use_dense - Use a dense linear solver within CVODE (serial only) 84616016685SHong Zhang 84795452b02SPatrick Sanan Notes: 84895452b02SPatrick Sanan This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply, 849897c9f78SHong Zhang only PETSc PC options. 8507cb94ee6SHong Zhang 8517cb94ee6SHong Zhang Level: beginner 8527cb94ee6SHong Zhang 853db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, `TSSundialsSetLinearTolerance()`, 854db781477SPatrick Sanan `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, `TSSundialsGetIterations()`, `TSSetExactFinalTime()` 8557cb94ee6SHong Zhang 8567cb94ee6SHong Zhang M*/ 8579371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts) { 8587cb94ee6SHong Zhang TS_Sundials *cvode; 85942528757SHong Zhang PC pc; 8607cb94ee6SHong Zhang 8617cb94ee6SHong Zhang PetscFunctionBegin; 862277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Sundials; 86328adb3f7SHong Zhang ts->ops->destroy = TSDestroy_Sundials; 86428adb3f7SHong Zhang ts->ops->view = TSView_Sundials; 865214bc6a2SJed Brown ts->ops->setup = TSSetUp_Sundials; 866b4eba00bSSean Farley ts->ops->step = TSStep_Sundials; 867b4eba00bSSean Farley ts->ops->interpolate = TSInterpolate_Sundials; 868214bc6a2SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Sundials; 8692ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 8707cb94ee6SHong Zhang 871*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&cvode)); 872bbd56ea5SKarl Rupp 873efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 874efd4aadfSBarry Smith 8757cb94ee6SHong Zhang ts->data = (void *)cvode; 8766fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF; 8776fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS; 878f61b2b6aSHong Zhang cvode->maxl = 5; 879c4e80e11SFlorian cvode->maxord = PETSC_DEFAULT; 8807cb94ee6SHong Zhang cvode->linear_tol = .05; 881b4eba00bSSean Farley cvode->monitorstep = PETSC_TRUE; 882918687b8SPatrick Sanan cvode->use_dense = PETSC_FALSE; 8837cb94ee6SHong Zhang 8849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)ts), &(cvode->comm_sundials))); 885f1cd61daSJed Brown 886f1cd61daSJed Brown cvode->mindt = -1.; 887f1cd61daSJed Brown cvode->maxdt = -1.; 888f1cd61daSJed Brown 8897cb94ee6SHong Zhang /* set tolerance for Sundials */ 8907cb94ee6SHong Zhang cvode->reltol = 1e-6; 8912c823083SHong Zhang cvode->abstol = 1e-6; 8927cb94ee6SHong Zhang 89342528757SHong Zhang /* set PCNONE as default pctype */ 8949566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC_Sundials(ts, &pc)); 8959566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 89642528757SHong Zhang 8979566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetType_C", TSSundialsSetType_Sundials)); 8989566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxl_C", TSSundialsSetMaxl_Sundials)); 8999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetLinearTolerance_C", TSSundialsSetLinearTolerance_Sundials)); 9009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetGramSchmidtType_C", TSSundialsSetGramSchmidtType_Sundials)); 9019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetTolerance_C", TSSundialsSetTolerance_Sundials)); 9029566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMinTimeStep_C", TSSundialsSetMinTimeStep_Sundials)); 9039566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxTimeStep_C", TSSundialsSetMaxTimeStep_Sundials)); 9049566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetPC_C", TSSundialsGetPC_Sundials)); 9059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetIterations_C", TSSundialsGetIterations_Sundials)); 9069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsMonitorInternalSteps_C", TSSundialsMonitorInternalSteps_Sundials)); 9079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetUseDense_C", TSSundialsSetUseDense_Sundials)); 9087cb94ee6SHong Zhang PetscFunctionReturn(0); 9097cb94ee6SHong Zhang } 910