1 static char help[] = "Solves the trivial ODE 2 du/dt = 1, u(0) = 0. \n\n"; 2 3 #include <petscts.h> 4 #include <petscpc.h> 5 6 PetscErrorCode IFunction(TS, PetscReal, Vec, Vec, Vec, void *); 7 PetscErrorCode IJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 8 9 int main(int argc, char **argv) { 10 TS ts; 11 Vec x; 12 Vec f; 13 Mat A; 14 15 PetscFunctionBeginUser; 16 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 17 18 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 19 PetscCall(TSSetEquationType(ts, TS_EQ_ODE_IMPLICIT)); 20 PetscCall(VecCreate(PETSC_COMM_WORLD, &f)); 21 PetscCall(VecSetSizes(f, 1, PETSC_DECIDE)); 22 PetscCall(VecSetFromOptions(f)); 23 PetscCall(VecSetUp(f)); 24 PetscCall(TSSetIFunction(ts, f, IFunction, NULL)); 25 PetscCall(VecDestroy(&f)); 26 27 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 28 PetscCall(MatSetSizes(A, 1, 1, PETSC_DECIDE, PETSC_DECIDE)); 29 PetscCall(MatSetFromOptions(A)); 30 PetscCall(MatSetUp(A)); 31 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 32 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 33 /* ensure that the Jacobian matrix has diagonal entries since that is required by TS */ 34 PetscCall(MatShift(A, (PetscReal)1)); 35 PetscCall(MatShift(A, (PetscReal)-1)); 36 PetscCall(TSSetIJacobian(ts, A, A, IJacobian, NULL)); 37 PetscCall(MatDestroy(&A)); 38 39 PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 40 PetscCall(VecSetSizes(x, 1, PETSC_DECIDE)); 41 PetscCall(VecSetFromOptions(x)); 42 PetscCall(VecSetUp(x)); 43 PetscCall(TSSetSolution(ts, x)); 44 PetscCall(VecDestroy(&x)); 45 PetscCall(TSSetFromOptions(ts)); 46 47 PetscCall(TSSetStepNumber(ts, 0)); 48 PetscCall(TSSetTimeStep(ts, 1)); 49 PetscCall(TSSetTime(ts, 0)); 50 PetscCall(TSSetMaxTime(ts, PETSC_MAX_REAL)); 51 PetscCall(TSSetMaxSteps(ts, 3)); 52 53 /* 54 When an ARKIMEX scheme with an explicit stage is used this will error with a message informing the user it is not possible to use 55 a non-trivial mass matrix with ARKIMEX schemes with explicit stages. 56 */ 57 PetscCall(TSSolve(ts, NULL)); 58 59 PetscCall(TSDestroy(&ts)); 60 PetscCall(PetscFinalize()); 61 return 0; 62 } 63 64 PetscErrorCode IFunction(TS ts, PetscReal t, Vec x, Vec xdot, Vec f, void *ctx) { 65 PetscFunctionBeginUser; 66 PetscCall(VecCopy(xdot, f)); 67 PetscCall(VecScale(f, 2.0)); 68 PetscCall(VecShift(f, -1.0)); 69 PetscFunctionReturn(0); 70 } 71 72 PetscErrorCode IJacobian(TS ts, PetscReal t, Vec x, Vec xdot, PetscReal shift, Mat A, Mat B, void *ctx) { 73 PetscScalar j; 74 75 PetscFunctionBeginUser; 76 j = shift * 2.0; 77 PetscCall(MatSetValue(B, 0, 0, j, INSERT_VALUES)); 78 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 79 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 80 PetscFunctionReturn(0); 81 } 82 83 /*TEST 84 85 test: 86 suffix: arkimex_explicit_stage 87 requires: !defined(PETSCTEST_VALGRIND) defined(PETSC_USE_DEBUG) 88 args: -ts_type arkimex -petsc_ci_portable_error_output -error_output_stdout 89 filter: egrep -v "(options_left|memory block|leaked context|is not freed before MPI_Finalize|Could be the program crashed)" 90 91 test: 92 suffix: arkimex_implicit_stage 93 args: -ts_type arkimex -ts_arkimex_type l2 -ts_monitor_solution -ts_monitor 94 95 TEST*/ 96