xref: /petsc/src/ts/tests/ex9.c (revision 9566063d113dddea24716c546802770db7481bc0)
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 
9c4762a1bSJed Brown     Same as ex6.c and ex7.c except calls the ARKIMEX integrator on the entire DAE
10c4762a1bSJed Brown */
11c4762a1bSJed Brown 
12c4762a1bSJed Brown /*
13c4762a1bSJed Brown    f(U,V) = U + V
14c4762a1bSJed Brown 
15c4762a1bSJed Brown */
16c4762a1bSJed Brown PetscErrorCode f(PetscReal t,Vec U,Vec V,Vec F)
17c4762a1bSJed Brown {
18c4762a1bSJed Brown   PetscErrorCode ierr;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   PetscFunctionBeginUser;
21*9566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(F,1.0,U,V));
22c4762a1bSJed Brown   PetscFunctionReturn(0);
23c4762a1bSJed Brown }
24c4762a1bSJed Brown 
25c4762a1bSJed Brown /*
26c4762a1bSJed Brown    F(U,V) = U - V
27c4762a1bSJed Brown 
28c4762a1bSJed Brown */
29c4762a1bSJed Brown PetscErrorCode F(PetscReal t,Vec U,Vec V,Vec F)
30c4762a1bSJed Brown {
31c4762a1bSJed Brown   PetscErrorCode ierr;
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   PetscFunctionBeginUser;
34*9566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(F,-1.0,V,U));
35c4762a1bSJed Brown   PetscFunctionReturn(0);
36c4762a1bSJed Brown }
37c4762a1bSJed Brown 
38c4762a1bSJed Brown typedef struct {
39c4762a1bSJed Brown   Vec            U,V;
40c4762a1bSJed Brown   Vec            UF,VF;
41c4762a1bSJed Brown   VecScatter     scatterU,scatterV;
42c4762a1bSJed Brown   PetscErrorCode (*f)(PetscReal,Vec,Vec,Vec);
43c4762a1bSJed Brown   PetscErrorCode (*F)(PetscReal,Vec,Vec,Vec);
44c4762a1bSJed Brown } AppCtx;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown extern PetscErrorCode TSFunctionRHS(TS,PetscReal,Vec,Vec,void*);
47c4762a1bSJed Brown extern PetscErrorCode TSFunctionI(TS,PetscReal,Vec,Vec,Vec,void*);
48c4762a1bSJed Brown 
49c4762a1bSJed Brown int main(int argc,char **argv)
50c4762a1bSJed Brown {
51c4762a1bSJed Brown   PetscErrorCode ierr;
52c4762a1bSJed Brown   AppCtx         ctx;
53c4762a1bSJed Brown   TS             ts;
54c4762a1bSJed Brown   Vec            tsrhs,UV;
55c4762a1bSJed Brown   IS             is;
56c4762a1bSJed Brown   PetscInt       I;
57c4762a1bSJed Brown   PetscMPIInt    rank;
58c4762a1bSJed Brown 
59*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
60*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
61*9566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
62*9566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
63*9566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSROSW));
64*9566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
65*9566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PETSC_COMM_WORLD,2,PETSC_DETERMINE,&tsrhs));
66*9566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PETSC_COMM_WORLD,2,PETSC_DETERMINE,&UV));
67*9566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts,tsrhs,TSFunctionRHS,&ctx));
68*9566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts,NULL,TSFunctionI,&ctx));
69c4762a1bSJed Brown   ctx.f = f;
70c4762a1bSJed Brown   ctx.F = F;
71c4762a1bSJed Brown 
72*9566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.U));
73*9566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.V));
74*9566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.UF));
75*9566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.VF));
76c4762a1bSJed Brown   I    = 2*rank;
77*9566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,1,&I,PETSC_COPY_VALUES,&is));
78*9566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(ctx.U,NULL,UV,is,&ctx.scatterU));
79*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
80c4762a1bSJed Brown   I    = 2*rank + 1;
81*9566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,1,&I,PETSC_COPY_VALUES,&is));
82*9566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(ctx.V,NULL,UV,is,&ctx.scatterV));
83*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
84c4762a1bSJed Brown 
85*9566063dSJacob Faibussowitsch   PetscCall(VecSet(UV,1.0));
86*9566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,UV));
87*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tsrhs));
88*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&UV));
89*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx.U));
90*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx.V));
91*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx.UF));
92*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx.VF));
93*9566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&ctx.scatterU));
94*9566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&ctx.scatterV));
95*9566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
96*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
97b122ec5aSJacob Faibussowitsch   return 0;
98c4762a1bSJed Brown }
99c4762a1bSJed Brown 
100c4762a1bSJed Brown /*
101c4762a1bSJed Brown    Defines the RHS function that is passed to the time-integrator.
102c4762a1bSJed Brown 
103c4762a1bSJed Brown */
104c4762a1bSJed Brown PetscErrorCode TSFunctionRHS(TS ts,PetscReal t,Vec UV,Vec F,void *actx)
105c4762a1bSJed Brown {
106c4762a1bSJed Brown   AppCtx         *ctx = (AppCtx*)actx;
107c4762a1bSJed Brown   PetscErrorCode ierr;
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   PetscFunctionBeginUser;
110*9566063dSJacob Faibussowitsch   PetscCall(VecSet(F,0.0));
111*9566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ctx->scatterU,UV,ctx->U,INSERT_VALUES,SCATTER_REVERSE));
112*9566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ctx->scatterU,UV,ctx->U,INSERT_VALUES,SCATTER_REVERSE));
113*9566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ctx->scatterV,UV,ctx->V,INSERT_VALUES,SCATTER_REVERSE));
114*9566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ctx->scatterV,UV,ctx->V,INSERT_VALUES,SCATTER_REVERSE));
115*9566063dSJacob Faibussowitsch   PetscCall((*ctx->f)(t,ctx->U,ctx->V,ctx->UF));
116*9566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ctx->scatterU,ctx->UF,F,INSERT_VALUES,SCATTER_FORWARD));
117*9566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ctx->scatterU,ctx->UF,F,INSERT_VALUES,SCATTER_FORWARD));
118c4762a1bSJed Brown   PetscFunctionReturn(0);
119c4762a1bSJed Brown }
120c4762a1bSJed Brown 
121c4762a1bSJed Brown /*
122c4762a1bSJed Brown    Defines the nonlinear function that is passed to the time-integrator
123c4762a1bSJed Brown 
124c4762a1bSJed Brown */
125c4762a1bSJed Brown PetscErrorCode TSFunctionI(TS ts,PetscReal t,Vec UV,Vec UVdot,Vec F,void *actx)
126c4762a1bSJed Brown {
127c4762a1bSJed Brown   AppCtx         *ctx = (AppCtx*)actx;
128c4762a1bSJed Brown   PetscErrorCode ierr;
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   PetscFunctionBeginUser;
131*9566063dSJacob Faibussowitsch   PetscCall(VecCopy(UVdot,F));
132*9566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ctx->scatterU,UV,ctx->U,INSERT_VALUES,SCATTER_REVERSE));
133*9566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ctx->scatterU,UV,ctx->U,INSERT_VALUES,SCATTER_REVERSE));
134*9566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ctx->scatterV,UV,ctx->V,INSERT_VALUES,SCATTER_REVERSE));
135*9566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ctx->scatterV,UV,ctx->V,INSERT_VALUES,SCATTER_REVERSE));
136*9566063dSJacob Faibussowitsch   PetscCall((*ctx->F)(t,ctx->U,ctx->V,ctx->VF));
137*9566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ctx->scatterV,ctx->VF,F,INSERT_VALUES,SCATTER_FORWARD));
138*9566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ctx->scatterV,ctx->VF,F,INSERT_VALUES,SCATTER_FORWARD));
139c4762a1bSJed Brown   PetscFunctionReturn(0);
140c4762a1bSJed Brown }
141