xref: /petsc/src/ts/tutorials/ex17.c (revision d0609ced746bc51b019815ca91d747429db24893)
1c4762a1bSJed Brown static const char help[] = "Time-dependent PDE in 1d. Simplified from ex15.c for illustrating how to solve DAEs. \n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown    u_t = uxx
4c4762a1bSJed Brown    0 < x < 1;
5c4762a1bSJed Brown    At t=0: u(x) = exp(c*r*r*r), if r=PetscSqrtReal((x-.5)*(x-.5)) < .125
6c4762a1bSJed Brown            u(x) = 0.0           if r >= .125
7c4762a1bSJed Brown 
8c4762a1bSJed Brown    Boundary conditions:
9c4762a1bSJed Brown    Dirichlet BC:
10c4762a1bSJed Brown    At x=0, x=1, u = 0.0
11c4762a1bSJed Brown 
12c4762a1bSJed Brown    Neumann BC:
13c4762a1bSJed Brown    At x=0, x=1: du(x,t)/dx = 0
14c4762a1bSJed Brown 
15c4762a1bSJed Brown    mpiexec -n 2 ./ex17 -da_grid_x 40 -ts_max_steps 2 -snes_monitor -ksp_monitor
16c4762a1bSJed Brown          ./ex17 -da_grid_x 40 -monitor_solution
17c4762a1bSJed Brown          ./ex17 -da_grid_x 100  -ts_type theta -ts_theta_theta 0.5 # Midpoint is not L-stable
18c4762a1bSJed Brown          ./ex17 -jac_type fd_coloring  -da_grid_x 500 -boundary 1
19c4762a1bSJed Brown          ./ex17 -da_grid_x 100  -ts_type gl -ts_adapt_type none -ts_max_steps 2
20c4762a1bSJed Brown */
21c4762a1bSJed Brown 
22c4762a1bSJed Brown #include <petscdm.h>
23c4762a1bSJed Brown #include <petscdmda.h>
24c4762a1bSJed Brown #include <petscts.h>
25c4762a1bSJed Brown 
26c4762a1bSJed Brown typedef enum {JACOBIAN_ANALYTIC,JACOBIAN_FD_COLORING,JACOBIAN_FD_FULL} JacobianType;
27c4762a1bSJed Brown static const char *const JacobianTypes[] = {"analytic","fd_coloring","fd_full","JacobianType","fd_",0};
28c4762a1bSJed Brown 
29c4762a1bSJed Brown /*
30c4762a1bSJed Brown    User-defined data structures and routines
31c4762a1bSJed Brown */
32c4762a1bSJed Brown typedef struct {
33c4762a1bSJed Brown   PetscReal c;
34c4762a1bSJed Brown   PetscInt  boundary;            /* Type of boundary condition */
35c4762a1bSJed Brown   PetscBool viewJacobian;
36c4762a1bSJed Brown } AppCtx;
37c4762a1bSJed Brown 
38c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
39c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
40c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*);
41c4762a1bSJed Brown 
42c4762a1bSJed Brown int main(int argc,char **argv)
43c4762a1bSJed Brown {
44c4762a1bSJed Brown   TS             ts;                   /* nonlinear solver */
45c4762a1bSJed Brown   Vec            u;                    /* solution, residual vectors */
46c4762a1bSJed Brown   Mat            J;                    /* Jacobian matrix */
47c4762a1bSJed Brown   PetscInt       nsteps;
48c4762a1bSJed Brown   PetscReal      vmin,vmax,norm;
49c4762a1bSJed Brown   DM             da;
50c4762a1bSJed Brown   PetscReal      ftime,dt;
51c4762a1bSJed Brown   AppCtx         user;              /* user-defined work context */
52c4762a1bSJed Brown   JacobianType   jacType;
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
58c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
599566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,1,1,NULL,&da));
609566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
619566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64c4762a1bSJed Brown      Extract global vectors from DMDA;
65c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
669566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&u));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* Initialize user application context */
69c4762a1bSJed Brown   user.c            = -30.0;
70c4762a1bSJed Brown   user.boundary     = 0;  /* 0: Dirichlet BC; 1: Neumann BC */
71c4762a1bSJed Brown   user.viewJacobian = PETSC_FALSE;
72c4762a1bSJed Brown 
739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-boundary",&user.boundary,NULL));
749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-viewJacobian",&user.viewJacobian));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77c4762a1bSJed Brown      Create timestepping solver context
78c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
799566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
809566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
819566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSTHETA));
829566063dSJacob Faibussowitsch   PetscCall(TSThetaSetTheta(ts,1.0)); /* Make the Theta method behave like backward Euler */
839566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts,NULL,FormIFunction,&user));
84c4762a1bSJed Brown 
859566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(da,MATAIJ));
869566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(da,&J));
87c4762a1bSJed Brown   jacType = JACOBIAN_ANALYTIC; /* use user-provide Jacobian */
88c4762a1bSJed Brown 
899566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts,da)); /* Use TSGetDM() to access. Setting here allows easy use of geometric multigrid. */
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   ftime = 1.0;
929566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,ftime));
939566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96c4762a1bSJed Brown      Set initial conditions
97c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
989566063dSJacob Faibussowitsch   PetscCall(FormInitialSolution(ts,u,&user));
999566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,u));
100c4762a1bSJed Brown   dt   = .01;
1019566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,dt));
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   /* Use slow fd Jacobian or fast fd Jacobian with colorings.
104c4762a1bSJed Brown      Note: this requirs snes which is not created until TSSetUp()/TSSetFromOptions() is called */
105*d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Options for Jacobian evaluation",NULL);
1069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-jac_type","Type of Jacobian","",JacobianTypes,(PetscEnum)jacType,(PetscEnum*)&jacType,0));
107*d0609cedSBarry Smith   PetscOptionsEnd();
108c4762a1bSJed Brown   if (jacType == JACOBIAN_ANALYTIC) {
1099566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(ts,J,J,FormIJacobian,&user));
110c4762a1bSJed Brown   } else if (jacType == JACOBIAN_FD_COLORING) {
111c4762a1bSJed Brown     SNES snes;
1129566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts,&snes));
1139566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,0));
114c4762a1bSJed Brown   } else if (jacType == JACOBIAN_FD_FULL) {
115c4762a1bSJed Brown     SNES snes;
1169566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts,&snes));
1179566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefault,&user));
118c4762a1bSJed Brown   }
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121c4762a1bSJed Brown      Set runtime options
122c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1239566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126c4762a1bSJed Brown      Integrate ODE system
127c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1289566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,u));
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131c4762a1bSJed Brown    Compute diagnostics of the solution
132c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1339566063dSJacob Faibussowitsch   PetscCall(VecNorm(u,NORM_1,&norm));
1349566063dSJacob Faibussowitsch   PetscCall(VecMax(u,NULL,&vmax));
1359566063dSJacob Faibussowitsch   PetscCall(VecMin(u,NULL,&vmin));
1369566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&nsteps));
1379566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts,&ftime));
1389566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"timestep %D: time %g, solution norm %g, max %g, min %g\n",nsteps,(double)ftime,(double)norm,(double)vmax,(double)vmin));
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141c4762a1bSJed Brown      Free work space.
142c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1439566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1449566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
1459566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1469566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1479566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
148b122ec5aSJacob Faibussowitsch   return 0;
149c4762a1bSJed Brown }
150c4762a1bSJed Brown /* ------------------------------------------------------------------- */
151c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr)
152c4762a1bSJed Brown {
153c4762a1bSJed Brown   AppCtx         *user=(AppCtx*)ptr;
154c4762a1bSJed Brown   DM             da;
155c4762a1bSJed Brown   PetscInt       i,Mx,xs,xm;
156c4762a1bSJed Brown   PetscReal      hx,sx;
157c4762a1bSJed Brown   PetscScalar    *u,*udot,*f;
158c4762a1bSJed Brown   Vec            localU;
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   PetscFunctionBeginUser;
1619566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
1629566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localU));
1639566063dSJacob 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));
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /*
168c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
169c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
170c4762a1bSJed Brown      By placing code between these two statements, computations can be
171c4762a1bSJed Brown      done while messages are in transition.
172c4762a1bSJed Brown   */
1739566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
1749566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   /* Get pointers to vector data */
1779566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,localU,&u));
1789566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,Udot,&udot));
1799566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,F,&f));
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   /* Get local grid boundaries */
1829566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
185c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
186c4762a1bSJed Brown     if (user->boundary == 0) { /* Dirichlet BC */
187c4762a1bSJed Brown       if (i == 0 || i == Mx-1) f[i] = u[i]; /* F = U */
188c4762a1bSJed Brown       else                     f[i] = udot[i] + (2.*u[i] - u[i-1] - u[i+1])*sx;
189c4762a1bSJed Brown     } else { /* Neumann BC */
190c4762a1bSJed Brown       if (i == 0)         f[i] = u[0] - u[1];
191c4762a1bSJed Brown       else if (i == Mx-1) f[i] = u[i] - u[i-1];
192c4762a1bSJed Brown       else                f[i] = udot[i] + (2.*u[i] - u[i-1] - u[i+1])*sx;
193c4762a1bSJed Brown     }
194c4762a1bSJed Brown   }
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   /* Restore vectors */
1979566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,localU,&u));
1989566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,Udot,&udot));
1999566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,F,&f));
2009566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localU));
201c4762a1bSJed Brown   PetscFunctionReturn(0);
202c4762a1bSJed Brown }
203c4762a1bSJed Brown 
204c4762a1bSJed Brown /* --------------------------------------------------------------------- */
205c4762a1bSJed Brown /*
206c4762a1bSJed Brown   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
207c4762a1bSJed Brown */
208c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat J,Mat Jpre,void *ctx)
209c4762a1bSJed Brown {
210c4762a1bSJed Brown   PetscInt       i,rstart,rend,Mx;
211c4762a1bSJed Brown   PetscReal      hx,sx;
212c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;
213c4762a1bSJed Brown   DM             da;
214c4762a1bSJed Brown   MatStencil     col[3],row;
215c4762a1bSJed Brown   PetscInt       nc;
216c4762a1bSJed Brown   PetscScalar    vals[3];
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   PetscFunctionBeginUser;
2199566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
2209566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Jpre,&rstart,&rend));
2219566063dSJacob 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));
222c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
223c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
224c4762a1bSJed Brown     nc    = 0;
225c4762a1bSJed Brown     row.i = i;
226c4762a1bSJed Brown     if (user->boundary == 0 && (i == 0 || i == Mx-1)) {
227c4762a1bSJed Brown       col[nc].i = i; vals[nc++] = 1.0;
228c4762a1bSJed Brown     } else if (user->boundary > 0 && i == 0) { /* Left Neumann */
229c4762a1bSJed Brown       col[nc].i = i;   vals[nc++] = 1.0;
230c4762a1bSJed Brown       col[nc].i = i+1; vals[nc++] = -1.0;
231c4762a1bSJed Brown     } else if (user->boundary > 0 && i == Mx-1) { /* Right Neumann */
232c4762a1bSJed Brown       col[nc].i = i-1; vals[nc++] = -1.0;
233c4762a1bSJed Brown       col[nc].i = i;   vals[nc++] = 1.0;
234c4762a1bSJed Brown     } else {                    /* Interior */
235c4762a1bSJed Brown       col[nc].i = i-1; vals[nc++] = -1.0*sx;
236c4762a1bSJed Brown       col[nc].i = i;   vals[nc++] = 2.0*sx + a;
237c4762a1bSJed Brown       col[nc].i = i+1; vals[nc++] = -1.0*sx;
238c4762a1bSJed Brown     }
2399566063dSJacob Faibussowitsch     PetscCall(MatSetValuesStencil(Jpre,1,&row,nc,col,vals,INSERT_VALUES));
240c4762a1bSJed Brown   }
241c4762a1bSJed Brown 
2429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY));
2439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY));
244c4762a1bSJed Brown   if (J != Jpre) {
2459566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
2469566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
247c4762a1bSJed Brown   }
248c4762a1bSJed Brown   if (user->viewJacobian) {
2499566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Jpre:\n"));
2509566063dSJacob Faibussowitsch     PetscCall(MatView(Jpre,PETSC_VIEWER_STDOUT_WORLD));
251c4762a1bSJed Brown   }
252c4762a1bSJed Brown   PetscFunctionReturn(0);
253c4762a1bSJed Brown }
254c4762a1bSJed Brown 
255c4762a1bSJed Brown /* ------------------------------------------------------------------- */
256c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec U,void *ptr)
257c4762a1bSJed Brown {
258c4762a1bSJed Brown   AppCtx         *user=(AppCtx*)ptr;
259c4762a1bSJed Brown   PetscReal      c    =user->c;
260c4762a1bSJed Brown   DM             da;
261c4762a1bSJed Brown   PetscInt       i,xs,xm,Mx;
262c4762a1bSJed Brown   PetscScalar    *u;
263c4762a1bSJed Brown   PetscReal      hx,x,r;
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   PetscFunctionBeginUser;
2669566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
2679566063dSJacob 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));
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1);
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   /* Get pointers to vector data */
2729566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,U,&u));
273c4762a1bSJed Brown 
274c4762a1bSJed Brown   /* Get local grid boundaries */
2759566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
276c4762a1bSJed Brown 
277c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
278c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
279c4762a1bSJed Brown     x = i*hx;
280c4762a1bSJed Brown     r = PetscSqrtReal((x-.5)*(x-.5));
281c4762a1bSJed Brown     if (r < .125) u[i] = PetscExpReal(c*r*r*r);
282c4762a1bSJed Brown     else          u[i] = 0.0;
283c4762a1bSJed Brown   }
284c4762a1bSJed Brown 
285c4762a1bSJed Brown   /* Restore vectors */
2869566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,U,&u));
287c4762a1bSJed Brown   PetscFunctionReturn(0);
288c4762a1bSJed Brown }
289c4762a1bSJed Brown 
290c4762a1bSJed Brown /*TEST
291c4762a1bSJed Brown 
292c4762a1bSJed Brown     test:
293c4762a1bSJed Brown       requires: !single
294c4762a1bSJed Brown       args: -da_grid_x 40 -ts_max_steps 2 -snes_monitor_short -ksp_monitor_short -ts_monitor
295c4762a1bSJed Brown 
296c4762a1bSJed Brown     test:
297c4762a1bSJed Brown       suffix: 2
298c4762a1bSJed Brown       requires: !single
299c4762a1bSJed Brown       args: -da_grid_x 100 -ts_type theta -ts_theta_theta 0.5
300c4762a1bSJed Brown 
301c4762a1bSJed Brown TEST*/
302