xref: /petsc/src/ts/tests/ex8.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1c4762a1bSJed Brown static char help[] = "Solves DAE with integrator only on non-algebraic terms \n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscts.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /*
6c4762a1bSJed Brown         \dot{U} = f(U,V)
7c4762a1bSJed Brown         F(U,V)  = 0
8c4762a1bSJed Brown 
93d5a8a6aSBarry Smith     Same as ex6.c and ex7.c except calls the TSROSW integrator on the entire DAE
10c4762a1bSJed Brown */
11c4762a1bSJed Brown 
12c4762a1bSJed Brown /*
13c4762a1bSJed Brown    f(U,V) = U + V
14c4762a1bSJed Brown 
15c4762a1bSJed Brown */
f(PetscReal t,Vec UV,Vec F)16d71ae5a4SJacob Faibussowitsch PetscErrorCode f(PetscReal t, Vec UV, Vec F)
17d71ae5a4SJacob Faibussowitsch {
18c4762a1bSJed Brown   const PetscScalar *u, *v;
19c4762a1bSJed Brown   PetscScalar       *f;
20c4762a1bSJed Brown   PetscInt           n, i;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   PetscFunctionBeginUser;
239566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(UV, &n));
24c4762a1bSJed Brown   n = n / 2;
259566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(UV, &u));
26c4762a1bSJed Brown   v = u + n;
279566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(F, &f));
28c4762a1bSJed Brown   for (i = 0; i < n; i++) f[i] = u[i] + v[i];
299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(UV, &u));
309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(F, &f));
313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32c4762a1bSJed Brown }
33c4762a1bSJed Brown 
34c4762a1bSJed Brown /*
35c4762a1bSJed Brown    F(U,V) = U - V
36c4762a1bSJed Brown 
37c4762a1bSJed Brown */
F(PetscReal t,Vec UV,Vec F)38d71ae5a4SJacob Faibussowitsch PetscErrorCode F(PetscReal t, Vec UV, Vec F)
39d71ae5a4SJacob Faibussowitsch {
40c4762a1bSJed Brown   const PetscScalar *u, *v;
41c4762a1bSJed Brown   PetscScalar       *f;
42c4762a1bSJed Brown   PetscInt           n, i;
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   PetscFunctionBeginUser;
459566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(UV, &n));
46c4762a1bSJed Brown   n = n / 2;
479566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(UV, &u));
48c4762a1bSJed Brown   v = u + n;
499566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(F, &f));
50c4762a1bSJed Brown   f = f + n;
51c4762a1bSJed Brown   for (i = 0; i < n; i++) f[i] = u[i] - v[i];
52c4762a1bSJed Brown   f = f - n;
539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(UV, &u));
549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(F, &f));
553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
56c4762a1bSJed Brown }
57c4762a1bSJed Brown 
58c4762a1bSJed Brown typedef struct {
59c4762a1bSJed Brown   PetscErrorCode (*f)(PetscReal, Vec, Vec);
60c4762a1bSJed Brown   PetscErrorCode (*F)(PetscReal, Vec, Vec);
61c4762a1bSJed Brown } AppCtx;
62c4762a1bSJed Brown 
63c4762a1bSJed Brown extern PetscErrorCode TSFunctionRHS(TS, PetscReal, Vec, Vec, void *);
64c4762a1bSJed Brown extern PetscErrorCode TSFunctionI(TS, PetscReal, Vec, Vec, Vec, void *);
65c4762a1bSJed Brown 
main(int argc,char ** argv)66d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
67d71ae5a4SJacob Faibussowitsch {
68c4762a1bSJed Brown   AppCtx ctx;
69c4762a1bSJed Brown   TS     ts;
70c4762a1bSJed Brown   Vec    tsrhs, UV;
71c4762a1bSJed Brown 
72327415f7SBarry Smith   PetscFunctionBeginUser;
73*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
749566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
759566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
769566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSROSW));
779566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
789566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PETSC_COMM_WORLD, 2, PETSC_DETERMINE, &tsrhs));
799566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tsrhs, &UV));
809566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, tsrhs, TSFunctionRHS, &ctx));
819566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, TSFunctionI, &ctx));
829566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 1.0));
83c4762a1bSJed Brown   ctx.f = f;
84c4762a1bSJed Brown   ctx.F = F;
85c4762a1bSJed Brown 
869566063dSJacob Faibussowitsch   PetscCall(VecSet(UV, 1.0));
879566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, UV));
889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tsrhs));
899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&UV));
909566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
919566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
92b122ec5aSJacob Faibussowitsch   return 0;
93c4762a1bSJed Brown }
94c4762a1bSJed Brown 
95c4762a1bSJed Brown /*
96c4762a1bSJed Brown    Defines the RHS function that is passed to the time-integrator.
97c4762a1bSJed Brown */
TSFunctionRHS(TS ts,PetscReal t,Vec UV,Vec F,void * actx)98d71ae5a4SJacob Faibussowitsch PetscErrorCode TSFunctionRHS(TS ts, PetscReal t, Vec UV, Vec F, void *actx)
99d71ae5a4SJacob Faibussowitsch {
100c4762a1bSJed Brown   AppCtx *ctx = (AppCtx *)actx;
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   PetscFunctionBeginUser;
1039566063dSJacob Faibussowitsch   PetscCall(VecSet(F, 0.0));
1049566063dSJacob Faibussowitsch   PetscCall((*ctx->f)(t, UV, F));
1053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
106c4762a1bSJed Brown }
107c4762a1bSJed Brown 
108c4762a1bSJed Brown /*
109c4762a1bSJed Brown    Defines the nonlinear function that is passed to the time-integrator
110c4762a1bSJed Brown */
TSFunctionI(TS ts,PetscReal t,Vec UV,Vec UVdot,Vec F,void * actx)111d71ae5a4SJacob Faibussowitsch PetscErrorCode TSFunctionI(TS ts, PetscReal t, Vec UV, Vec UVdot, Vec F, void *actx)
112d71ae5a4SJacob Faibussowitsch {
113c4762a1bSJed Brown   AppCtx *ctx = (AppCtx *)actx;
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   PetscFunctionBeginUser;
1169566063dSJacob Faibussowitsch   PetscCall(VecCopy(UVdot, F));
1179566063dSJacob Faibussowitsch   PetscCall((*ctx->F)(t, UV, F));
1183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
119c4762a1bSJed Brown }
120f7f07198SBarry Smith 
121f7f07198SBarry Smith /*TEST
122f7f07198SBarry Smith 
123f7f07198SBarry Smith     test:
124f7f07198SBarry Smith       args: -ts_view
125f7f07198SBarry Smith 
126f7f07198SBarry Smith     test:
127f7f07198SBarry Smith       suffix: 2
128f7f07198SBarry Smith       args: -snes_lag_jacobian 2 -ts_view
129f7f07198SBarry Smith 
130f7f07198SBarry Smith TEST*/
131