xref: /petsc/src/ts/tutorials/phasefield/biharmonic2.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Solves biharmonic equation in 1d.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown   Solves the equation biharmonic equation in split form
6c4762a1bSJed Brown 
7c4762a1bSJed Brown     w = -kappa \Delta u
8c4762a1bSJed Brown     u_t =  \Delta w
9c4762a1bSJed Brown     -1  <= u <= 1
10c4762a1bSJed Brown     Periodic boundary conditions
11c4762a1bSJed Brown 
12c4762a1bSJed Brown Evolve the biharmonic heat equation with bounds:  (same as biharmonic)
13c4762a1bSJed Brown ---------------
14c4762a1bSJed Brown ./biharmonic2 -ts_monitor -snes_monitor -ts_monitor_draw_solution  -pc_type lu  -draw_pause .1 -snes_converged_reason  -ts_type beuler  -da_refine 5 -draw_fields 1 -ts_dt 9.53674e-9
15c4762a1bSJed Brown 
16c4762a1bSJed Brown     w = -kappa \Delta u  + u^3  - u
17c4762a1bSJed Brown     u_t =  \Delta w
18c4762a1bSJed Brown     -1  <= u <= 1
19c4762a1bSJed Brown     Periodic boundary conditions
20c4762a1bSJed Brown 
21c4762a1bSJed Brown Evolve the Cahn-Hillard equations: (this fails after a few timesteps 12/17/2017)
22c4762a1bSJed Brown ---------------
23c4762a1bSJed Brown ./biharmonic2 -ts_monitor -snes_monitor -ts_monitor_draw_solution  -pc_type lu  -draw_pause .1 -snes_converged_reason   -ts_type beuler    -da_refine 6  -draw_fields 1  -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard
24c4762a1bSJed Brown 
25c4762a1bSJed Brown */
26c4762a1bSJed Brown #include <petscdm.h>
27c4762a1bSJed Brown #include <petscdmda.h>
28c4762a1bSJed Brown #include <petscts.h>
29c4762a1bSJed Brown #include <petscdraw.h>
30c4762a1bSJed Brown 
31c4762a1bSJed Brown /*
32c4762a1bSJed Brown    User-defined routines
33c4762a1bSJed Brown */
34c4762a1bSJed Brown extern PetscErrorCode FormFunction(TS,PetscReal,Vec,Vec,Vec,void*),FormInitialSolution(DM,Vec,PetscReal);
35c4762a1bSJed Brown typedef struct {PetscBool cahnhillard;PetscReal kappa;PetscInt energy;PetscReal tol;PetscReal theta;PetscReal theta_c;} UserCtx;
36c4762a1bSJed Brown 
37c4762a1bSJed Brown int main(int argc,char **argv)
38c4762a1bSJed Brown {
39c4762a1bSJed Brown   TS             ts;                           /* nonlinear solver */
40c4762a1bSJed Brown   Vec            x,r;                          /* solution, residual vectors */
41c4762a1bSJed Brown   Mat            J;                            /* Jacobian matrix */
42c4762a1bSJed Brown   PetscInt       steps,Mx;
43c4762a1bSJed Brown   DM             da;
44c4762a1bSJed Brown   MatFDColoring  matfdcoloring;
45c4762a1bSJed Brown   ISColoring     iscoloring;
46c4762a1bSJed Brown   PetscReal      dt;
47c4762a1bSJed Brown   PetscReal      vbounds[] = {-100000,100000,-1.1,1.1};
48c4762a1bSJed Brown   SNES           snes;
49c4762a1bSJed Brown   UserCtx        ctx;
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52c4762a1bSJed Brown      Initialize program
53c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54*327415f7SBarry Smith   PetscFunctionBeginUser;
559566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
56c4762a1bSJed Brown   ctx.kappa = 1.0;
579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL));
58c4762a1bSJed Brown   ctx.cahnhillard = PETSC_FALSE;
59c4762a1bSJed Brown 
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-cahn-hillard",&ctx.cahnhillard,NULL));
619566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),2,vbounds));
629566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),600,600));
63c4762a1bSJed Brown   ctx.energy = 1;
649566063dSJacob Faibussowitsch   /*PetscCall(PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL));*/
659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL));
66c4762a1bSJed Brown   ctx.tol     = 1.0e-8;
679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL));
68c4762a1bSJed Brown   ctx.theta   = .001;
69c4762a1bSJed Brown   ctx.theta_c = 1.0;
709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL));
719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL));
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
75c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
769566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10,2,2,NULL,&da));
779566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
789566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
799566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da,0,"Biharmonic heat equation: w = -kappa*u_xx"));
809566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da,1,"Biharmonic heat equation: u"));
819566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0));
82c4762a1bSJed Brown   dt   = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx);
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
86c4762a1bSJed Brown      vectors that are the same types
87c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
889566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&x));
899566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&r));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92c4762a1bSJed Brown      Create timestepping solver context
93c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
949566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
959566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts,da));
969566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
979566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts,NULL,FormFunction,&ctx));
989566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,.02));
999566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_INTERPOLATE));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
103c4762a1bSJed Brown 
104c4762a1bSJed Brown <     Set Jacobian matrix data structure and default Jacobian evaluation
105c4762a1bSJed Brown      routine. User can override with:
106c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
107c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
108c4762a1bSJed Brown      -snes_mf_operator : form preconditioning matrix as set by the user,
109c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
110c4762a1bSJed Brown                          products within Newton-Krylov method
111c4762a1bSJed Brown 
112c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1139566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts,&snes));
1149566063dSJacob Faibussowitsch   PetscCall(DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring));
1159566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(da,MATAIJ));
1169566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(da,&J));
1179566063dSJacob Faibussowitsch   PetscCall(MatFDColoringCreate(J,iscoloring,&matfdcoloring));
1189566063dSJacob Faibussowitsch   PetscCall(MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts));
1199566063dSJacob Faibussowitsch   PetscCall(MatFDColoringSetFromOptions(matfdcoloring));
1209566063dSJacob Faibussowitsch   PetscCall(MatFDColoringSetUp(J,iscoloring,matfdcoloring));
1219566063dSJacob Faibussowitsch   PetscCall(ISColoringDestroy(&iscoloring));
1229566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring));
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125c4762a1bSJed Brown      Customize nonlinear solver
126c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1279566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSBEULER));
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130c4762a1bSJed Brown      Set initial conditions
131c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1329566063dSJacob Faibussowitsch   PetscCall(FormInitialSolution(da,x,ctx.kappa));
1339566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,dt));
1349566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,x));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137c4762a1bSJed Brown      Set runtime options
138c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1399566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142c4762a1bSJed Brown      Solve nonlinear system
143c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1449566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,x));
1459566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&steps));
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
149c4762a1bSJed Brown      are no longer needed.
150c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1519566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1529566063dSJacob Faibussowitsch   PetscCall(MatFDColoringDestroy(&matfdcoloring));
1539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
1559566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1569566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
157c4762a1bSJed Brown 
1589566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
159b122ec5aSJacob Faibussowitsch   return 0;
160c4762a1bSJed Brown }
161c4762a1bSJed Brown 
162c4762a1bSJed Brown typedef struct {PetscScalar w,u;} Field;
163c4762a1bSJed Brown /* ------------------------------------------------------------------- */
164c4762a1bSJed Brown /*
165c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
166c4762a1bSJed Brown 
167c4762a1bSJed Brown    Input Parameters:
168c4762a1bSJed Brown .  ts - the TS context
169c4762a1bSJed Brown .  X - input vector
170c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
171c4762a1bSJed Brown 
172c4762a1bSJed Brown    Output Parameter:
173c4762a1bSJed Brown .  F - function vector
174c4762a1bSJed Brown  */
175c4762a1bSJed Brown PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec Xdot,Vec F,void *ptr)
176c4762a1bSJed Brown {
177c4762a1bSJed Brown   DM             da;
178c4762a1bSJed Brown   PetscInt       i,Mx,xs,xm;
179c4762a1bSJed Brown   PetscReal      hx,sx;
180c4762a1bSJed Brown   Field          *x,*xdot,*f;
181c4762a1bSJed Brown   Vec            localX,localXdot;
182c4762a1bSJed Brown   UserCtx        *ctx = (UserCtx*)ptr;
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   PetscFunctionBegin;
1859566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
1869566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localX));
1879566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localXdot));
1889566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   /*
193c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
194c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
195c4762a1bSJed Brown      By placing code between these two statements, computations can be
196c4762a1bSJed Brown      done while messages are in transition.
197c4762a1bSJed Brown   */
1989566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
1999566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
2009566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,Xdot,INSERT_VALUES,localXdot));
2019566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,Xdot,INSERT_VALUES,localXdot));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /*
204c4762a1bSJed Brown      Get pointers to vector data
205c4762a1bSJed Brown   */
2069566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,localX,&x));
2079566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,localXdot,&xdot));
2089566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,F,&f));
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   /*
211c4762a1bSJed Brown      Get local grid boundaries
212c4762a1bSJed Brown   */
2139566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   /*
216c4762a1bSJed Brown      Compute function over the locally owned part of the grid
217c4762a1bSJed Brown   */
218c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
219c4762a1bSJed Brown     f[i].w =  x[i].w + ctx->kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx;
220c4762a1bSJed Brown     if (ctx->cahnhillard) {
221c4762a1bSJed Brown       switch (ctx->energy) {
222c4762a1bSJed Brown       case 1: /* double well */
223c4762a1bSJed Brown         f[i].w += -x[i].u*x[i].u*x[i].u + x[i].u;
224c4762a1bSJed Brown         break;
225c4762a1bSJed Brown       case 2: /* double obstacle */
226c4762a1bSJed Brown         f[i].w += x[i].u;
227c4762a1bSJed Brown         break;
228c4762a1bSJed Brown       case 3: /* logarithmic */
229c4762a1bSJed Brown         if (PetscRealPart(x[i].u) < -1.0 + 2.0*ctx->tol)     f[i].w += .5*ctx->theta*(-PetscLogReal(ctx->tol) + PetscLogScalar((1.0-x[i].u)/2.0)) + ctx->theta_c*x[i].u;
230c4762a1bSJed Brown         else if (PetscRealPart(x[i].u) > 1.0 - 2.0*ctx->tol) f[i].w += .5*ctx->theta*(-PetscLogScalar((1.0+x[i].u)/2.0) + PetscLogReal(ctx->tol)) + ctx->theta_c*x[i].u;
231c4762a1bSJed Brown         else                                                 f[i].w += .5*ctx->theta*(-PetscLogScalar((1.0+x[i].u)/2.0) + PetscLogScalar((1.0-x[i].u)/2.0)) + ctx->theta_c*x[i].u;
232c4762a1bSJed Brown         break;
233c4762a1bSJed Brown       }
234c4762a1bSJed Brown     }
235c4762a1bSJed Brown     f[i].u = xdot[i].u - (x[i-1].w + x[i+1].w - 2.0*x[i].w)*sx;
236c4762a1bSJed Brown   }
237c4762a1bSJed Brown 
238c4762a1bSJed Brown   /*
239c4762a1bSJed Brown      Restore vectors
240c4762a1bSJed Brown   */
2419566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,localXdot,&xdot));
2429566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,localX,&x));
2439566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,F,&f));
2449566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localX));
2459566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localXdot));
246c4762a1bSJed Brown   PetscFunctionReturn(0);
247c4762a1bSJed Brown }
248c4762a1bSJed Brown 
249c4762a1bSJed Brown /* ------------------------------------------------------------------- */
250c4762a1bSJed Brown PetscErrorCode FormInitialSolution(DM da,Vec X,PetscReal kappa)
251c4762a1bSJed Brown {
252c4762a1bSJed Brown   PetscInt       i,xs,xm,Mx,xgs,xgm;
253c4762a1bSJed Brown   Field          *x;
254c4762a1bSJed Brown   PetscReal      hx,xx,r,sx;
255c4762a1bSJed Brown   Vec            Xg;
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   PetscFunctionBegin;
2589566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   hx = 1.0/(PetscReal)Mx;
261c4762a1bSJed Brown   sx = 1.0/(hx*hx);
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   /*
264c4762a1bSJed Brown      Get pointers to vector data
265c4762a1bSJed Brown   */
2669566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(da,&Xg));
2679566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,Xg,&x));
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   /*
270c4762a1bSJed Brown      Get local grid boundaries
271c4762a1bSJed Brown   */
2729566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
2739566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da,&xgs,NULL,NULL,&xgm,NULL,NULL));
274c4762a1bSJed Brown 
275c4762a1bSJed Brown   /*
276c4762a1bSJed Brown      Compute u function over the locally owned part of the grid including ghost points
277c4762a1bSJed Brown   */
278c4762a1bSJed Brown   for (i=xgs; i<xgs+xgm; i++) {
279c4762a1bSJed Brown     xx = i*hx;
280c4762a1bSJed Brown     r = PetscSqrtReal((xx-.5)*(xx-.5));
281c4762a1bSJed Brown     if (r < .125) x[i].u = 1.0;
282c4762a1bSJed Brown     else          x[i].u = -.50;
283c4762a1bSJed Brown     /* fill in x[i].w so that valgrind doesn't detect use of uninitialized memory */
284c4762a1bSJed Brown     x[i].w = 0;
285c4762a1bSJed Brown   }
286c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) x[i].w = -kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx;
287c4762a1bSJed Brown 
288c4762a1bSJed Brown   /*
289c4762a1bSJed Brown      Restore vectors
290c4762a1bSJed Brown   */
2919566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,Xg,&x));
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   /* Grab only the global part of the vector */
2949566063dSJacob Faibussowitsch   PetscCall(VecSet(X,0));
2959566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da,Xg,ADD_VALUES,X));
2969566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da,Xg,ADD_VALUES,X));
2979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Xg));
298c4762a1bSJed Brown   PetscFunctionReturn(0);
299c4762a1bSJed Brown }
300c4762a1bSJed Brown 
301c4762a1bSJed Brown /*TEST
302c4762a1bSJed Brown 
303c4762a1bSJed Brown    build:
304c4762a1bSJed Brown      requires: !complex !single
305c4762a1bSJed Brown 
306c4762a1bSJed Brown    test:
307c4762a1bSJed Brown      args: -ts_monitor -snes_monitor  -pc_type lu   -snes_converged_reason  -ts_type beuler  -da_refine 5 -ts_dt 9.53674e-9 -ts_max_steps 50
308c4762a1bSJed Brown      requires: x
309c4762a1bSJed Brown 
310c4762a1bSJed Brown TEST*/
311