Actual source code: ex7.c
petsc-3.14.6 2021-03-30
1: static char help[] = "Solves DAE with integrator only on non-algebraic terms \n";
3: #include <petscts.h>
5: /*
6: \dot{U} = f(U,V)
7: F(U,V) = 0
9: Same as ex6.c except the user provided functions take input values as a single vector instead of two vectors
10: */
13: /*
14: f(U,V) = U + V
16: */
17: PetscErrorCode f(PetscReal t,Vec UV,Vec F)
18: {
19: PetscErrorCode ierr;
20: const PetscScalar *u,*v;
21: PetscScalar *f;
22: PetscInt n,i;
25: VecGetLocalSize(UV,&n);
26: n = n/2;
27: VecGetArrayRead(UV,&u);
28: v = u + n;
29: VecGetArrayWrite(F,&f);
30: for (i=0; i<n; i++) f[i] = u[i] + v[i];
31: VecRestoreArrayRead(UV,&u);
32: VecRestoreArrayWrite(F,&f);
33: return(0);
34: }
36: /*
37: F(U,V) = U - V
38: */
39: PetscErrorCode F(PetscReal t,Vec UV,Vec F)
40: {
41: PetscErrorCode ierr;
42: const PetscScalar *u,*v;
43: PetscScalar *f;
44: PetscInt n,i;
47: VecGetLocalSize(UV,&n);
48: n = n/2;
49: VecGetArrayRead(UV,&u);
50: v = u + n;
51: VecGetArrayWrite(F,&f);
52: for (i=0; i<n; i++) f[i] = u[i] - v[i];
53: VecRestoreArrayRead(UV,&u);
54: VecRestoreArrayWrite(F,&f);
55: return(0);
56: }
58: typedef struct {
59: PetscReal t;
60: SNES snes;
61: Vec UV,V;
62: VecScatter scatterU,scatterV;
63: PetscErrorCode (*f)(PetscReal,Vec,Vec);
64: PetscErrorCode (*F)(PetscReal,Vec,Vec);
65: } AppCtx;
67: extern PetscErrorCode TSFunction(TS,PetscReal,Vec,Vec,void*);
68: extern PetscErrorCode SNESFunction(SNES,Vec,Vec,void*);
70: int main(int argc,char **argv)
71: {
73: AppCtx ctx;
74: TS ts;
75: Vec tsrhs,U;
76: IS is;
77: PetscInt i;
78: PetscMPIInt rank;
80: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
81: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
82: TSCreate(PETSC_COMM_WORLD,&ts);
83: TSSetProblemType(ts,TS_NONLINEAR);
84: TSSetType(ts,TSEULER);
85: TSSetFromOptions(ts);
86: VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&tsrhs);
87: VecDuplicate(tsrhs,&U);
88: TSSetRHSFunction(ts,tsrhs,TSFunction,&ctx);
89: TSSetMaxTime(ts,1.0);
90: ctx.f = f;
92: SNESCreate(PETSC_COMM_WORLD,&ctx.snes);
93: SNESSetFromOptions(ctx.snes);
94: SNESSetFunction(ctx.snes,NULL,SNESFunction,&ctx);
95: SNESSetJacobian(ctx.snes,NULL,NULL,SNESComputeJacobianDefault,&ctx);
96: ctx.F = F;
97: VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.V);
99: /* Create scatters to move between separate U and V representation and UV representation of solution */
100: VecCreateMPI(PETSC_COMM_WORLD,2,PETSC_DETERMINE,&ctx.UV);
101: i = 2*rank;
102: ISCreateGeneral(PETSC_COMM_WORLD,1,&i,PETSC_COPY_VALUES,&is);
103: VecScatterCreate(U,NULL,ctx.UV,is,&ctx.scatterU);
104: ISDestroy(&is);
105: i = 2*rank + 1;
106: ISCreateGeneral(PETSC_COMM_WORLD,1,&i,PETSC_COPY_VALUES,&is);
107: VecScatterCreate(ctx.V,NULL,ctx.UV,is,&ctx.scatterV);
108: ISDestroy(&is);
110: VecSet(U,1.0);
111: TSSolve(ts,U);
113: VecDestroy(&ctx.V);
114: VecDestroy(&ctx.UV);
115: VecScatterDestroy(&ctx.scatterU);
116: VecScatterDestroy(&ctx.scatterV);
117: VecDestroy(&tsrhs);
118: VecDestroy(&U);
119: SNESDestroy(&ctx.snes);
120: TSDestroy(&ts);
121: PetscFinalize();
122: return ierr;
123: }
125: /*
126: Defines the RHS function that is passed to the time-integrator.
128: Solves F(U,V) for V and then computes f(U,V)
129: */
130: PetscErrorCode TSFunction(TS ts,PetscReal t,Vec U,Vec F,void *actx)
131: {
132: AppCtx *ctx = (AppCtx*)actx;
136: ctx->t = t;
137: VecScatterBegin(ctx->scatterU,U,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);
138: VecScatterEnd(ctx->scatterU,U,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);
139: SNESSolve(ctx->snes,NULL,ctx->V);
140: VecScatterBegin(ctx->scatterV,ctx->V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);
141: VecScatterEnd(ctx->scatterV,ctx->V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);
142: (*ctx->f)(t,ctx->UV,F);
143: return(0);
144: }
146: /*
147: Defines the nonlinear function that is passed to the nonlinear solver
149: */
150: PetscErrorCode SNESFunction(SNES snes,Vec V,Vec F,void *actx)
151: {
152: AppCtx *ctx = (AppCtx*)actx;
156: VecScatterBegin(ctx->scatterV,V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);
157: VecScatterEnd(ctx->scatterV,V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);
158: (*ctx->F)(ctx->t,ctx->UV,F);
159: return(0);
160: }
162: /*TEST
164: test:
165: args: -ts_monitor -ts_view
167: TEST*/