1c4762a1bSJed Brown static char help[] = "Nonlinear DAE benchmark problems.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this 5c4762a1bSJed Brown file automatically includes: 6c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 7c4762a1bSJed Brown petscmat.h - matrices 8c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 9c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 10c4762a1bSJed Brown petscksp.h - linear solvers 11c4762a1bSJed Brown */ 12c4762a1bSJed Brown #include <petscts.h> 13c4762a1bSJed Brown 14c4762a1bSJed Brown typedef struct _Problem *Problem; 15c4762a1bSJed Brown struct _Problem { 16c4762a1bSJed Brown PetscErrorCode (*destroy)(Problem); 17*8434afd1SBarry Smith TSIFunctionFn *function; 18*8434afd1SBarry Smith TSIJacobianFn *jacobian; 19c4762a1bSJed Brown PetscErrorCode (*solution)(PetscReal, Vec, void *); 20c4762a1bSJed Brown MPI_Comm comm; 21c4762a1bSJed Brown PetscReal final_time; 22c4762a1bSJed Brown PetscInt n; 23c4762a1bSJed Brown PetscBool hasexact; 24c4762a1bSJed Brown void *data; 25c4762a1bSJed Brown }; 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* 28c4762a1bSJed Brown Stiff 3-variable system from chemical reactions, due to Robertson (1966), problem ROBER in Hairer&Wanner, ODE 2, 1996 29c4762a1bSJed Brown */ 30d71ae5a4SJacob Faibussowitsch static PetscErrorCode RoberFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) 31d71ae5a4SJacob Faibussowitsch { 32c4762a1bSJed Brown PetscScalar *f; 33c4762a1bSJed Brown const PetscScalar *x, *xdot; 34c4762a1bSJed Brown 35c4762a1bSJed Brown PetscFunctionBeginUser; 369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot, &xdot)); 389566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 39c4762a1bSJed Brown f[0] = xdot[0] + 0.04 * x[0] - 1e4 * x[1] * x[2]; 40c4762a1bSJed Brown f[1] = xdot[1] - 0.04 * x[0] + 1e4 * x[1] * x[2] + 3e7 * PetscSqr(x[1]); 41c4762a1bSJed Brown f[2] = xdot[2] - 3e7 * PetscSqr(x[1]); 429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot, &xdot)); 449566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 48d71ae5a4SJacob Faibussowitsch static PetscErrorCode RoberJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) 49d71ae5a4SJacob Faibussowitsch { 50c4762a1bSJed Brown PetscInt rowcol[] = {0, 1, 2}; 51c4762a1bSJed Brown PetscScalar J[3][3]; 52c4762a1bSJed Brown const PetscScalar *x, *xdot; 53c4762a1bSJed Brown 54c4762a1bSJed Brown PetscFunctionBeginUser; 559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot, &xdot)); 579371c9d4SSatish Balay J[0][0] = a + 0.04; 589371c9d4SSatish Balay J[0][1] = -1e4 * x[2]; 599371c9d4SSatish Balay J[0][2] = -1e4 * x[1]; 609371c9d4SSatish Balay J[1][0] = -0.04; 619371c9d4SSatish Balay J[1][1] = a + 1e4 * x[2] + 3e7 * 2 * x[1]; 629371c9d4SSatish Balay J[1][2] = 1e4 * x[1]; 639371c9d4SSatish Balay J[2][0] = 0; 649371c9d4SSatish Balay J[2][1] = -3e7 * 2 * x[1]; 659371c9d4SSatish Balay J[2][2] = a; 669566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 3, rowcol, 3, rowcol, &J[0][0], INSERT_VALUES)); 679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 689566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot, &xdot)); 69c4762a1bSJed Brown 709566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 719566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 72c4762a1bSJed Brown if (A != B) { 739566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 749566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 75c4762a1bSJed Brown } 763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown 79d71ae5a4SJacob Faibussowitsch static PetscErrorCode RoberSolution(PetscReal t, Vec X, void *ctx) 80d71ae5a4SJacob Faibussowitsch { 81c4762a1bSJed Brown PetscScalar *x; 82c4762a1bSJed Brown 83c4762a1bSJed Brown PetscFunctionBeginUser; 843c633725SBarry Smith PetscCheck(t == 0, PETSC_COMM_WORLD, PETSC_ERR_SUP, "not implemented"); 859566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 86c4762a1bSJed Brown x[0] = 1; 87c4762a1bSJed Brown x[1] = 0; 88c4762a1bSJed Brown x[2] = 0; 899566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 91c4762a1bSJed Brown } 92c4762a1bSJed Brown 93d71ae5a4SJacob Faibussowitsch static PetscErrorCode RoberCreate(Problem p) 94d71ae5a4SJacob Faibussowitsch { 95c4762a1bSJed Brown PetscFunctionBeginUser; 96c4762a1bSJed Brown p->destroy = 0; 97c4762a1bSJed Brown p->function = &RoberFunction; 98c4762a1bSJed Brown p->jacobian = &RoberJacobian; 99c4762a1bSJed Brown p->solution = &RoberSolution; 100c4762a1bSJed Brown p->final_time = 1e11; 101c4762a1bSJed Brown p->n = 3; 1023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* 106c4762a1bSJed Brown Stiff scalar valued problem 107c4762a1bSJed Brown */ 108c4762a1bSJed Brown 109c4762a1bSJed Brown typedef struct { 110c4762a1bSJed Brown PetscReal lambda; 111c4762a1bSJed Brown } CECtx; 112c4762a1bSJed Brown 113d71ae5a4SJacob Faibussowitsch static PetscErrorCode CEDestroy(Problem p) 114d71ae5a4SJacob Faibussowitsch { 115c4762a1bSJed Brown PetscFunctionBeginUser; 1169566063dSJacob Faibussowitsch PetscCall(PetscFree(p->data)); 1173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 118c4762a1bSJed Brown } 119c4762a1bSJed Brown 120d71ae5a4SJacob Faibussowitsch static PetscErrorCode CEFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) 121d71ae5a4SJacob Faibussowitsch { 122c4762a1bSJed Brown PetscReal l = ((CECtx *)ctx)->lambda; 123c4762a1bSJed Brown PetscScalar *f; 124c4762a1bSJed Brown const PetscScalar *x, *xdot; 125c4762a1bSJed Brown 126c4762a1bSJed Brown PetscFunctionBeginUser; 1279566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 1289566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot, &xdot)); 1299566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 130c4762a1bSJed Brown f[0] = xdot[0] + l * (x[0] - PetscCosReal(t)); 131c4762a1bSJed Brown #if 0 1329566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," f(t=%g,x=%g,xdot=%g) = %g\n",(double)t,(double)x[0],(double)xdot[0],(double)f[0])); 133c4762a1bSJed Brown #endif 1349566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 1359566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot, &xdot)); 1369566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 1373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 138c4762a1bSJed Brown } 139c4762a1bSJed Brown 140d71ae5a4SJacob Faibussowitsch static PetscErrorCode CEJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) 141d71ae5a4SJacob Faibussowitsch { 142c4762a1bSJed Brown PetscReal l = ((CECtx *)ctx)->lambda; 143c4762a1bSJed Brown PetscInt rowcol[] = {0}; 144c4762a1bSJed Brown PetscScalar J[1][1]; 145c4762a1bSJed Brown const PetscScalar *x, *xdot; 146c4762a1bSJed Brown 147c4762a1bSJed Brown PetscFunctionBeginUser; 1489566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 1499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot, &xdot)); 150c4762a1bSJed Brown J[0][0] = a + l; 1519566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, rowcol, 1, rowcol, &J[0][0], INSERT_VALUES)); 1529566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 1539566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot, &xdot)); 154c4762a1bSJed Brown 1559566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1569566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 157c4762a1bSJed Brown if (A != B) { 1589566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1599566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 160c4762a1bSJed Brown } 1613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 162c4762a1bSJed Brown } 163c4762a1bSJed Brown 164d71ae5a4SJacob Faibussowitsch static PetscErrorCode CESolution(PetscReal t, Vec X, void *ctx) 165d71ae5a4SJacob Faibussowitsch { 166c4762a1bSJed Brown PetscReal l = ((CECtx *)ctx)->lambda; 167c4762a1bSJed Brown PetscScalar *x; 168c4762a1bSJed Brown 169c4762a1bSJed Brown PetscFunctionBeginUser; 1709566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 171c4762a1bSJed Brown x[0] = l / (l * l + 1) * (l * PetscCosReal(t) + PetscSinReal(t)) - l * l / (l * l + 1) * PetscExpReal(-l * t); 1729566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 1733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176d71ae5a4SJacob Faibussowitsch static PetscErrorCode CECreate(Problem p) 177d71ae5a4SJacob Faibussowitsch { 178c4762a1bSJed Brown CECtx *ce; 179c4762a1bSJed Brown 180c4762a1bSJed Brown PetscFunctionBeginUser; 1819566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(CECtx), &ce)); 182c4762a1bSJed Brown p->data = (void *)ce; 183c4762a1bSJed Brown 184c4762a1bSJed Brown p->destroy = &CEDestroy; 185c4762a1bSJed Brown p->function = &CEFunction; 186c4762a1bSJed Brown p->jacobian = &CEJacobian; 187c4762a1bSJed Brown p->solution = &CESolution; 188c4762a1bSJed Brown p->final_time = 10; 189c4762a1bSJed Brown p->n = 1; 190c4762a1bSJed Brown p->hasexact = PETSC_TRUE; 191c4762a1bSJed Brown 192c4762a1bSJed Brown ce->lambda = 10; 193d0609cedSBarry Smith PetscOptionsBegin(p->comm, NULL, "CE options", ""); 194d71ae5a4SJacob Faibussowitsch { 195d71ae5a4SJacob Faibussowitsch PetscCall(PetscOptionsReal("-problem_ce_lambda", "Parameter controlling stiffness: xdot + lambda*(x - cos(t))", "", ce->lambda, &ce->lambda, NULL)); 196d71ae5a4SJacob Faibussowitsch } 197d0609cedSBarry Smith PetscOptionsEnd(); 1983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 199c4762a1bSJed Brown } 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* 2020e3d61c9SBarry Smith Stiff 3-variable oscillatory system from chemical reactions. problem OREGO in Hairer&Wanner 203c4762a1bSJed Brown */ 204d71ae5a4SJacob Faibussowitsch static PetscErrorCode OregoFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) 205d71ae5a4SJacob Faibussowitsch { 206c4762a1bSJed Brown PetscScalar *f; 207c4762a1bSJed Brown const PetscScalar *x, *xdot; 208c4762a1bSJed Brown 209c4762a1bSJed Brown PetscFunctionBeginUser; 2109566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 2119566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot, &xdot)); 2129566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 213c4762a1bSJed Brown f[0] = xdot[0] - 77.27 * (x[1] + x[0] * (1. - 8.375e-6 * x[0] - x[1])); 214c4762a1bSJed Brown f[1] = xdot[1] - 1 / 77.27 * (x[2] - (1. + x[0]) * x[1]); 215c4762a1bSJed Brown f[2] = xdot[2] - 0.161 * (x[0] - x[2]); 2169566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 2179566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot, &xdot)); 2189566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 2193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 220c4762a1bSJed Brown } 221c4762a1bSJed Brown 222d71ae5a4SJacob Faibussowitsch static PetscErrorCode OregoJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) 223d71ae5a4SJacob Faibussowitsch { 224c4762a1bSJed Brown PetscInt rowcol[] = {0, 1, 2}; 225c4762a1bSJed Brown PetscScalar J[3][3]; 226c4762a1bSJed Brown const PetscScalar *x, *xdot; 227c4762a1bSJed Brown 228c4762a1bSJed Brown PetscFunctionBeginUser; 2299566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 2309566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot, &xdot)); 231c4762a1bSJed Brown J[0][0] = a - 77.27 * ((1. - 8.375e-6 * x[0] - x[1]) - 8.375e-6 * x[0]); 232c4762a1bSJed Brown J[0][1] = -77.27 * (1. - x[0]); 233c4762a1bSJed Brown J[0][2] = 0; 234c4762a1bSJed Brown J[1][0] = 1. / 77.27 * x[1]; 235c4762a1bSJed Brown J[1][1] = a + 1. / 77.27 * (1. + x[0]); 236c4762a1bSJed Brown J[1][2] = -1. / 77.27; 237c4762a1bSJed Brown J[2][0] = -0.161; 238c4762a1bSJed Brown J[2][1] = 0; 239c4762a1bSJed Brown J[2][2] = a + 0.161; 2409566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 3, rowcol, 3, rowcol, &J[0][0], INSERT_VALUES)); 2419566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 2429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot, &xdot)); 243c4762a1bSJed Brown 2449566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 2459566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 246c4762a1bSJed Brown if (A != B) { 2479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 2489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 249c4762a1bSJed Brown } 2503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 251c4762a1bSJed Brown } 252c4762a1bSJed Brown 253d71ae5a4SJacob Faibussowitsch static PetscErrorCode OregoSolution(PetscReal t, Vec X, void *ctx) 254d71ae5a4SJacob Faibussowitsch { 255c4762a1bSJed Brown PetscScalar *x; 256c4762a1bSJed Brown 257c4762a1bSJed Brown PetscFunctionBeginUser; 2583c633725SBarry Smith PetscCheck(t == 0, PETSC_COMM_WORLD, PETSC_ERR_SUP, "not implemented"); 2599566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 260c4762a1bSJed Brown x[0] = 1; 261c4762a1bSJed Brown x[1] = 2; 262c4762a1bSJed Brown x[2] = 3; 2639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 2643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 265c4762a1bSJed Brown } 266c4762a1bSJed Brown 267d71ae5a4SJacob Faibussowitsch static PetscErrorCode OregoCreate(Problem p) 268d71ae5a4SJacob Faibussowitsch { 269c4762a1bSJed Brown PetscFunctionBeginUser; 270c4762a1bSJed Brown p->destroy = 0; 271c4762a1bSJed Brown p->function = &OregoFunction; 272c4762a1bSJed Brown p->jacobian = &OregoJacobian; 273c4762a1bSJed Brown p->solution = &OregoSolution; 274c4762a1bSJed Brown p->final_time = 360; 275c4762a1bSJed Brown p->n = 3; 2763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 277c4762a1bSJed Brown } 278c4762a1bSJed Brown 279c4762a1bSJed Brown /* 2800e3d61c9SBarry Smith User-defined monitor for comparing to exact solutions when possible 281c4762a1bSJed Brown */ 282c4762a1bSJed Brown typedef struct { 283c4762a1bSJed Brown MPI_Comm comm; 284c4762a1bSJed Brown Problem problem; 285c4762a1bSJed Brown Vec x; 286c4762a1bSJed Brown } MonitorCtx; 287c4762a1bSJed Brown 288d71ae5a4SJacob Faibussowitsch static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal t, Vec x, void *ctx) 289d71ae5a4SJacob Faibussowitsch { 290c4762a1bSJed Brown MonitorCtx *mon = (MonitorCtx *)ctx; 291c4762a1bSJed Brown PetscReal h, nrm_x, nrm_exact, nrm_diff; 292c4762a1bSJed Brown 293c4762a1bSJed Brown PetscFunctionBeginUser; 2943ba16761SJacob Faibussowitsch if (!mon->problem->solution) PetscFunctionReturn(PETSC_SUCCESS); 2959566063dSJacob Faibussowitsch PetscCall((*mon->problem->solution)(t, mon->x, mon->problem->data)); 2969566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &nrm_x)); 2979566063dSJacob Faibussowitsch PetscCall(VecNorm(mon->x, NORM_2, &nrm_exact)); 2989566063dSJacob Faibussowitsch PetscCall(VecAYPX(mon->x, -1, x)); 2999566063dSJacob Faibussowitsch PetscCall(VecNorm(mon->x, NORM_2, &nrm_diff)); 3009566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &h)); 30148a46eb9SPierre Jolivet if (step < 0) PetscCall(PetscPrintf(mon->comm, "Interpolated final solution ")); 30263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(mon->comm, "step %4" PetscInt_FMT " t=%12.8e h=% 8.2e |x|=%9.2e |x_e|=%9.2e |x-x_e|=%9.2e\n", step, (double)t, (double)h, (double)nrm_x, (double)nrm_exact, (double)nrm_diff)); 3033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 304c4762a1bSJed Brown } 305c4762a1bSJed Brown 306d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 307d71ae5a4SJacob Faibussowitsch { 308c4762a1bSJed Brown PetscFunctionList plist = NULL; 309c4762a1bSJed Brown char pname[256]; 310c4762a1bSJed Brown TS ts; /* nonlinear solver */ 311c4762a1bSJed Brown Vec x, r; /* solution, residual vectors */ 312c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 313c4762a1bSJed Brown Problem problem; 3143d5a8a6aSBarry Smith PetscBool use_monitor = PETSC_FALSE; 3153d5a8a6aSBarry Smith PetscBool use_result = PETSC_FALSE; 316c4762a1bSJed Brown PetscInt steps, nonlinits, linits, snesfails, rejects; 317c4762a1bSJed Brown PetscReal ftime; 318c4762a1bSJed Brown MonitorCtx mon; 319c4762a1bSJed Brown PetscMPIInt size; 320c4762a1bSJed Brown 321c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 322c4762a1bSJed Brown Initialize program 323c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 324327415f7SBarry Smith PetscFunctionBeginUser; 3259566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 3269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 3273c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 328c4762a1bSJed Brown 329c4762a1bSJed Brown /* Register the available problems */ 3309566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&plist, "rober", &RoberCreate)); 3319566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&plist, "ce", &CECreate)); 3329566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&plist, "orego", &OregoCreate)); 333c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(pname, "ce", sizeof(pname))); 334c4762a1bSJed Brown 335c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 336c4762a1bSJed Brown Set runtime options 337c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 338d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Timestepping benchmark options", ""); 339c4762a1bSJed Brown { 3409566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-problem_type", "Name of problem to run", "", plist, pname, pname, sizeof(pname), NULL)); 341c4762a1bSJed Brown use_monitor = PETSC_FALSE; 3429566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-monitor_error", "Display errors relative to exact solutions", "", use_monitor, &use_monitor, NULL)); 3439566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-monitor_result", "Display result", "", use_result, &use_result, NULL)); 344c4762a1bSJed Brown } 345d0609cedSBarry Smith PetscOptionsEnd(); 346c4762a1bSJed Brown 347c4762a1bSJed Brown /* Create the new problem */ 3489566063dSJacob Faibussowitsch PetscCall(PetscNew(&problem)); 349c4762a1bSJed Brown problem->comm = MPI_COMM_WORLD; 350c4762a1bSJed Brown { 351c4762a1bSJed Brown PetscErrorCode (*pcreate)(Problem); 352c4762a1bSJed Brown 3539566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(plist, pname, &pcreate)); 3543c633725SBarry Smith PetscCheck(pcreate, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "No problem '%s'", pname); 3559566063dSJacob Faibussowitsch PetscCall((*pcreate)(problem)); 356c4762a1bSJed Brown } 357c4762a1bSJed Brown 358c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 359c4762a1bSJed Brown Create necessary matrix and vectors 360c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3619566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 3629566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, problem->n, problem->n, PETSC_DETERMINE, PETSC_DETERMINE)); 3639566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 3649566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 365c4762a1bSJed Brown 3669566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, NULL)); 3679566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 368c4762a1bSJed Brown 369c4762a1bSJed Brown mon.comm = PETSC_COMM_WORLD; 370c4762a1bSJed Brown mon.problem = problem; 3719566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &mon.x)); 372c4762a1bSJed Brown 373c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 374c4762a1bSJed Brown Create timestepping solver context 375c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3769566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 3779566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 3789566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSROSW)); /* Rosenbrock-W */ 3799566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, problem->function, problem->data)); 3809566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, problem->jacobian, problem->data)); 3819566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, problem->final_time)); 3829566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 3839566063dSJacob Faibussowitsch PetscCall(TSSetMaxStepRejections(ts, 10)); 3849566063dSJacob Faibussowitsch PetscCall(TSSetMaxSNESFailures(ts, -1)); /* unlimited */ 38548a46eb9SPierre Jolivet if (use_monitor) PetscCall(TSMonitorSet(ts, &MonitorError, &mon, NULL)); 386c4762a1bSJed Brown 387c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 388c4762a1bSJed Brown Set initial conditions 389c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3909566063dSJacob Faibussowitsch PetscCall((*problem->solution)(0, x, problem->data)); 3919566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .001)); 3929566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x)); 393c4762a1bSJed Brown 394c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 395c4762a1bSJed Brown Set runtime options 396c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3979566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 398c4762a1bSJed Brown 399c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 400c4762a1bSJed Brown Solve nonlinear system 401c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 4029566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 4039566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 4049566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 4059566063dSJacob Faibussowitsch PetscCall(TSGetSNESFailures(ts, &snesfails)); 4069566063dSJacob Faibussowitsch PetscCall(TSGetStepRejections(ts, &rejects)); 4079566063dSJacob Faibussowitsch PetscCall(TSGetSNESIterations(ts, &nonlinits)); 4089566063dSJacob Faibussowitsch PetscCall(TSGetKSPIterations(ts, &linits)); 4093d5a8a6aSBarry Smith if (use_result) { 41063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "steps %" PetscInt_FMT " (%" PetscInt_FMT " rejected, %" PetscInt_FMT " SNES fails), ftime %g, nonlinits %" PetscInt_FMT ", linits %" PetscInt_FMT "\n", steps, rejects, snesfails, (double)ftime, nonlinits, linits)); 4113d5a8a6aSBarry Smith } 412c4762a1bSJed Brown 413c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 414c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 415c4762a1bSJed Brown are no longer needed. 416c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 4179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 4189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 4199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 4209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mon.x)); 4219566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 4221baa6e33SBarry Smith if (problem->destroy) PetscCall((*problem->destroy)(problem)); 4239566063dSJacob Faibussowitsch PetscCall(PetscFree(problem)); 4249566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&plist)); 425c4762a1bSJed Brown 4269566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 427b122ec5aSJacob Faibussowitsch return 0; 428c4762a1bSJed Brown } 429c4762a1bSJed Brown 430c4762a1bSJed Brown /*TEST 431c4762a1bSJed Brown 432c4762a1bSJed Brown test: 433c4762a1bSJed Brown requires: !complex 4343d5a8a6aSBarry Smith args: -monitor_result -monitor_error -ts_atol 1e-2 -ts_rtol 1e-2 -ts_exact_final_time interpolate -ts_type arkimex 435c4762a1bSJed Brown 436c4762a1bSJed Brown test: 437c4762a1bSJed Brown suffix: 2 438c4762a1bSJed Brown requires: !single !complex 4393d5a8a6aSBarry Smith args: -monitor_result -ts_atol 1e-2 -ts_rtol 1e-2 -ts_max_time 15 -ts_type arkimex -ts_arkimex_type 2e -problem_type orego -ts_arkimex_initial_guess_extrapolate 0 -ts_adapt_time_step_increase_delay 4 440c4762a1bSJed Brown 441c4762a1bSJed Brown test: 442c4762a1bSJed Brown suffix: 3 443c4762a1bSJed Brown requires: !single !complex 4443d5a8a6aSBarry Smith args: -monitor_result -ts_atol 1e-2 -ts_rtol 1e-2 -ts_max_time 15 -ts_type arkimex -ts_arkimex_type 2e -problem_type orego -ts_arkimex_initial_guess_extrapolate 1 4453d5a8a6aSBarry Smith 4463d5a8a6aSBarry Smith test: 4473d5a8a6aSBarry Smith suffix: 4 4483d5a8a6aSBarry Smith 4493d5a8a6aSBarry Smith test: 4503d5a8a6aSBarry Smith suffix: 5 4513d5a8a6aSBarry Smith args: -snes_lag_jacobian 20 -snes_lag_jacobian_persists 452c4762a1bSJed Brown 453c4762a1bSJed Brown TEST*/ 454