xref: /petsc/src/ts/tests/ex7.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 except the user provided functions take input values as a single vector instead of two vectors
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 UV,Vec F)
18*c4762a1bSJed Brown {
19*c4762a1bSJed Brown   PetscErrorCode    ierr;
20*c4762a1bSJed Brown   const PetscScalar *u,*v;
21*c4762a1bSJed Brown   PetscScalar       *f;
22*c4762a1bSJed Brown   PetscInt          n,i;
23*c4762a1bSJed Brown 
24*c4762a1bSJed Brown   PetscFunctionBeginUser;
25*c4762a1bSJed Brown   ierr = VecGetLocalSize(UV,&n);CHKERRQ(ierr);
26*c4762a1bSJed Brown   n    = n/2;
27*c4762a1bSJed Brown   ierr = VecGetArrayRead(UV,&u);CHKERRQ(ierr);
28*c4762a1bSJed Brown   v    = u + n;
29*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
30*c4762a1bSJed Brown   for (i=0; i<n; i++) f[i] = u[i] + v[i];
31*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(UV,&u);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
33*c4762a1bSJed Brown   PetscFunctionReturn(0);
34*c4762a1bSJed Brown }
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown /*
37*c4762a1bSJed Brown    F(U,V) = U - V
38*c4762a1bSJed Brown 
39*c4762a1bSJed Brown */
40*c4762a1bSJed Brown PetscErrorCode F(PetscReal t,Vec UV,Vec F)
41*c4762a1bSJed Brown {
42*c4762a1bSJed Brown   PetscErrorCode    ierr;
43*c4762a1bSJed Brown   const PetscScalar *u,*v;
44*c4762a1bSJed Brown   PetscScalar       *f;
45*c4762a1bSJed Brown   PetscInt          n,i;
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown   PetscFunctionBeginUser;
48*c4762a1bSJed Brown   ierr = VecGetLocalSize(UV,&n);CHKERRQ(ierr);
49*c4762a1bSJed Brown   n    = n/2;
50*c4762a1bSJed Brown   ierr = VecGetArrayRead(UV,&u);CHKERRQ(ierr);
51*c4762a1bSJed Brown   v    = u + n;
52*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
53*c4762a1bSJed Brown   for (i=0; i<n; i++) f[i] = u[i] - v[i];
54*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(UV,&u);CHKERRQ(ierr);
55*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
56*c4762a1bSJed Brown   PetscFunctionReturn(0);
57*c4762a1bSJed Brown }
58*c4762a1bSJed Brown 
59*c4762a1bSJed Brown typedef struct {
60*c4762a1bSJed Brown   PetscReal      t;
61*c4762a1bSJed Brown   SNES           snes;
62*c4762a1bSJed Brown   Vec            UV,V;
63*c4762a1bSJed Brown   VecScatter     scatterU,scatterV;
64*c4762a1bSJed Brown   PetscErrorCode (*f)(PetscReal,Vec,Vec);
65*c4762a1bSJed Brown   PetscErrorCode (*F)(PetscReal,Vec,Vec);
66*c4762a1bSJed Brown } AppCtx;
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown extern PetscErrorCode TSFunction(TS,PetscReal,Vec,Vec,void*);
69*c4762a1bSJed Brown extern PetscErrorCode SNESFunction(SNES,Vec,Vec,void*);
70*c4762a1bSJed Brown 
71*c4762a1bSJed Brown int main(int argc,char **argv)
72*c4762a1bSJed Brown {
73*c4762a1bSJed Brown   PetscErrorCode ierr;
74*c4762a1bSJed Brown   AppCtx         ctx;
75*c4762a1bSJed Brown   TS             ts;
76*c4762a1bSJed Brown   Vec            tsrhs,U;
77*c4762a1bSJed Brown   IS             is;
78*c4762a1bSJed Brown   PetscInt       I;
79*c4762a1bSJed Brown   PetscMPIInt    rank;
80*c4762a1bSJed Brown 
81*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
82*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
84*c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
85*c4762a1bSJed Brown   ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr);
86*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&tsrhs);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&U);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,tsrhs,TSFunction,&ctx);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ctx.f = f;
91*c4762a1bSJed Brown 
92*c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD,&ctx.snes);CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr = SNESSetFromOptions(ctx.snes);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = SNESSetFunction(ctx.snes,NULL,SNESFunction,&ctx);CHKERRQ(ierr);
95*c4762a1bSJed Brown   ierr = SNESSetJacobian(ctx.snes,NULL,NULL,SNESComputeJacobianDefault,&ctx);CHKERRQ(ierr);
96*c4762a1bSJed Brown   ctx.F = F;
97*c4762a1bSJed Brown   ierr = VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.V);CHKERRQ(ierr);
98*c4762a1bSJed Brown 
99*c4762a1bSJed Brown   /* Create scatters to move between separate U and V representation and UV representation of solution */
100*c4762a1bSJed Brown   ierr = VecCreateMPI(PETSC_COMM_WORLD,2,PETSC_DETERMINE,&ctx.UV);CHKERRQ(ierr);
101*c4762a1bSJed Brown   I    = 2*rank;
102*c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_WORLD,1,&I,PETSC_COPY_VALUES,&is);CHKERRQ(ierr);
103*c4762a1bSJed Brown   ierr = VecScatterCreate(U,NULL,ctx.UV,is,&ctx.scatterU);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = ISDestroy(&is);CHKERRQ(ierr);
105*c4762a1bSJed Brown   I    = 2*rank + 1;
106*c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_WORLD,1,&I,PETSC_COPY_VALUES,&is);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = VecScatterCreate(ctx.V,NULL,ctx.UV,is,&ctx.scatterV);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = ISDestroy(&is);CHKERRQ(ierr);
109*c4762a1bSJed Brown 
110*c4762a1bSJed Brown   ierr = VecSet(U,1.0);CHKERRQ(ierr);
111*c4762a1bSJed Brown   ierr = TSSolve(ts,U);CHKERRQ(ierr);
112*c4762a1bSJed Brown 
113*c4762a1bSJed Brown   ierr = VecDestroy(&ctx.V);CHKERRQ(ierr);
114*c4762a1bSJed Brown   ierr = VecDestroy(&ctx.UV);CHKERRQ(ierr);
115*c4762a1bSJed Brown   ierr = VecScatterDestroy(&ctx.scatterU);CHKERRQ(ierr);
116*c4762a1bSJed Brown   ierr = VecScatterDestroy(&ctx.scatterV);CHKERRQ(ierr);
117*c4762a1bSJed Brown   ierr = VecDestroy(&tsrhs);CHKERRQ(ierr);
118*c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);
119*c4762a1bSJed Brown   ierr = SNESDestroy(&ctx.snes);CHKERRQ(ierr);
120*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = PetscFinalize();
122*c4762a1bSJed Brown   return ierr;
123*c4762a1bSJed Brown }
124*c4762a1bSJed Brown 
125*c4762a1bSJed Brown /*
126*c4762a1bSJed Brown    Defines the RHS function that is passed to the time-integrator.
127*c4762a1bSJed Brown 
128*c4762a1bSJed Brown    Solves F(U,V) for V and then computes f(U,V)
129*c4762a1bSJed Brown 
130*c4762a1bSJed Brown */
131*c4762a1bSJed Brown PetscErrorCode TSFunction(TS ts,PetscReal t,Vec U,Vec F,void *actx)
132*c4762a1bSJed Brown {
133*c4762a1bSJed Brown   AppCtx         *ctx = (AppCtx*)actx;
134*c4762a1bSJed Brown   PetscErrorCode ierr;
135*c4762a1bSJed Brown 
136*c4762a1bSJed Brown   PetscFunctionBeginUser;
137*c4762a1bSJed Brown   ctx->t = t;
138*c4762a1bSJed Brown   ierr = VecScatterBegin(ctx->scatterU,U,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
139*c4762a1bSJed Brown   ierr = VecScatterEnd(ctx->scatterU,U,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
140*c4762a1bSJed Brown   ierr = SNESSolve(ctx->snes,NULL,ctx->V);CHKERRQ(ierr);
141*c4762a1bSJed Brown   ierr = VecScatterBegin(ctx->scatterV,ctx->V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
142*c4762a1bSJed Brown   ierr = VecScatterEnd(ctx->scatterV,ctx->V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
143*c4762a1bSJed Brown   ierr = (*ctx->f)(t,ctx->UV,F);CHKERRQ(ierr);
144*c4762a1bSJed Brown   PetscFunctionReturn(0);
145*c4762a1bSJed Brown }
146*c4762a1bSJed Brown 
147*c4762a1bSJed Brown /*
148*c4762a1bSJed Brown    Defines the nonlinear function that is passed to the nonlinear solver
149*c4762a1bSJed Brown 
150*c4762a1bSJed Brown */
151*c4762a1bSJed Brown PetscErrorCode SNESFunction(SNES snes,Vec V,Vec F,void *actx)
152*c4762a1bSJed Brown {
153*c4762a1bSJed Brown   AppCtx         *ctx = (AppCtx*)actx;
154*c4762a1bSJed Brown   PetscErrorCode ierr;
155*c4762a1bSJed Brown 
156*c4762a1bSJed Brown   PetscFunctionBeginUser;
157*c4762a1bSJed Brown   ierr = VecScatterBegin(ctx->scatterV,V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
158*c4762a1bSJed Brown   ierr = VecScatterEnd(ctx->scatterV,V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
159*c4762a1bSJed Brown   ierr = (*ctx->F)(ctx->t,ctx->UV,F);CHKERRQ(ierr);
160*c4762a1bSJed Brown   PetscFunctionReturn(0);
161*c4762a1bSJed Brown }
162*c4762a1bSJed Brown 
163*c4762a1bSJed Brown 
164