xref: /petsc/src/ts/tutorials/ex17.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static const char help[] = "Time-dependent PDE in 1d. Simplified from ex15.c for illustrating how to solve DAEs. \n";
2*c4762a1bSJed Brown /*
3*c4762a1bSJed Brown    u_t = uxx
4*c4762a1bSJed Brown    0 < x < 1;
5*c4762a1bSJed Brown    At t=0: u(x) = exp(c*r*r*r), if r=PetscSqrtReal((x-.5)*(x-.5)) < .125
6*c4762a1bSJed Brown            u(x) = 0.0           if r >= .125
7*c4762a1bSJed Brown 
8*c4762a1bSJed Brown 
9*c4762a1bSJed Brown    Boundary conditions:
10*c4762a1bSJed Brown    Dirichlet BC:
11*c4762a1bSJed Brown    At x=0, x=1, u = 0.0
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown    Neumann BC:
14*c4762a1bSJed Brown    At x=0, x=1: du(x,t)/dx = 0
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown    mpiexec -n 2 ./ex17 -da_grid_x 40 -ts_max_steps 2 -snes_monitor -ksp_monitor
17*c4762a1bSJed Brown          ./ex17 -da_grid_x 40 -monitor_solution
18*c4762a1bSJed Brown          ./ex17 -da_grid_x 100  -ts_type theta -ts_theta_theta 0.5 # Midpoint is not L-stable
19*c4762a1bSJed Brown          ./ex17 -jac_type fd_coloring  -da_grid_x 500 -boundary 1
20*c4762a1bSJed Brown          ./ex17 -da_grid_x 100  -ts_type gl -ts_adapt_type none -ts_max_steps 2
21*c4762a1bSJed Brown */
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown #include <petscdm.h>
24*c4762a1bSJed Brown #include <petscdmda.h>
25*c4762a1bSJed Brown #include <petscts.h>
26*c4762a1bSJed Brown 
27*c4762a1bSJed Brown typedef enum {JACOBIAN_ANALYTIC,JACOBIAN_FD_COLORING,JACOBIAN_FD_FULL} JacobianType;
28*c4762a1bSJed Brown static const char *const JacobianTypes[] = {"analytic","fd_coloring","fd_full","JacobianType","fd_",0};
29*c4762a1bSJed Brown 
30*c4762a1bSJed Brown /*
31*c4762a1bSJed Brown    User-defined data structures and routines
32*c4762a1bSJed Brown */
33*c4762a1bSJed Brown typedef struct {
34*c4762a1bSJed Brown   PetscReal c;
35*c4762a1bSJed Brown   PetscInt  boundary;            /* Type of boundary condition */
36*c4762a1bSJed Brown   PetscBool viewJacobian;
37*c4762a1bSJed Brown } AppCtx;
38*c4762a1bSJed Brown 
39*c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
40*c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
41*c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*);
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown int main(int argc,char **argv)
44*c4762a1bSJed Brown {
45*c4762a1bSJed Brown   TS             ts;                   /* nonlinear solver */
46*c4762a1bSJed Brown   Vec            u;                    /* solution, residual vectors */
47*c4762a1bSJed Brown   Mat            J;                    /* Jacobian matrix */
48*c4762a1bSJed Brown   PetscInt       nsteps;
49*c4762a1bSJed Brown   PetscReal      vmin,vmax,norm;
50*c4762a1bSJed Brown   PetscErrorCode ierr;
51*c4762a1bSJed Brown   DM             da;
52*c4762a1bSJed Brown   PetscReal      ftime,dt;
53*c4762a1bSJed Brown   AppCtx         user;              /* user-defined work context */
54*c4762a1bSJed Brown   JacobianType   jacType;
55*c4762a1bSJed Brown 
56*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
57*c4762a1bSJed Brown 
58*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59*c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
60*c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61*c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,1,1,NULL,&da);CHKERRQ(ierr);
62*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
64*c4762a1bSJed Brown 
65*c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66*c4762a1bSJed Brown      Extract global vectors from DMDA;
67*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&u);CHKERRQ(ierr);
69*c4762a1bSJed Brown 
70*c4762a1bSJed Brown   /* Initialize user application context */
71*c4762a1bSJed Brown   user.c            = -30.0;
72*c4762a1bSJed Brown   user.boundary     = 0;  /* 0: Dirichlet BC; 1: Neumann BC */
73*c4762a1bSJed Brown   user.viewJacobian = PETSC_FALSE;
74*c4762a1bSJed Brown 
75*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-boundary",&user.boundary,NULL);CHKERRQ(ierr);
76*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-viewJacobian",&user.viewJacobian);CHKERRQ(ierr);
77*c4762a1bSJed Brown 
78*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79*c4762a1bSJed Brown      Create timestepping solver context
80*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
81*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
82*c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr = TSSetType(ts,TSTHETA);CHKERRQ(ierr);
84*c4762a1bSJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); /* Make the Theta method behave like backward Euler */
85*c4762a1bSJed Brown   ierr = TSSetIFunction(ts,NULL,FormIFunction,&user);CHKERRQ(ierr);
86*c4762a1bSJed Brown 
87*c4762a1bSJed Brown   ierr = DMSetMatType(da,MATAIJ);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = DMCreateMatrix(da,&J);CHKERRQ(ierr);
89*c4762a1bSJed Brown   jacType = JACOBIAN_ANALYTIC; /* use user-provide Jacobian */
90*c4762a1bSJed Brown 
91*c4762a1bSJed Brown   ierr = TSSetDM(ts,da);CHKERRQ(ierr); /* Use TSGetDM() to access. Setting here allows easy use of geometric multigrid. */
92*c4762a1bSJed Brown 
93*c4762a1bSJed Brown   ftime = 1.0;
94*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr);
95*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
96*c4762a1bSJed Brown 
97*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98*c4762a1bSJed Brown      Set initial conditions
99*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100*c4762a1bSJed Brown   ierr = FormInitialSolution(ts,u,&user);CHKERRQ(ierr);
101*c4762a1bSJed Brown   ierr = TSSetSolution(ts,u);CHKERRQ(ierr);
102*c4762a1bSJed Brown   dt   = .01;
103*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
104*c4762a1bSJed Brown 
105*c4762a1bSJed Brown 
106*c4762a1bSJed Brown   /* Use slow fd Jacobian or fast fd Jacobian with colorings.
107*c4762a1bSJed Brown      Note: this requirs snes which is not created until TSSetUp()/TSSetFromOptions() is called */
108*c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Options for Jacobian evaluation",NULL);CHKERRQ(ierr);
109*c4762a1bSJed Brown   ierr = PetscOptionsEnum("-jac_type","Type of Jacobian","",JacobianTypes,(PetscEnum)jacType,(PetscEnum*)&jacType,0);CHKERRQ(ierr);
110*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
111*c4762a1bSJed Brown   if (jacType == JACOBIAN_ANALYTIC) {
112*c4762a1bSJed Brown     ierr = TSSetIJacobian(ts,J,J,FormIJacobian,&user);CHKERRQ(ierr);
113*c4762a1bSJed Brown   } else if (jacType == JACOBIAN_FD_COLORING) {
114*c4762a1bSJed Brown     SNES snes;
115*c4762a1bSJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
116*c4762a1bSJed Brown     ierr = SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,0);CHKERRQ(ierr);
117*c4762a1bSJed Brown   } else if (jacType == JACOBIAN_FD_FULL) {
118*c4762a1bSJed Brown     SNES snes;
119*c4762a1bSJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
120*c4762a1bSJed Brown     ierr = SNESSetJacobian(snes,J,J,SNESComputeJacobianDefault,&user);CHKERRQ(ierr);
121*c4762a1bSJed Brown   }
122*c4762a1bSJed Brown 
123*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124*c4762a1bSJed Brown      Set runtime options
125*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
126*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
127*c4762a1bSJed Brown 
128*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129*c4762a1bSJed Brown      Integrate ODE system
130*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131*c4762a1bSJed Brown   ierr = TSSolve(ts,u);CHKERRQ(ierr);
132*c4762a1bSJed Brown 
133*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134*c4762a1bSJed Brown    Compute diagnostics of the solution
135*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136*c4762a1bSJed Brown   ierr = VecNorm(u,NORM_1,&norm);CHKERRQ(ierr);
137*c4762a1bSJed Brown   ierr = VecMax(u,NULL,&vmax);CHKERRQ(ierr);
138*c4762a1bSJed Brown   ierr = VecMin(u,NULL,&vmin);CHKERRQ(ierr);
139*c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&nsteps);CHKERRQ(ierr);
140*c4762a1bSJed Brown   ierr = TSGetTime(ts,&ftime);CHKERRQ(ierr);
141*c4762a1bSJed Brown   ierr = 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);CHKERRQ(ierr);
142*c4762a1bSJed Brown 
143*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144*c4762a1bSJed Brown      Free work space.
145*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146*c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
148*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
149*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
150*c4762a1bSJed Brown   ierr = PetscFinalize();
151*c4762a1bSJed Brown   return ierr;
152*c4762a1bSJed Brown }
153*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
154*c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr)
155*c4762a1bSJed Brown {
156*c4762a1bSJed Brown   AppCtx         *user=(AppCtx*)ptr;
157*c4762a1bSJed Brown   DM             da;
158*c4762a1bSJed Brown   PetscErrorCode ierr;
159*c4762a1bSJed Brown   PetscInt       i,Mx,xs,xm;
160*c4762a1bSJed Brown   PetscReal      hx,sx;
161*c4762a1bSJed Brown   PetscScalar    *u,*udot,*f;
162*c4762a1bSJed Brown   Vec            localU;
163*c4762a1bSJed Brown 
164*c4762a1bSJed Brown   PetscFunctionBeginUser;
165*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
166*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr);
167*c4762a1bSJed Brown   ierr = 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);CHKERRQ(ierr);
168*c4762a1bSJed Brown 
169*c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
170*c4762a1bSJed Brown 
171*c4762a1bSJed Brown   /*
172*c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
173*c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
174*c4762a1bSJed Brown      By placing code between these two statements, computations can be
175*c4762a1bSJed Brown      done while messages are in transition.
176*c4762a1bSJed Brown   */
177*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
178*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
179*c4762a1bSJed Brown 
180*c4762a1bSJed Brown   /* Get pointers to vector data */
181*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr);
182*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,Udot,&udot);CHKERRQ(ierr);
183*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
184*c4762a1bSJed Brown 
185*c4762a1bSJed Brown   /* Get local grid boundaries */
186*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
187*c4762a1bSJed Brown 
188*c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
189*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
190*c4762a1bSJed Brown     if (user->boundary == 0) { /* Dirichlet BC */
191*c4762a1bSJed Brown       if (i == 0 || i == Mx-1) f[i] = u[i]; /* F = U */
192*c4762a1bSJed Brown       else                     f[i] = udot[i] + (2.*u[i] - u[i-1] - u[i+1])*sx;
193*c4762a1bSJed Brown     } else { /* Neumann BC */
194*c4762a1bSJed Brown       if (i == 0)         f[i] = u[0] - u[1];
195*c4762a1bSJed Brown       else if (i == Mx-1) f[i] = u[i] - u[i-1];
196*c4762a1bSJed Brown       else                f[i] = udot[i] + (2.*u[i] - u[i-1] - u[i+1])*sx;
197*c4762a1bSJed Brown     }
198*c4762a1bSJed Brown   }
199*c4762a1bSJed Brown 
200*c4762a1bSJed Brown   /* Restore vectors */
201*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr);
202*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,Udot,&udot);CHKERRQ(ierr);
203*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
204*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr);
205*c4762a1bSJed Brown   PetscFunctionReturn(0);
206*c4762a1bSJed Brown }
207*c4762a1bSJed Brown 
208*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
209*c4762a1bSJed Brown /*
210*c4762a1bSJed Brown   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
211*c4762a1bSJed Brown */
212*c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat J,Mat Jpre,void *ctx)
213*c4762a1bSJed Brown {
214*c4762a1bSJed Brown   PetscErrorCode ierr;
215*c4762a1bSJed Brown   PetscInt       i,rstart,rend,Mx;
216*c4762a1bSJed Brown   PetscReal      hx,sx;
217*c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;
218*c4762a1bSJed Brown   DM             da;
219*c4762a1bSJed Brown   MatStencil     col[3],row;
220*c4762a1bSJed Brown   PetscInt       nc;
221*c4762a1bSJed Brown   PetscScalar    vals[3];
222*c4762a1bSJed Brown 
223*c4762a1bSJed Brown   PetscFunctionBeginUser;
224*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
225*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(Jpre,&rstart,&rend);CHKERRQ(ierr);
226*c4762a1bSJed Brown   ierr = 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);CHKERRQ(ierr);
227*c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
228*c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
229*c4762a1bSJed Brown     nc    = 0;
230*c4762a1bSJed Brown     row.i = i;
231*c4762a1bSJed Brown     if (user->boundary == 0 && (i == 0 || i == Mx-1)) {
232*c4762a1bSJed Brown       col[nc].i = i; vals[nc++] = 1.0;
233*c4762a1bSJed Brown     } else if (user->boundary > 0 && i == 0) { /* Left Neumann */
234*c4762a1bSJed Brown       col[nc].i = i;   vals[nc++] = 1.0;
235*c4762a1bSJed Brown       col[nc].i = i+1; vals[nc++] = -1.0;
236*c4762a1bSJed Brown     } else if (user->boundary > 0 && i == Mx-1) { /* Right Neumann */
237*c4762a1bSJed Brown       col[nc].i = i-1; vals[nc++] = -1.0;
238*c4762a1bSJed Brown       col[nc].i = i;   vals[nc++] = 1.0;
239*c4762a1bSJed Brown     } else {                    /* Interior */
240*c4762a1bSJed Brown       col[nc].i = i-1; vals[nc++] = -1.0*sx;
241*c4762a1bSJed Brown       col[nc].i = i;   vals[nc++] = 2.0*sx + a;
242*c4762a1bSJed Brown       col[nc].i = i+1; vals[nc++] = -1.0*sx;
243*c4762a1bSJed Brown     }
244*c4762a1bSJed Brown     ierr = MatSetValuesStencil(Jpre,1,&row,nc,col,vals,INSERT_VALUES);CHKERRQ(ierr);
245*c4762a1bSJed Brown   }
246*c4762a1bSJed Brown 
247*c4762a1bSJed Brown   ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
248*c4762a1bSJed Brown   ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
249*c4762a1bSJed Brown   if (J != Jpre) {
250*c4762a1bSJed Brown     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
251*c4762a1bSJed Brown     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
252*c4762a1bSJed Brown   }
253*c4762a1bSJed Brown   if (user->viewJacobian) {
254*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Jpre:\n");CHKERRQ(ierr);
255*c4762a1bSJed Brown     ierr = MatView(Jpre,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
256*c4762a1bSJed Brown   }
257*c4762a1bSJed Brown   PetscFunctionReturn(0);
258*c4762a1bSJed Brown }
259*c4762a1bSJed Brown 
260*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
261*c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec U,void *ptr)
262*c4762a1bSJed Brown {
263*c4762a1bSJed Brown   AppCtx         *user=(AppCtx*)ptr;
264*c4762a1bSJed Brown   PetscReal      c    =user->c;
265*c4762a1bSJed Brown   DM             da;
266*c4762a1bSJed Brown   PetscErrorCode ierr;
267*c4762a1bSJed Brown   PetscInt       i,xs,xm,Mx;
268*c4762a1bSJed Brown   PetscScalar    *u;
269*c4762a1bSJed Brown   PetscReal      hx,x,r;
270*c4762a1bSJed Brown 
271*c4762a1bSJed Brown   PetscFunctionBeginUser;
272*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
273*c4762a1bSJed Brown   ierr = 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);CHKERRQ(ierr);
274*c4762a1bSJed Brown 
275*c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1);
276*c4762a1bSJed Brown 
277*c4762a1bSJed Brown   /* Get pointers to vector data */
278*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
279*c4762a1bSJed Brown 
280*c4762a1bSJed Brown   /* Get local grid boundaries */
281*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
282*c4762a1bSJed Brown 
283*c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
284*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
285*c4762a1bSJed Brown     x = i*hx;
286*c4762a1bSJed Brown     r = PetscSqrtReal((x-.5)*(x-.5));
287*c4762a1bSJed Brown     if (r < .125) u[i] = PetscExpReal(c*r*r*r);
288*c4762a1bSJed Brown     else          u[i] = 0.0;
289*c4762a1bSJed Brown   }
290*c4762a1bSJed Brown 
291*c4762a1bSJed Brown   /* Restore vectors */
292*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
293*c4762a1bSJed Brown   PetscFunctionReturn(0);
294*c4762a1bSJed Brown }
295*c4762a1bSJed Brown 
296*c4762a1bSJed Brown /*TEST
297*c4762a1bSJed Brown 
298*c4762a1bSJed Brown     test:
299*c4762a1bSJed Brown       requires: !single
300*c4762a1bSJed Brown       args: -da_grid_x 40 -ts_max_steps 2 -snes_monitor_short -ksp_monitor_short -ts_monitor
301*c4762a1bSJed Brown 
302*c4762a1bSJed Brown     test:
303*c4762a1bSJed Brown       suffix: 2
304*c4762a1bSJed Brown       requires: !single
305*c4762a1bSJed Brown       args: -da_grid_x 100 -ts_type theta -ts_theta_theta 0.5
306*c4762a1bSJed Brown 
307*c4762a1bSJed Brown 
308*c4762a1bSJed Brown TEST*/
309*c4762a1bSJed Brown 
310