xref: /petsc/src/ts/tests/ex26.c (revision 76be6f4ff3bd4e251c19fc00ebbebfd58b6e7589)
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 {
11   TS              ts;
12   Vec             x;
13   Vec             f;
14   Mat             A;
15   PetscErrorCode  ierr;
16 
17   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
18 
19   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
20   PetscCall(TSSetEquationType(ts,TS_EQ_ODE_IMPLICIT));
21   PetscCall(VecCreate(PETSC_COMM_WORLD,&f));
22   PetscCall(VecSetSizes(f,1,PETSC_DECIDE));
23   PetscCall(VecSetFromOptions(f));
24   PetscCall(VecSetUp(f));
25   PetscCall(TSSetIFunction(ts,f,IFunction,NULL));
26   PetscCall(VecDestroy(&f));
27 
28   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
29   PetscCall(MatSetSizes(A,1,1,PETSC_DECIDE,PETSC_DECIDE));
30   PetscCall(MatSetFromOptions(A));
31   PetscCall(MatSetUp(A));
32   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
33   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
34   /* ensure that the Jacobian matrix has diagonal entries since that is required by TS */
35   PetscCall(MatShift(A,(PetscReal)1));
36   PetscCall(MatShift(A,(PetscReal)-1));
37   PetscCall(TSSetIJacobian(ts,A,A,IJacobian,NULL));
38   PetscCall(MatDestroy(&A));
39 
40   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
41   PetscCall(VecSetSizes(x,1,PETSC_DECIDE));
42   PetscCall(VecSetFromOptions(x));
43   PetscCall(VecSetUp(x));
44   PetscCall(TSSetSolution(ts,x));
45   PetscCall(VecDestroy(&x));
46   PetscCall(TSSetFromOptions(ts));
47 
48   PetscCall(TSSetStepNumber(ts,0));
49   PetscCall(TSSetTimeStep(ts,1));
50   PetscCall(TSSetTime(ts,0));
51   PetscCall(TSSetMaxTime(ts,PETSC_MAX_REAL));
52   PetscCall(TSSetMaxSteps(ts,3));
53 
54   /*
55       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
56       a non-trivial mass matrix with ARKIMEX schemes with explicit stages.
57   */
58   ierr = TSSolve(ts,NULL);
59   if (ierr != PETSC_ERR_ARG_INCOMP) PetscCall(ierr);
60 
61   PetscCall(TSDestroy(&ts));
62   PetscCall(PetscFinalize());
63   return 0;
64 }
65 
66 PetscErrorCode IFunction(TS ts,PetscReal t,Vec x,Vec xdot,Vec f,void *ctx)
67 {
68   PetscFunctionBegin;
69   PetscCall(VecCopy(xdot,f));
70   PetscCall(VecScale(f,2.0));
71   PetscCall(VecShift(f,-1.0));
72   PetscFunctionReturn(0);
73 }
74 
75 PetscErrorCode IJacobian(TS ts,PetscReal t,Vec x,Vec xdot,PetscReal shift,Mat A,Mat B,void *ctx)
76 {
77   PetscScalar    j;
78 
79   PetscFunctionBegin;
80   j = shift*2.0;
81   PetscCall(MatSetValue(B,0,0,j,INSERT_VALUES));
82   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
83   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
84   PetscFunctionReturn(0);
85 }
86 
87 /*TEST
88 
89     test:
90       suffix: arkimex_explicit_stage
91       requires: defined(PETSC_USE_DEBUG)
92       args: -ts_type arkimex -error_output_stdout
93       filter:  egrep -v "(Petsc|on a| at |Configure)"
94 
95     test:
96       suffix: arkimex_implicit_stage
97       args: -ts_type arkimex -ts_arkimex_type l2 -ts_monitor_solution -ts_monitor
98 
99 TEST*/
100