xref: /petsc/src/ts/event/tests/ex16.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1fa9584fbSIlya Fursov static char help[] = "Solves the trivial ODE du/dt = 1, u(0) = 0. \n\n";
2fa9584fbSIlya Fursov /*
3fa9584fbSIlya Fursov   This example tests TSEvent's capability to handle complicated cases.
4fa9584fbSIlya Fursov   Test 1: an event close to endpoint of a time step should not be detected twice.
5fa9584fbSIlya Fursov   Test 2: two events in the same time step should be detected in the correct order.
6fa9584fbSIlya Fursov   Test 3: three events in the same time step should be detected in the correct order.
7fa9584fbSIlya Fursov */
8fa9584fbSIlya Fursov 
9fa9584fbSIlya Fursov #include <petscts.h>
10fa9584fbSIlya Fursov 
11fa9584fbSIlya Fursov static PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *);
12fa9584fbSIlya Fursov static PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *);
13fa9584fbSIlya Fursov 
14fa9584fbSIlya Fursov static PetscErrorCode Event(TS, PetscReal, Vec, PetscReal *, void *);
15fa9584fbSIlya Fursov static PetscErrorCode PostEvent(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *);
16fa9584fbSIlya Fursov 
main(int argc,char ** argv)17fa9584fbSIlya Fursov int main(int argc, char **argv)
18fa9584fbSIlya Fursov {
19fa9584fbSIlya Fursov   TS              ts;
20fa9584fbSIlya Fursov   const PetscReal t_end = 2.0;
21fa9584fbSIlya Fursov   Vec             x;
22fa9584fbSIlya Fursov   Vec             f;
23fa9584fbSIlya Fursov   Mat             A;
24fa9584fbSIlya Fursov 
25fa9584fbSIlya Fursov   PetscFunctionBeginUser;
26c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
27fa9584fbSIlya Fursov 
28fa9584fbSIlya Fursov   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
29fa9584fbSIlya Fursov 
30fa9584fbSIlya Fursov   PetscCall(VecCreate(PETSC_COMM_WORLD, &f));
31fa9584fbSIlya Fursov   PetscCall(VecSetSizes(f, 1, PETSC_DECIDE));
32fa9584fbSIlya Fursov   PetscCall(VecSetFromOptions(f));
33fa9584fbSIlya Fursov   PetscCall(VecSetUp(f));
34fa9584fbSIlya Fursov   PetscCall(TSSetRHSFunction(ts, f, RHSFunction, NULL));
35fa9584fbSIlya Fursov   PetscCall(VecDestroy(&f));
36fa9584fbSIlya Fursov 
37fa9584fbSIlya Fursov   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
38fa9584fbSIlya Fursov   PetscCall(MatSetSizes(A, 1, 1, PETSC_DECIDE, PETSC_DECIDE));
39fa9584fbSIlya Fursov   PetscCall(MatSetFromOptions(A));
40fa9584fbSIlya Fursov   PetscCall(MatSetUp(A));
41fa9584fbSIlya Fursov   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
42fa9584fbSIlya Fursov   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
43fa9584fbSIlya Fursov   /* ensure that the Jacobian matrix has diagonal entries since that is required by TS */
44fa9584fbSIlya Fursov   PetscCall(MatShift(A, (PetscReal)1));
45fa9584fbSIlya Fursov   PetscCall(MatShift(A, (PetscReal)-1));
46fa9584fbSIlya Fursov   PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, NULL));
47fa9584fbSIlya Fursov   PetscCall(MatDestroy(&A));
48fa9584fbSIlya Fursov 
49fa9584fbSIlya Fursov   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
50fa9584fbSIlya Fursov   PetscCall(VecSetSizes(x, 1, PETSC_DECIDE));
51fa9584fbSIlya Fursov   PetscCall(VecSetFromOptions(x));
52fa9584fbSIlya Fursov   PetscCall(VecSetUp(x));
53fa9584fbSIlya Fursov   PetscCall(TSSetSolution(ts, x));
54fa9584fbSIlya Fursov   PetscCall(VecDestroy(&x));
55fa9584fbSIlya Fursov 
56fa9584fbSIlya Fursov   {
57fa9584fbSIlya Fursov     PetscInt  direction[3];
58fa9584fbSIlya Fursov     PetscBool terminate[3];
59fa9584fbSIlya Fursov     direction[0] = +1;
60fa9584fbSIlya Fursov     terminate[0] = PETSC_FALSE;
61fa9584fbSIlya Fursov     direction[1] = -1;
62fa9584fbSIlya Fursov     terminate[1] = PETSC_FALSE;
63fa9584fbSIlya Fursov     direction[2] = 0;
64fa9584fbSIlya Fursov     terminate[2] = PETSC_FALSE;
65fa9584fbSIlya Fursov     PetscCall(TSSetEventHandler(ts, 3, direction, terminate, Event, PostEvent, NULL));
66fa9584fbSIlya Fursov   }
67fa9584fbSIlya Fursov   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
68fa9584fbSIlya Fursov   PetscCall(TSSetMaxTime(ts, t_end));
69fa9584fbSIlya Fursov   PetscCall(TSSetFromOptions(ts));
70fa9584fbSIlya Fursov 
71fa9584fbSIlya Fursov   PetscCall(TSSolve(ts, NULL));
72fa9584fbSIlya Fursov 
73fa9584fbSIlya Fursov   PetscCall(TSDestroy(&ts));
74fa9584fbSIlya Fursov   PetscCall(PetscFinalize());
75fa9584fbSIlya Fursov   return 0;
76fa9584fbSIlya Fursov }
77fa9584fbSIlya Fursov 
RHSFunction(TS ts,PetscReal t,Vec x,Vec f,PetscCtx ctx)78*2a8381b2SBarry Smith PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec x, Vec f, PetscCtx ctx)
79fa9584fbSIlya Fursov {
80fa9584fbSIlya Fursov   PetscFunctionBeginUser;
81fa9584fbSIlya Fursov   PetscCall(VecSet(f, (PetscReal)1));
82fa9584fbSIlya Fursov   PetscFunctionReturn(PETSC_SUCCESS);
83fa9584fbSIlya Fursov }
84fa9584fbSIlya Fursov 
RHSJacobian(TS ts,PetscReal t,Vec x,Mat A,Mat B,PetscCtx ctx)85*2a8381b2SBarry Smith PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec x, Mat A, Mat B, PetscCtx ctx)
86fa9584fbSIlya Fursov {
87fa9584fbSIlya Fursov   PetscFunctionBeginUser;
88fa9584fbSIlya Fursov   PetscCall(MatZeroEntries(B));
89fa9584fbSIlya Fursov   if (B != A) PetscCall(MatZeroEntries(A));
90fa9584fbSIlya Fursov   PetscFunctionReturn(PETSC_SUCCESS);
91fa9584fbSIlya Fursov }
92fa9584fbSIlya Fursov 
Event(TS ts,PetscReal t,Vec x,PetscReal * fvalue,PetscCtx ctx)93*2a8381b2SBarry Smith PetscErrorCode Event(TS ts, PetscReal t, Vec x, PetscReal *fvalue, PetscCtx ctx)
94fa9584fbSIlya Fursov {
95fa9584fbSIlya Fursov   PetscFunctionBeginUser;
96fa9584fbSIlya Fursov   fvalue[0] = t - 1.1;
97fa9584fbSIlya Fursov   fvalue[1] = 1.2 - t;
98fa9584fbSIlya Fursov   fvalue[2] = t - 1.3;
99fa9584fbSIlya Fursov   PetscFunctionReturn(PETSC_SUCCESS);
100fa9584fbSIlya Fursov }
101fa9584fbSIlya Fursov 
PostEvent(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec x,PetscBool forwardsolve,PetscCtx ctx)102*2a8381b2SBarry Smith PetscErrorCode PostEvent(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec x, PetscBool forwardsolve, PetscCtx ctx)
103fa9584fbSIlya Fursov {
104fa9584fbSIlya Fursov   PetscInt           i;
105fa9584fbSIlya Fursov   const PetscScalar *a;
106fa9584fbSIlya Fursov 
107fa9584fbSIlya Fursov   PetscFunctionBeginUser;
108fa9584fbSIlya Fursov   PetscCall(TSGetStepNumber(ts, &i));
109fa9584fbSIlya Fursov   PetscCall(VecGetArrayRead(x, &a));
110fa9584fbSIlya Fursov   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "%-10s-> step %" PetscInt_FMT " time %g value %g\n", PETSC_FUNCTION_NAME, i, (double)t, (double)PetscRealPart(a[0])));
111fa9584fbSIlya Fursov   PetscCall(VecRestoreArrayRead(x, &a));
112fa9584fbSIlya Fursov   PetscFunctionReturn(PETSC_SUCCESS);
113fa9584fbSIlya Fursov }
114fa9584fbSIlya Fursov 
115fa9584fbSIlya Fursov /*TEST
116fa9584fbSIlya Fursov 
117fa9584fbSIlya Fursov     test:
118188af4bfSBarry Smith       args: -ts_type beuler -ts_time_step 0.1 -ts_event_monitor
119fa9584fbSIlya Fursov 
120fa9584fbSIlya Fursov     test:
121fa9584fbSIlya Fursov       suffix: 2
122188af4bfSBarry Smith       args: -ts_type beuler -ts_time_step 0.2 -ts_event_monitor
123fa9584fbSIlya Fursov 
124fa9584fbSIlya Fursov     test:
125fa9584fbSIlya Fursov       suffix: 3
126188af4bfSBarry Smith       args: -ts_type beuler -ts_time_step 0.5 -ts_event_monitor
127fa9584fbSIlya Fursov TEST*/
128