xref: /petsc/src/ts/tests/ex8.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 static char help[] = "Solves DAE with integrator only on non-algebraic terms \n";
2 
3 #include <petscts.h>
4 
5 /*
6         \dot{U} = f(U,V)
7         F(U,V)  = 0
8 
9     Same as ex6.c and ex7.c except calls the TSROSW integrator on the entire DAE
10 */
11 
12 /*
13    f(U,V) = U + V
14 
15 */
16 PetscErrorCode f(PetscReal t,Vec UV,Vec F)
17 {
18   const PetscScalar *u,*v;
19   PetscScalar       *f;
20   PetscInt          n,i;
21 
22   PetscFunctionBeginUser;
23   CHKERRQ(VecGetLocalSize(UV,&n));
24   n    = n/2;
25   CHKERRQ(VecGetArrayRead(UV,&u));
26   v    = u + n;
27   CHKERRQ(VecGetArrayWrite(F,&f));
28   for (i=0; i<n; i++) f[i] = u[i] + v[i];
29   CHKERRQ(VecRestoreArrayRead(UV,&u));
30   CHKERRQ(VecRestoreArrayWrite(F,&f));
31   PetscFunctionReturn(0);
32 }
33 
34 /*
35    F(U,V) = U - V
36 
37 */
38 PetscErrorCode F(PetscReal t,Vec UV,Vec F)
39 {
40   const PetscScalar *u,*v;
41   PetscScalar       *f;
42   PetscInt          n,i;
43 
44   PetscFunctionBeginUser;
45   CHKERRQ(VecGetLocalSize(UV,&n));
46   n    = n/2;
47   CHKERRQ(VecGetArrayRead(UV,&u));
48   v    = u + n;
49   CHKERRQ(VecGetArrayWrite(F,&f));
50   f    = f + n;
51   for (i=0; i<n; i++) f[i] = u[i] - v[i];
52   f    = f - n;
53   CHKERRQ(VecRestoreArrayRead(UV,&u));
54   CHKERRQ(VecRestoreArrayWrite(F,&f));
55   PetscFunctionReturn(0);
56 }
57 
58 typedef struct {
59   PetscErrorCode (*f)(PetscReal,Vec,Vec);
60   PetscErrorCode (*F)(PetscReal,Vec,Vec);
61 } AppCtx;
62 
63 extern PetscErrorCode TSFunctionRHS(TS,PetscReal,Vec,Vec,void*);
64 extern PetscErrorCode TSFunctionI(TS,PetscReal,Vec,Vec,Vec,void*);
65 
66 int main(int argc,char **argv)
67 {
68   PetscErrorCode ierr;
69   AppCtx         ctx;
70   TS             ts;
71   Vec            tsrhs,UV;
72 
73   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
74   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
75   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
76   CHKERRQ(TSSetType(ts,TSROSW));
77   CHKERRQ(TSSetFromOptions(ts));
78   CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,2,PETSC_DETERMINE,&tsrhs));
79   CHKERRQ(VecDuplicate(tsrhs,&UV));
80   CHKERRQ(TSSetRHSFunction(ts,tsrhs,TSFunctionRHS,&ctx));
81   CHKERRQ(TSSetIFunction(ts,NULL,TSFunctionI,&ctx));
82   CHKERRQ(TSSetMaxTime(ts,1.0));
83   ctx.f = f;
84   ctx.F = F;
85 
86   CHKERRQ(VecSet(UV,1.0));
87   CHKERRQ(TSSolve(ts,UV));
88   CHKERRQ(VecDestroy(&tsrhs));
89   CHKERRQ(VecDestroy(&UV));
90   CHKERRQ(TSDestroy(&ts));
91   ierr = PetscFinalize();
92   return ierr;
93 }
94 
95 /*
96    Defines the RHS function that is passed to the time-integrator.
97 */
98 PetscErrorCode TSFunctionRHS(TS ts,PetscReal t,Vec UV,Vec F,void *actx)
99 {
100   AppCtx         *ctx = (AppCtx*)actx;
101 
102   PetscFunctionBeginUser;
103   CHKERRQ(VecSet(F,0.0));
104   CHKERRQ((*ctx->f)(t,UV,F));
105   PetscFunctionReturn(0);
106 }
107 
108 /*
109    Defines the nonlinear function that is passed to the time-integrator
110 */
111 PetscErrorCode TSFunctionI(TS ts,PetscReal t,Vec UV,Vec UVdot,Vec F,void *actx)
112 {
113   AppCtx         *ctx = (AppCtx*)actx;
114 
115   PetscFunctionBeginUser;
116   CHKERRQ(VecCopy(UVdot,F));
117   CHKERRQ((*ctx->F)(t,UV,F));
118   PetscFunctionReturn(0);
119 }
120 
121 /*TEST
122 
123     test:
124       args:  -ts_view
125 
126     test:
127       suffix: 2
128       args: -snes_lag_jacobian 2 -ts_view
129 
130 TEST*/
131