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