xref: /petsc/src/ts/tutorials/ex21.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] ="Solves a time-dependent nonlinear PDE with lower and upper bounds on the interior grid points. Uses implicit\n\
3c4762a1bSJed Brown timestepping.  Runtime options include:\n\
4c4762a1bSJed Brown   -M <xg>, where <xg> = number of grid points\n\
5c4762a1bSJed Brown   -debug : Activate debugging printouts\n\
6c4762a1bSJed Brown   -nox   : Deactivate x-window graphics\n\
7c4762a1bSJed Brown   -ul   : lower bound\n\
8c4762a1bSJed Brown   -uh  : upper bound\n\n";
9c4762a1bSJed Brown 
10c4762a1bSJed Brown /* ------------------------------------------------------------------------
11c4762a1bSJed Brown 
12c4762a1bSJed Brown    This is a variation of ex2.c to solve the PDE
13c4762a1bSJed Brown 
14c4762a1bSJed Brown                u * u_xx
15c4762a1bSJed Brown          u_t = ---------
16c4762a1bSJed Brown                2*(t+1)^2
17c4762a1bSJed Brown 
18c4762a1bSJed Brown     with box constraints on the interior grid points
19c4762a1bSJed Brown     ul <= u(t,x) <= uh with x != 0,1
20c4762a1bSJed Brown     on the domain 0 <= x <= 1, with boundary conditions
21c4762a1bSJed Brown          u(t,0) = t + 1,  u(t,1) = 2*t + 2,
22c4762a1bSJed Brown     and initial condition
23c4762a1bSJed Brown          u(0,x) = 1 + x*x.
24c4762a1bSJed Brown 
25c4762a1bSJed Brown     The exact solution is:
26c4762a1bSJed Brown          u(t,x) = (1 + x*x) * (1 + t)
27c4762a1bSJed Brown 
28c4762a1bSJed Brown     We use by default the backward Euler method.
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   ------------------------------------------------------------------------- */
31c4762a1bSJed Brown 
32c4762a1bSJed Brown /*
33c4762a1bSJed Brown    Include "petscts.h" to use the PETSc timestepping routines. Note that
34c4762a1bSJed Brown    this file automatically includes "petscsys.h" and other lower-level
35c4762a1bSJed Brown    PETSc include files.
36c4762a1bSJed Brown 
37c4762a1bSJed Brown    Include the "petscdmda.h" to allow us to use the distributed array data
38c4762a1bSJed Brown    structures to manage the parallel grid.
39c4762a1bSJed Brown */
40c4762a1bSJed Brown #include <petscts.h>
41c4762a1bSJed Brown #include <petscdm.h>
42c4762a1bSJed Brown #include <petscdmda.h>
43c4762a1bSJed Brown #include <petscdraw.h>
44c4762a1bSJed Brown 
45c4762a1bSJed Brown /*
46c4762a1bSJed Brown    User-defined application context - contains data needed by the
47c4762a1bSJed Brown    application-provided callback routines.
48c4762a1bSJed Brown */
49c4762a1bSJed Brown typedef struct {
50c4762a1bSJed Brown   MPI_Comm  comm;           /* communicator */
51c4762a1bSJed Brown   DM        da;             /* distributed array data structure */
52c4762a1bSJed Brown   Vec       localwork;      /* local ghosted work vector */
53c4762a1bSJed Brown   Vec       u_local;        /* local ghosted approximate solution vector */
54c4762a1bSJed Brown   Vec       solution;       /* global exact solution vector */
55c4762a1bSJed Brown   PetscInt  m;              /* total number of grid points */
56c4762a1bSJed Brown   PetscReal h;              /* mesh width: h = 1/(m-1) */
57c4762a1bSJed Brown   PetscBool debug;          /* flag (1 indicates activation of debugging printouts) */
58c4762a1bSJed Brown } AppCtx;
59c4762a1bSJed Brown 
60c4762a1bSJed Brown /*
61c4762a1bSJed Brown    User-defined routines, provided below.
62c4762a1bSJed Brown */
63c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec,AppCtx*);
64c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
65c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
66c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
67c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);
68c4762a1bSJed Brown extern PetscErrorCode SetBounds(Vec,Vec,PetscScalar,PetscScalar,AppCtx*);
69c4762a1bSJed Brown 
70c4762a1bSJed Brown int main(int argc,char **argv)
71c4762a1bSJed Brown {
72c4762a1bSJed Brown   AppCtx         appctx;                 /* user-defined application context */
73c4762a1bSJed Brown   TS             ts;                     /* timestepping context */
74c4762a1bSJed Brown   Mat            A;                      /* Jacobian matrix data structure */
75c4762a1bSJed Brown   Vec            u;                      /* approximate solution vector */
76c4762a1bSJed Brown   Vec            r;                      /* residual vector */
77c4762a1bSJed Brown   PetscInt       time_steps_max = 1000;  /* default max timesteps */
78c4762a1bSJed Brown   PetscReal      dt;
79c4762a1bSJed Brown   PetscReal      time_total_max = 100.0; /* default max total time */
80c4762a1bSJed Brown   Vec            xl,xu; /* Lower and upper bounds on variables */
81c4762a1bSJed Brown   PetscScalar    ul=0.0,uh = 3.0;
82c4762a1bSJed Brown   PetscBool      mymonitor;
83c4762a1bSJed Brown   PetscReal      bounds[] = {1.0, 3.3};
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86c4762a1bSJed Brown      Initialize program and set problem parameters
87c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88c4762a1bSJed Brown 
89*327415f7SBarry Smith   PetscFunctionBeginUser;
909566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
919566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,bounds));
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   appctx.comm = PETSC_COMM_WORLD;
94c4762a1bSJed Brown   appctx.m    = 60;
959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&appctx.m,NULL));
969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-ul",&ul,NULL));
979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-uh",&uh,NULL));
989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug));
999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-mymonitor",&mymonitor));
100c4762a1bSJed Brown   appctx.h    = 1.0/(appctx.m-1.0);
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103c4762a1bSJed Brown      Create vector data structures
104c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /*
107c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
108c4762a1bSJed Brown      and to set up the ghost point communication pattern.  There are M
109c4762a1bSJed Brown      total grid values spread equally among all the processors.
110c4762a1bSJed Brown   */
1119566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,appctx.m,1,1,NULL,&appctx.da));
1129566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(appctx.da));
1139566063dSJacob Faibussowitsch   PetscCall(DMSetUp(appctx.da));
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /*
116c4762a1bSJed Brown      Extract global and local vectors from DMDA; we use these to store the
117c4762a1bSJed Brown      approximate solution.  Then duplicate these for remaining vectors that
118c4762a1bSJed Brown      have the same types.
119c4762a1bSJed Brown   */
1209566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(appctx.da,&u));
1219566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(appctx.da,&appctx.u_local));
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   /*
124c4762a1bSJed Brown      Create local work vector for use in evaluating right-hand-side function;
125c4762a1bSJed Brown      create global work vector for storing exact solution.
126c4762a1bSJed Brown   */
1279566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(appctx.u_local,&appctx.localwork));
1289566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&appctx.solution));
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* Create residual vector */
1319566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&r));
132c4762a1bSJed Brown   /* Create lower and upper bound vectors */
1339566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&xl));
1349566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&xu));
1359566063dSJacob Faibussowitsch   PetscCall(SetBounds(xl,xu,ul,uh,&appctx));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138c4762a1bSJed Brown      Create timestepping solver context; set callback routine for
139c4762a1bSJed Brown      right-hand-side function evaluation.
140c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141c4762a1bSJed Brown 
1429566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
1439566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
1449566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts,r,RHSFunction,&appctx));
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147c4762a1bSJed Brown      Set optional user-defined monitoring routine
148c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   if (mymonitor) {
1519566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts,Monitor,&appctx,NULL));
152c4762a1bSJed Brown   }
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155c4762a1bSJed Brown      For nonlinear problems, the user can provide a Jacobian evaluation
156c4762a1bSJed Brown      routine (or use a finite differencing approximation).
157c4762a1bSJed Brown 
158c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine.
159c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
160c4762a1bSJed Brown 
1619566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
1629566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,appctx.m,appctx.m));
1639566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1649566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
1659566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts,A,A,RHSJacobian,&appctx));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168c4762a1bSJed Brown      Set solution vector and initial timestep
169c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   dt   = appctx.h/2.0;
1729566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,dt));
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175c4762a1bSJed Brown      Customize timestepping solver:
176c4762a1bSJed Brown        - Set the solution method to be the Backward Euler method.
177c4762a1bSJed Brown        - Set timestepping duration info
178c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
179c4762a1bSJed Brown      For example,
180c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
181c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
182c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
183c4762a1bSJed Brown 
1849566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSBEULER));
1859566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts,time_steps_max));
1869566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,time_total_max));
1879566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
188c4762a1bSJed Brown   /* Set lower and upper bound on the solution vector for each time step */
1899566063dSJacob Faibussowitsch   PetscCall(TSVISetVariableBounds(ts,xl,xu));
1909566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193c4762a1bSJed Brown      Solve the problem
194c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   /*
197c4762a1bSJed Brown      Evaluate initial conditions
198c4762a1bSJed Brown   */
1999566063dSJacob Faibussowitsch   PetscCall(InitialConditions(u,&appctx));
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   /*
202c4762a1bSJed Brown      Run the timestepping solver
203c4762a1bSJed Brown   */
2049566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,u));
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
208c4762a1bSJed Brown      are no longer needed.
209c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
210c4762a1bSJed Brown 
2119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
2129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xl));
2139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xu));
2149566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
2169566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2179566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&appctx.da));
2189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.localwork));
2199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.solution));
2209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.u_local));
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   /*
223c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
224c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
225c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
226c4762a1bSJed Brown          options are chosen (e.g., -log_view).
227c4762a1bSJed Brown   */
2289566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
229b122ec5aSJacob Faibussowitsch   return 0;
230c4762a1bSJed Brown }
231c4762a1bSJed Brown /* --------------------------------------------------------------------- */
232c4762a1bSJed Brown /*
233c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
234c4762a1bSJed Brown 
235c4762a1bSJed Brown    Input Parameters:
236c4762a1bSJed Brown    u - uninitialized solution vector (global)
237c4762a1bSJed Brown    appctx - user-defined application context
238c4762a1bSJed Brown 
239c4762a1bSJed Brown    Output Parameter:
240c4762a1bSJed Brown    u - vector with solution at initial time (global)
241c4762a1bSJed Brown */
242c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
243c4762a1bSJed Brown {
244c4762a1bSJed Brown   PetscScalar    *u_localptr,h = appctx->h,x;
245c4762a1bSJed Brown   PetscInt       i,mybase,myend;
246c4762a1bSJed Brown 
247c4762a1bSJed Brown   /*
248c4762a1bSJed Brown      Determine starting point of each processor's range of
249c4762a1bSJed Brown      grid values.
250c4762a1bSJed Brown   */
2519566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(u,&mybase,&myend));
252c4762a1bSJed Brown 
253c4762a1bSJed Brown   /*
254c4762a1bSJed Brown     Get a pointer to vector data.
255c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
256c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
257c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
258c4762a1bSJed Brown       the array.
259c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
260c4762a1bSJed Brown       C version.  See the users manual for details.
261c4762a1bSJed Brown   */
2629566063dSJacob Faibussowitsch   PetscCall(VecGetArray(u,&u_localptr));
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   /*
265c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
266c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
267c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
268c4762a1bSJed Brown   */
269c4762a1bSJed Brown   for (i=mybase; i<myend; i++) {
270c4762a1bSJed Brown     x = h*(PetscReal)i; /* current location in global grid */
271c4762a1bSJed Brown     u_localptr[i-mybase] = 1.0 + x*x;
272c4762a1bSJed Brown   }
273c4762a1bSJed Brown 
274c4762a1bSJed Brown   /*
275c4762a1bSJed Brown      Restore vector
276c4762a1bSJed Brown   */
2779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(u,&u_localptr));
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   /*
280c4762a1bSJed Brown      Print debugging information if desired
281c4762a1bSJed Brown   */
282c4762a1bSJed Brown   if (appctx->debug) {
2839566063dSJacob Faibussowitsch      PetscCall(PetscPrintf(appctx->comm,"initial guess vector\n"));
2849566063dSJacob Faibussowitsch      PetscCall(VecView(u,PETSC_VIEWER_STDOUT_WORLD));
285c4762a1bSJed Brown   }
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   return 0;
288c4762a1bSJed Brown }
289c4762a1bSJed Brown 
290c4762a1bSJed Brown /* --------------------------------------------------------------------- */
291c4762a1bSJed Brown /*
292c4762a1bSJed Brown   SetBounds - Sets the lower and uper bounds on the interior points
293c4762a1bSJed Brown 
294c4762a1bSJed Brown   Input parameters:
295c4762a1bSJed Brown   xl - vector of lower bounds
296c4762a1bSJed Brown   xu - vector of upper bounds
297c4762a1bSJed Brown   ul - constant lower bound for all variables
298c4762a1bSJed Brown   uh - constant upper bound for all variables
299c4762a1bSJed Brown   appctx - Application context
300c4762a1bSJed Brown  */
301c4762a1bSJed Brown PetscErrorCode SetBounds(Vec xl, Vec xu, PetscScalar ul, PetscScalar uh,AppCtx *appctx)
302c4762a1bSJed Brown {
303c4762a1bSJed Brown   PetscScalar       *l,*u;
304c4762a1bSJed Brown   PetscMPIInt       rank,size;
305c4762a1bSJed Brown   PetscInt          localsize;
306c4762a1bSJed Brown 
307c4762a1bSJed Brown   PetscFunctionBeginUser;
3089566063dSJacob Faibussowitsch   PetscCall(VecSet(xl,ul));
3099566063dSJacob Faibussowitsch   PetscCall(VecSet(xu,uh));
3109566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xl,&localsize));
3119566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xl,&l));
3129566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xu,&u));
313c4762a1bSJed Brown 
3149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(appctx->comm,&rank));
3159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(appctx->comm,&size));
316dd400576SPatrick Sanan   if (rank == 0) {
317c4762a1bSJed Brown     l[0] = -PETSC_INFINITY;
318c4762a1bSJed Brown     u[0] =  PETSC_INFINITY;
319c4762a1bSJed Brown   }
320c4762a1bSJed Brown   if (rank == size-1) {
321c4762a1bSJed Brown     l[localsize-1] = -PETSC_INFINITY;
322c4762a1bSJed Brown     u[localsize-1] = PETSC_INFINITY;
323c4762a1bSJed Brown   }
3249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xl,&l));
3259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xu,&u));
326c4762a1bSJed Brown   PetscFunctionReturn(0);
327c4762a1bSJed Brown }
328c4762a1bSJed Brown 
329c4762a1bSJed Brown /* --------------------------------------------------------------------- */
330c4762a1bSJed Brown /*
331c4762a1bSJed Brown    ExactSolution - Computes the exact solution at a given time.
332c4762a1bSJed Brown 
333c4762a1bSJed Brown    Input Parameters:
334c4762a1bSJed Brown    t - current time
335c4762a1bSJed Brown    solution - vector in which exact solution will be computed
336c4762a1bSJed Brown    appctx - user-defined application context
337c4762a1bSJed Brown 
338c4762a1bSJed Brown    Output Parameter:
339c4762a1bSJed Brown    solution - vector with the newly computed exact solution
340c4762a1bSJed Brown */
341c4762a1bSJed Brown PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
342c4762a1bSJed Brown {
343c4762a1bSJed Brown   PetscScalar    *s_localptr,h = appctx->h,x;
344c4762a1bSJed Brown   PetscInt       i,mybase,myend;
345c4762a1bSJed Brown 
346c4762a1bSJed Brown   /*
347c4762a1bSJed Brown      Determine starting and ending points of each processor's
348c4762a1bSJed Brown      range of grid values
349c4762a1bSJed Brown   */
3509566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(solution,&mybase,&myend));
351c4762a1bSJed Brown 
352c4762a1bSJed Brown   /*
353c4762a1bSJed Brown      Get a pointer to vector data.
354c4762a1bSJed Brown   */
3559566063dSJacob Faibussowitsch   PetscCall(VecGetArray(solution,&s_localptr));
356c4762a1bSJed Brown 
357c4762a1bSJed Brown   /*
358c4762a1bSJed Brown      Simply write the solution directly into the array locations.
359c4762a1bSJed Brown      Alternatively, we could use VecSetValues() or VecSetValuesLocal().
360c4762a1bSJed Brown   */
361c4762a1bSJed Brown   for (i=mybase; i<myend; i++) {
362c4762a1bSJed Brown     x = h*(PetscReal)i;
363c4762a1bSJed Brown     s_localptr[i-mybase] = (t + 1.0)*(1.0 + x*x);
364c4762a1bSJed Brown   }
365c4762a1bSJed Brown 
366c4762a1bSJed Brown   /*
367c4762a1bSJed Brown      Restore vector
368c4762a1bSJed Brown   */
3699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(solution,&s_localptr));
370c4762a1bSJed Brown   return 0;
371c4762a1bSJed Brown }
372c4762a1bSJed Brown /* --------------------------------------------------------------------- */
373c4762a1bSJed Brown /*
374c4762a1bSJed Brown    Monitor - User-provided routine to monitor the solution computed at
375c4762a1bSJed Brown    each timestep.  This example plots the solution and computes the
376c4762a1bSJed Brown    error in two different norms.
377c4762a1bSJed Brown 
378c4762a1bSJed Brown    Input Parameters:
379c4762a1bSJed Brown    ts     - the timestep context
380c4762a1bSJed Brown    step   - the count of the current step (with 0 meaning the
381c4762a1bSJed Brown             initial condition)
382c4762a1bSJed Brown    time   - the current time
383c4762a1bSJed Brown    u      - the solution at this timestep
384c4762a1bSJed Brown    ctx    - the user-provided context for this monitoring routine.
385c4762a1bSJed Brown             In this case we use the application context which contains
386c4762a1bSJed Brown             information about the problem size, workspace and the exact
387c4762a1bSJed Brown             solution.
388c4762a1bSJed Brown */
389c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx)
390c4762a1bSJed Brown {
391c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*) ctx;   /* user-defined application context */
392c4762a1bSJed Brown   PetscReal      en2,en2s,enmax;
393c4762a1bSJed Brown   PetscDraw      draw;
394c4762a1bSJed Brown 
395c4762a1bSJed Brown   /*
396c4762a1bSJed Brown      We use the default X windows viewer
397c4762a1bSJed Brown              PETSC_VIEWER_DRAW_(appctx->comm)
398c4762a1bSJed Brown      that is associated with the current communicator. This saves
399c4762a1bSJed Brown      the effort of calling PetscViewerDrawOpen() to create the window.
400c4762a1bSJed Brown      Note that if we wished to plot several items in separate windows we
401c4762a1bSJed Brown      would create each viewer with PetscViewerDrawOpen() and store them in
402c4762a1bSJed Brown      the application context, appctx.
403c4762a1bSJed Brown 
404c4762a1bSJed Brown      PetscReal buffering makes graphics look better.
405c4762a1bSJed Brown   */
4069566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(PETSC_VIEWER_DRAW_(appctx->comm),0,&draw));
4079566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
4089566063dSJacob Faibussowitsch   PetscCall(VecView(u,PETSC_VIEWER_DRAW_(appctx->comm)));
409c4762a1bSJed Brown 
410c4762a1bSJed Brown   /*
411c4762a1bSJed Brown      Compute the exact solution at this timestep
412c4762a1bSJed Brown   */
4139566063dSJacob Faibussowitsch   PetscCall(ExactSolution(time,appctx->solution,appctx));
414c4762a1bSJed Brown 
415c4762a1bSJed Brown   /*
416c4762a1bSJed Brown      Print debugging information if desired
417c4762a1bSJed Brown   */
418c4762a1bSJed Brown   if (appctx->debug) {
4199566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm,"Computed solution vector\n"));
4209566063dSJacob Faibussowitsch     PetscCall(VecView(u,PETSC_VIEWER_STDOUT_WORLD));
4219566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm,"Exact solution vector\n"));
4229566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD));
423c4762a1bSJed Brown   }
424c4762a1bSJed Brown 
425c4762a1bSJed Brown   /*
426c4762a1bSJed Brown      Compute the 2-norm and max-norm of the error
427c4762a1bSJed Brown   */
4289566063dSJacob Faibussowitsch   PetscCall(VecAXPY(appctx->solution,-1.0,u));
4299566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution,NORM_2,&en2));
430c4762a1bSJed Brown   en2s = PetscSqrtReal(appctx->h)*en2;  /* scale the 2-norm by the grid spacing */
4319566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution,NORM_MAX,&enmax));
432c4762a1bSJed Brown 
433c4762a1bSJed Brown   /*
434c4762a1bSJed Brown      PetscPrintf() causes only the first processor in this
435c4762a1bSJed Brown      communicator to print the timestep information.
436c4762a1bSJed Brown   */
43763a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(appctx->comm,"Timestep %" PetscInt_FMT ": time = %g,2-norm error = %g, max norm error = %g\n",step,(double)time,(double)en2s,(double)enmax));
438c4762a1bSJed Brown 
439c4762a1bSJed Brown   /*
440c4762a1bSJed Brown      Print debugging information if desired
441c4762a1bSJed Brown    */
442c4762a1bSJed Brown   /*  if (appctx->debug) {
4439566063dSJacob Faibussowitsch      PetscCall(PetscPrintf(appctx->comm,"Error vector\n"));
4449566063dSJacob Faibussowitsch      PetscCall(VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD));
445c4762a1bSJed Brown    } */
446c4762a1bSJed Brown   return 0;
447c4762a1bSJed Brown }
448c4762a1bSJed Brown /* --------------------------------------------------------------------- */
449c4762a1bSJed Brown /*
450c4762a1bSJed Brown    RHSFunction - User-provided routine that evalues the right-hand-side
451c4762a1bSJed Brown    function of the ODE.  This routine is set in the main program by
452c4762a1bSJed Brown    calling TSSetRHSFunction().  We compute:
453c4762a1bSJed Brown           global_out = F(global_in)
454c4762a1bSJed Brown 
455c4762a1bSJed Brown    Input Parameters:
456c4762a1bSJed Brown    ts         - timesteping context
457c4762a1bSJed Brown    t          - current time
458c4762a1bSJed Brown    global_in  - vector containing the current iterate
459c4762a1bSJed Brown    ctx        - (optional) user-provided context for function evaluation.
460c4762a1bSJed Brown                 In this case we use the appctx defined above.
461c4762a1bSJed Brown 
462c4762a1bSJed Brown    Output Parameter:
463c4762a1bSJed Brown    global_out - vector containing the newly evaluated function
464c4762a1bSJed Brown */
465c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec global_in,Vec global_out,void *ctx)
466c4762a1bSJed Brown {
467c4762a1bSJed Brown   AppCtx            *appctx   = (AppCtx*) ctx;     /* user-defined application context */
468c4762a1bSJed Brown   DM                da        = appctx->da;        /* distributed array */
469c4762a1bSJed Brown   Vec               local_in  = appctx->u_local;   /* local ghosted input vector */
470c4762a1bSJed Brown   Vec               localwork = appctx->localwork; /* local ghosted work vector */
471c4762a1bSJed Brown   PetscInt          i,localsize;
472c4762a1bSJed Brown   PetscMPIInt       rank,size;
473c4762a1bSJed Brown   PetscScalar       *copyptr,sc;
474c4762a1bSJed Brown   const PetscScalar *localptr;
475c4762a1bSJed Brown 
476c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
477c4762a1bSJed Brown      Get ready for local function computations
478c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
479c4762a1bSJed Brown   /*
480c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
481c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
482c4762a1bSJed Brown      By placing code between these two statements, computations can be
483c4762a1bSJed Brown      done while messages are in transition.
484c4762a1bSJed Brown   */
4859566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in));
4869566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in));
487c4762a1bSJed Brown 
488c4762a1bSJed Brown   /*
489c4762a1bSJed Brown       Access directly the values in our local INPUT work array
490c4762a1bSJed Brown   */
4919566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(local_in,&localptr));
492c4762a1bSJed Brown 
493c4762a1bSJed Brown   /*
494c4762a1bSJed Brown       Access directly the values in our local OUTPUT work array
495c4762a1bSJed Brown   */
4969566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localwork,&copyptr));
497c4762a1bSJed Brown 
498c4762a1bSJed Brown   sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t));
499c4762a1bSJed Brown 
500c4762a1bSJed Brown   /*
501c4762a1bSJed Brown       Evaluate our function on the nodes owned by this processor
502c4762a1bSJed Brown   */
5039566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(local_in,&localsize));
504c4762a1bSJed Brown 
505c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
506c4762a1bSJed Brown      Compute entries for the locally owned part
507c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
508c4762a1bSJed Brown 
509c4762a1bSJed Brown   /*
510c4762a1bSJed Brown      Handle boundary conditions: This is done by using the boundary condition
511c4762a1bSJed Brown         u(t,boundary) = g(t,boundary)
512c4762a1bSJed Brown      for some function g. Now take the derivative with respect to t to obtain
513c4762a1bSJed Brown         u_{t}(t,boundary) = g_{t}(t,boundary)
514c4762a1bSJed Brown 
515c4762a1bSJed Brown      In our case, u(t,0) = t + 1, so that u_{t}(t,0) = 1
516c4762a1bSJed Brown              and  u(t,1) = 2t+ 2, so that u_{t}(t,1) = 2
517c4762a1bSJed Brown   */
5189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(appctx->comm,&rank));
5199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(appctx->comm,&size));
520dd400576SPatrick Sanan   if (rank == 0) copyptr[0] = 1.0;
521c4762a1bSJed Brown   if (rank == size-1) copyptr[localsize-1] = (t < .5) ? 2.0 : 0.0;
522c4762a1bSJed Brown 
523c4762a1bSJed Brown   /*
524c4762a1bSJed Brown      Handle the interior nodes where the PDE is replace by finite
525c4762a1bSJed Brown      difference operators.
526c4762a1bSJed Brown   */
527c4762a1bSJed Brown   for (i=1; i<localsize-1; i++) copyptr[i] =  localptr[i] * sc * (localptr[i+1] + localptr[i-1] - 2.0*localptr[i]);
528c4762a1bSJed Brown 
529c4762a1bSJed Brown   /*
530c4762a1bSJed Brown      Restore vectors
531c4762a1bSJed Brown   */
5329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(local_in,&localptr));
5339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localwork,&copyptr));
534c4762a1bSJed Brown 
535c4762a1bSJed Brown   /*
536c4762a1bSJed Brown      Insert values from the local OUTPUT vector into the global
537c4762a1bSJed Brown      output vector
538c4762a1bSJed Brown   */
5399566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da,localwork,INSERT_VALUES,global_out));
5409566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da,localwork,INSERT_VALUES,global_out));
541c4762a1bSJed Brown 
542c4762a1bSJed Brown   /* Print debugging information if desired */
543c4762a1bSJed Brown   /*  if (appctx->debug) {
5449566063dSJacob Faibussowitsch      PetscCall(PetscPrintf(appctx->comm,"RHS function vector\n"));
5459566063dSJacob Faibussowitsch      PetscCall(VecView(global_out,PETSC_VIEWER_STDOUT_WORLD));
546c4762a1bSJed Brown    } */
547c4762a1bSJed Brown 
548c4762a1bSJed Brown   return 0;
549c4762a1bSJed Brown }
550c4762a1bSJed Brown /* --------------------------------------------------------------------- */
551c4762a1bSJed Brown /*
552c4762a1bSJed Brown    RHSJacobian - User-provided routine to compute the Jacobian of
553c4762a1bSJed Brown    the nonlinear right-hand-side function of the ODE.
554c4762a1bSJed Brown 
555c4762a1bSJed Brown    Input Parameters:
556c4762a1bSJed Brown    ts - the TS context
557c4762a1bSJed Brown    t - current time
558c4762a1bSJed Brown    global_in - global input vector
559c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
560c4762a1bSJed Brown 
561c4762a1bSJed Brown    Output Parameters:
562c4762a1bSJed Brown    AA - Jacobian matrix
563c4762a1bSJed Brown    BB - optionally different preconditioning matrix
564c4762a1bSJed Brown    str - flag indicating matrix structure
565c4762a1bSJed Brown 
566c4762a1bSJed Brown   Notes:
567c4762a1bSJed Brown   RHSJacobian computes entries for the locally owned part of the Jacobian.
568c4762a1bSJed Brown    - Currently, all PETSc parallel matrix formats are partitioned by
569c4762a1bSJed Brown      contiguous chunks of rows across the processors.
570c4762a1bSJed Brown    - Each processor needs to insert only elements that it owns
571c4762a1bSJed Brown      locally (but any non-local elements will be sent to the
572c4762a1bSJed Brown      appropriate processor during matrix assembly).
573c4762a1bSJed Brown    - Always specify global row and columns of matrix entries when
574c4762a1bSJed Brown      using MatSetValues().
575c4762a1bSJed Brown    - Here, we set all entries for a particular row at once.
576c4762a1bSJed Brown    - Note that MatSetValues() uses 0-based row and column numbers
577c4762a1bSJed Brown      in Fortran as well as in C.
578c4762a1bSJed Brown */
579c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec global_in,Mat AA,Mat B,void *ctx)
580c4762a1bSJed Brown {
581c4762a1bSJed Brown   AppCtx            *appctx  = (AppCtx*)ctx;    /* user-defined application context */
582c4762a1bSJed Brown   Vec               local_in = appctx->u_local;   /* local ghosted input vector */
583c4762a1bSJed Brown   DM                da       = appctx->da;        /* distributed array */
584c4762a1bSJed Brown   PetscScalar       v[3],sc;
585c4762a1bSJed Brown   const PetscScalar *localptr;
586c4762a1bSJed Brown   PetscInt          i,mstart,mend,mstarts,mends,idx[3],is;
587c4762a1bSJed Brown 
588c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
589c4762a1bSJed Brown      Get ready for local Jacobian computations
590c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
591c4762a1bSJed Brown   /*
592c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
593c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
594c4762a1bSJed Brown      By placing code between these two statements, computations can be
595c4762a1bSJed Brown      done while messages are in transition.
596c4762a1bSJed Brown   */
5979566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in));
5989566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in));
599c4762a1bSJed Brown 
600c4762a1bSJed Brown   /*
601c4762a1bSJed Brown      Get pointer to vector data
602c4762a1bSJed Brown   */
6039566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(local_in,&localptr));
604c4762a1bSJed Brown 
605c4762a1bSJed Brown   /*
606c4762a1bSJed Brown      Get starting and ending locally owned rows of the matrix
607c4762a1bSJed Brown   */
6089566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(B,&mstarts,&mends));
609c4762a1bSJed Brown   mstart = mstarts; mend = mends;
610c4762a1bSJed Brown 
611c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
612c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
613c4762a1bSJed Brown       - Currently, all PETSc parallel matrix formats are partitioned by
614c4762a1bSJed Brown         contiguous chunks of rows across the processors.
615c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
616c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
617c4762a1bSJed Brown         appropriate processor during matrix assembly).
618c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
619c4762a1bSJed Brown       - We can set matrix entries either using either
620c4762a1bSJed Brown         MatSetValuesLocal() or MatSetValues().
621c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
622c4762a1bSJed Brown 
623c4762a1bSJed Brown   /*
624c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
625c4762a1bSJed Brown   */
626c4762a1bSJed Brown   if (mstart == 0) {
627c4762a1bSJed Brown     v[0] = 0.0;
6289566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B,1,&mstart,1,&mstart,v,INSERT_VALUES));
629c4762a1bSJed Brown     mstart++;
630c4762a1bSJed Brown   }
631c4762a1bSJed Brown   if (mend == appctx->m) {
632c4762a1bSJed Brown     mend--;
633c4762a1bSJed Brown     v[0] = 0.0;
6349566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B,1,&mend,1,&mend,v,INSERT_VALUES));
635c4762a1bSJed Brown   }
636c4762a1bSJed Brown 
637c4762a1bSJed Brown   /*
638c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
639c4762a1bSJed Brown      matrix one row at a time.
640c4762a1bSJed Brown   */
641c4762a1bSJed Brown   sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t));
642c4762a1bSJed Brown   for (i=mstart; i<mend; i++) {
643c4762a1bSJed Brown     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
644c4762a1bSJed Brown     is     = i - mstart + 1;
645c4762a1bSJed Brown     v[0]   = sc*localptr[is];
646c4762a1bSJed Brown     v[1]   = sc*(localptr[is+1] + localptr[is-1] - 4.0*localptr[is]);
647c4762a1bSJed Brown     v[2]   = sc*localptr[is];
6489566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B,1,&i,3,idx,v,INSERT_VALUES));
649c4762a1bSJed Brown   }
650c4762a1bSJed Brown 
651c4762a1bSJed Brown   /*
652c4762a1bSJed Brown      Restore vector
653c4762a1bSJed Brown   */
6549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(local_in,&localptr));
655c4762a1bSJed Brown 
656c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
657c4762a1bSJed Brown      Complete the matrix assembly process and set some options
658c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
659c4762a1bSJed Brown   /*
660c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
661c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
662c4762a1bSJed Brown      Computations can be done while messages are in transition
663c4762a1bSJed Brown      by placing code between these two statements.
664c4762a1bSJed Brown   */
6659566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
6669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
667c4762a1bSJed Brown 
668c4762a1bSJed Brown   /*
669c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
670c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
671c4762a1bSJed Brown   */
6729566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
673c4762a1bSJed Brown 
674c4762a1bSJed Brown   return 0;
675c4762a1bSJed Brown }
676c4762a1bSJed Brown 
677c4762a1bSJed Brown /*TEST
678c4762a1bSJed Brown 
679c4762a1bSJed Brown     test:
680c4762a1bSJed Brown       args: -snes_type vinewtonrsls -ts_type glee -mymonitor -ts_max_steps 10 -nox
681c4762a1bSJed Brown       requires: !single
682c4762a1bSJed Brown 
683c4762a1bSJed Brown TEST*/
684