xref: /petsc/src/ts/tests/ex7.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 except the user provided functions take input values as a single vector instead of two vectors
10c4762a1bSJed Brown */
11c4762a1bSJed Brown 
12c4762a1bSJed Brown /*
13c4762a1bSJed Brown    f(U,V) = U + V
14c4762a1bSJed Brown 
15c4762a1bSJed Brown */
16c4762a1bSJed Brown PetscErrorCode f(PetscReal t,Vec UV,Vec F)
17c4762a1bSJed Brown {
18c4762a1bSJed Brown   const PetscScalar *u,*v;
19c4762a1bSJed Brown   PetscScalar       *f;
20c4762a1bSJed Brown   PetscInt          n,i;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   PetscFunctionBeginUser;
235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(UV,&n));
24c4762a1bSJed Brown   n    = n/2;
255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(UV,&u));
26c4762a1bSJed Brown   v    = u + n;
275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayWrite(F,&f));
28c4762a1bSJed Brown   for (i=0; i<n; i++) f[i] = u[i] + v[i];
295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(UV,&u));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayWrite(F,&f));
31c4762a1bSJed Brown   PetscFunctionReturn(0);
32c4762a1bSJed Brown }
33c4762a1bSJed Brown 
34c4762a1bSJed Brown /*
35c4762a1bSJed Brown    F(U,V) = U - V
36c4762a1bSJed Brown */
37c4762a1bSJed Brown PetscErrorCode F(PetscReal t,Vec UV,Vec F)
38c4762a1bSJed Brown {
39c4762a1bSJed Brown   const PetscScalar *u,*v;
40c4762a1bSJed Brown   PetscScalar       *f;
41c4762a1bSJed Brown   PetscInt          n,i;
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   PetscFunctionBeginUser;
445f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(UV,&n));
45c4762a1bSJed Brown   n    = n/2;
465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(UV,&u));
47c4762a1bSJed Brown   v    = u + n;
485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayWrite(F,&f));
49c4762a1bSJed Brown   for (i=0; i<n; i++) f[i] = u[i] - v[i];
505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(UV,&u));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayWrite(F,&f));
52c4762a1bSJed Brown   PetscFunctionReturn(0);
53c4762a1bSJed Brown }
54c4762a1bSJed Brown 
55c4762a1bSJed Brown typedef struct {
56c4762a1bSJed Brown   PetscReal      t;
57c4762a1bSJed Brown   SNES           snes;
58c4762a1bSJed Brown   Vec            UV,V;
59c4762a1bSJed Brown   VecScatter     scatterU,scatterV;
60c4762a1bSJed Brown   PetscErrorCode (*f)(PetscReal,Vec,Vec);
61c4762a1bSJed Brown   PetscErrorCode (*F)(PetscReal,Vec,Vec);
62c4762a1bSJed Brown } AppCtx;
63c4762a1bSJed Brown 
64c4762a1bSJed Brown extern PetscErrorCode TSFunction(TS,PetscReal,Vec,Vec,void*);
65c4762a1bSJed Brown extern PetscErrorCode SNESFunction(SNES,Vec,Vec,void*);
66c4762a1bSJed Brown 
67c4762a1bSJed Brown int main(int argc,char **argv)
68c4762a1bSJed Brown {
69c4762a1bSJed Brown   AppCtx         ctx;
70c4762a1bSJed Brown   TS             ts;
71c4762a1bSJed Brown   Vec            tsrhs,U;
72c4762a1bSJed Brown   IS             is;
73303a5415SBarry Smith   PetscInt       i;
74c4762a1bSJed Brown   PetscMPIInt    rank;
75c4762a1bSJed Brown 
76*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
775f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSEULER));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&tsrhs));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(tsrhs,&U));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,tsrhs,TSFunction,&ctx));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,1.0));
86c4762a1bSJed Brown   ctx.f = f;
87c4762a1bSJed Brown 
885f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&ctx.snes));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(ctx.snes));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFunction(ctx.snes,NULL,SNESFunction,&ctx));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetJacobian(ctx.snes,NULL,NULL,SNESComputeJacobianDefault,&ctx));
92c4762a1bSJed Brown   ctx.F = F;
935f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.V));
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* Create scatters to move between separate U and V representation and UV representation of solution */
965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,2,PETSC_DETERMINE,&ctx.UV));
97303a5415SBarry Smith   i    = 2*rank;
985f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,1,&i,PETSC_COPY_VALUES,&is));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(U,NULL,ctx.UV,is,&ctx.scatterU));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is));
101303a5415SBarry Smith   i    = 2*rank + 1;
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,1,&i,PETSC_COPY_VALUES,&is));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(ctx.V,NULL,ctx.UV,is,&ctx.scatterV));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is));
105c4762a1bSJed Brown 
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(U,1.0));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,U));
108c4762a1bSJed Brown 
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&ctx.V));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&ctx.UV));
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&ctx.scatterU));
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&ctx.scatterV));
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&tsrhs));
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&U));
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&ctx.snes));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
117*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
118*b122ec5aSJacob Faibussowitsch   return 0;
119c4762a1bSJed Brown }
120c4762a1bSJed Brown 
121c4762a1bSJed Brown /*
122c4762a1bSJed Brown    Defines the RHS function that is passed to the time-integrator.
123c4762a1bSJed Brown 
124c4762a1bSJed Brown    Solves F(U,V) for V and then computes f(U,V)
125c4762a1bSJed Brown */
126c4762a1bSJed Brown PetscErrorCode TSFunction(TS ts,PetscReal t,Vec U,Vec F,void *actx)
127c4762a1bSJed Brown {
128c4762a1bSJed Brown   AppCtx         *ctx = (AppCtx*)actx;
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   PetscFunctionBeginUser;
131c4762a1bSJed Brown   ctx->t = t;
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(ctx->scatterU,U,ctx->UV,INSERT_VALUES,SCATTER_FORWARD));
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(ctx->scatterU,U,ctx->UV,INSERT_VALUES,SCATTER_FORWARD));
1345f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(ctx->snes,NULL,ctx->V));
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(ctx->scatterV,ctx->V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD));
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(ctx->scatterV,ctx->V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD));
1375f80ce2aSJacob Faibussowitsch   CHKERRQ((*ctx->f)(t,ctx->UV,F));
138c4762a1bSJed Brown   PetscFunctionReturn(0);
139c4762a1bSJed Brown }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown /*
142c4762a1bSJed Brown    Defines the nonlinear function that is passed to the nonlinear solver
143c4762a1bSJed Brown 
144c4762a1bSJed Brown */
145c4762a1bSJed Brown PetscErrorCode SNESFunction(SNES snes,Vec V,Vec F,void *actx)
146c4762a1bSJed Brown {
147c4762a1bSJed Brown   AppCtx         *ctx = (AppCtx*)actx;
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   PetscFunctionBeginUser;
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(ctx->scatterV,V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD));
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(ctx->scatterV,V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ((*ctx->F)(ctx->t,ctx->UV,F));
153c4762a1bSJed Brown   PetscFunctionReturn(0);
154c4762a1bSJed Brown }
155c4762a1bSJed Brown 
156303a5415SBarry Smith /*TEST
157303a5415SBarry Smith 
158303a5415SBarry Smith     test:
159303a5415SBarry Smith       args: -ts_monitor -ts_view
160303a5415SBarry Smith 
161303a5415SBarry Smith TEST*/
162