xref: /petsc/src/ts/tests/ex26.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
16aad120cSJose E. Roman static char help[] = "Solves the trivial ODE 2 du/dt = 1, u(0) = 0. \n\n";
2c79dcfbdSBarry Smith 
3c79dcfbdSBarry Smith #include <petscts.h>
4c79dcfbdSBarry Smith #include <petscpc.h>
5c79dcfbdSBarry Smith 
691de4d3cSBarry Smith PetscErrorCode IFunction(TS, PetscReal, Vec, Vec, Vec, void *);
791de4d3cSBarry Smith PetscErrorCode IJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
8c79dcfbdSBarry Smith 
9d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
10d71ae5a4SJacob Faibussowitsch {
11c79dcfbdSBarry Smith   TS  ts;
12c79dcfbdSBarry Smith   Vec x;
13c79dcfbdSBarry Smith   Vec f;
14c79dcfbdSBarry Smith   Mat A;
15c79dcfbdSBarry Smith 
16327415f7SBarry Smith   PetscFunctionBeginUser;
179566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
18c79dcfbdSBarry Smith 
199566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
209566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(ts, TS_EQ_ODE_IMPLICIT));
219566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &f));
229566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(f, 1, PETSC_DECIDE));
239566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(f));
249566063dSJacob Faibussowitsch   PetscCall(VecSetUp(f));
259566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, f, IFunction, NULL));
269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&f));
27c79dcfbdSBarry Smith 
289566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
299566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, 1, 1, PETSC_DECIDE, PETSC_DECIDE));
309566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
319566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
339566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
34c79dcfbdSBarry Smith   /* ensure that the Jacobian matrix has diagonal entries since that is required by TS */
359566063dSJacob Faibussowitsch   PetscCall(MatShift(A, (PetscReal)1));
369566063dSJacob Faibussowitsch   PetscCall(MatShift(A, (PetscReal)-1));
379566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, A, A, IJacobian, NULL));
389566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
39c79dcfbdSBarry Smith 
409566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
419566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, 1, PETSC_DECIDE));
429566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
439566063dSJacob Faibussowitsch   PetscCall(VecSetUp(x));
449566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, x));
459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
469566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
47c79dcfbdSBarry Smith 
489566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
499566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 1));
509566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0));
519566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, PETSC_MAX_REAL));
529566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 3));
53c79dcfbdSBarry Smith 
54c79dcfbdSBarry Smith   /*
55c79dcfbdSBarry Smith       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
56c79dcfbdSBarry Smith       a non-trivial mass matrix with ARKIMEX schemes with explicit stages.
57c79dcfbdSBarry Smith   */
58cab2b192SBarry Smith   PetscCall(TSSolve(ts, NULL));
59c79dcfbdSBarry Smith 
609566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
619566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
62b122ec5aSJacob Faibussowitsch   return 0;
63c79dcfbdSBarry Smith }
64c79dcfbdSBarry Smith 
65d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction(TS ts, PetscReal t, Vec x, Vec xdot, Vec f, void *ctx)
66d71ae5a4SJacob Faibussowitsch {
677510d9b0SBarry Smith   PetscFunctionBeginUser;
689566063dSJacob Faibussowitsch   PetscCall(VecCopy(xdot, f));
699566063dSJacob Faibussowitsch   PetscCall(VecScale(f, 2.0));
709566063dSJacob Faibussowitsch   PetscCall(VecShift(f, -1.0));
71*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
72c79dcfbdSBarry Smith }
73c79dcfbdSBarry Smith 
74d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian(TS ts, PetscReal t, Vec x, Vec xdot, PetscReal shift, Mat A, Mat B, void *ctx)
75d71ae5a4SJacob Faibussowitsch {
76c79dcfbdSBarry Smith   PetscScalar j;
77c79dcfbdSBarry Smith 
787510d9b0SBarry Smith   PetscFunctionBeginUser;
79c79dcfbdSBarry Smith   j = shift * 2.0;
809566063dSJacob Faibussowitsch   PetscCall(MatSetValue(B, 0, 0, j, INSERT_VALUES));
819566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
829566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
83*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
84c79dcfbdSBarry Smith }
85c79dcfbdSBarry Smith 
86c79dcfbdSBarry Smith /*TEST
87c79dcfbdSBarry Smith 
88c79dcfbdSBarry Smith     test:
89c79dcfbdSBarry Smith       suffix: arkimex_explicit_stage
90cab2b192SBarry Smith       requires: !defined(PETSCTEST_VALGRIND) defined(PETSC_USE_DEBUG)
91cab2b192SBarry Smith       args: -ts_type arkimex -petsc_ci_portable_error_output -error_output_stdout
92f53b81b6SPierre Jolivet       filter: grep -E -v "(options_left|memory block|leaked context|not freed before MPI_Finalize|Could be the program crashed)"
93c79dcfbdSBarry Smith 
94c79dcfbdSBarry Smith     test:
95c79dcfbdSBarry Smith       suffix: arkimex_implicit_stage
96c79dcfbdSBarry Smith       args: -ts_type arkimex -ts_arkimex_type l2 -ts_monitor_solution -ts_monitor
97c79dcfbdSBarry Smith 
98c79dcfbdSBarry Smith TEST*/
99