xref: /petsc/src/ts/tests/ex12.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] ="Tests PetscObjectSetOptions() for TS object\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown    Concepts: TS^time-dependent nonlinear problems
6c4762a1bSJed Brown    Processors: n
7c4762a1bSJed Brown */
8c4762a1bSJed Brown 
9c4762a1bSJed Brown /* ------------------------------------------------------------------------
10c4762a1bSJed Brown 
11c4762a1bSJed Brown    This program solves the PDE
12c4762a1bSJed Brown 
13c4762a1bSJed Brown                u * u_xx
14c4762a1bSJed Brown          u_t = ---------
15c4762a1bSJed Brown                2*(t+1)^2
16c4762a1bSJed Brown 
17c4762a1bSJed Brown     on the domain 0 <= x <= 1, with boundary conditions
18c4762a1bSJed Brown          u(t,0) = t + 1,  u(t,1) = 2*t + 2,
19c4762a1bSJed Brown     and initial condition
20c4762a1bSJed Brown          u(0,x) = 1 + x*x.
21c4762a1bSJed Brown 
22c4762a1bSJed Brown     The exact solution is:
23c4762a1bSJed Brown          u(t,x) = (1 + x*x) * (1 + t)
24c4762a1bSJed Brown 
25c4762a1bSJed Brown     Note that since the solution is linear in time and quadratic in x,
26c4762a1bSJed Brown     the finite difference scheme actually computes the "exact" solution.
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 } AppCtx;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown /*
60c4762a1bSJed Brown    User-defined routines, provided below.
61c4762a1bSJed Brown */
62c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec,AppCtx*);
63c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
64c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
65c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);
66c4762a1bSJed Brown 
67c4762a1bSJed Brown int main(int argc,char **argv)
68c4762a1bSJed Brown {
69c4762a1bSJed Brown   AppCtx         appctx;                 /* user-defined application context */
70c4762a1bSJed Brown   TS             ts;                     /* timestepping context */
71c4762a1bSJed Brown   Mat            A;                      /* Jacobian matrix data structure */
72c4762a1bSJed Brown   Vec            u;                      /* approximate solution vector */
73c4762a1bSJed Brown   PetscInt       time_steps_max = 100;  /* default max timesteps */
74c4762a1bSJed Brown   PetscErrorCode ierr;
75c4762a1bSJed Brown   PetscReal      dt;
76c4762a1bSJed Brown   PetscReal      time_total_max = 100.0; /* default max total time */
77c4762a1bSJed Brown   PetscOptions   options,optionscopy;
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80c4762a1bSJed Brown      Initialize program and set problem parameters
81c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
84c4762a1bSJed Brown 
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsCreate(&options));
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsSetValue(options,"-ts_monitor","ascii"));
87*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsSetValue(options,"-snes_monitor","ascii"));
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsSetValue(options,"-ksp_monitor","ascii"));
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   appctx.comm = PETSC_COMM_WORLD;
91c4762a1bSJed Brown   appctx.m    = 60;
92c4762a1bSJed Brown 
93*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(options,NULL,"-M",&appctx.m,NULL));
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   appctx.h    = 1.0/(appctx.m-1.0);
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98c4762a1bSJed Brown      Create vector data structures
99c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /*
102c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
103c4762a1bSJed Brown      and to set up the ghost point communication pattern.  There are M
104c4762a1bSJed Brown      total grid values spread equally among all the processors.
105c4762a1bSJed Brown   */
106*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,appctx.m,1,1,NULL,&appctx.da));
107*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptions((PetscObject)appctx.da,options));
108*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(appctx.da));
109*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(appctx.da));
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /*
112c4762a1bSJed Brown      Extract global and local vectors from DMDA; we use these to store the
113c4762a1bSJed Brown      approximate solution.  Then duplicate these for remaining vectors that
114c4762a1bSJed Brown      have the same types.
115c4762a1bSJed Brown   */
116*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(appctx.da,&u));
117*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateLocalVector(appctx.da,&appctx.u_local));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /*
120c4762a1bSJed Brown      Create local work vector for use in evaluating right-hand-side function;
121c4762a1bSJed Brown      create global work vector for storing exact solution.
122c4762a1bSJed Brown   */
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(appctx.u_local,&appctx.localwork));
124*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(u,&appctx.solution));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127c4762a1bSJed Brown      Create timestepping solver context; set callback routine for
128c4762a1bSJed Brown      right-hand-side function evaluation.
129c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130c4762a1bSJed Brown 
131*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
132*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptions((PetscObject)ts,options));
133*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
134*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunction,&appctx));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137c4762a1bSJed Brown      For nonlinear problems, the user can provide a Jacobian evaluation
138c4762a1bSJed Brown      routine (or use a finite differencing approximation).
139c4762a1bSJed Brown 
140c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine.
141c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142c4762a1bSJed Brown 
143*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
144*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptions((PetscObject)A,options));
145*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,appctx.m,appctx.m));
146*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
148*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobian(ts,A,A,RHSJacobian,&appctx));
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151c4762a1bSJed Brown      Set solution vector and initial timestep
152c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   dt   = appctx.h/2.0;
155*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,dt));
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158c4762a1bSJed Brown      Customize timestepping solver:
159c4762a1bSJed Brown        - Set the solution method to be the Backward Euler method.
160c4762a1bSJed Brown        - Set timestepping duration info
161c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
162c4762a1bSJed Brown      For example,
163c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
164c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
165c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166c4762a1bSJed Brown 
167*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSBEULER));
168*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxSteps(ts,time_steps_max));
169*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,time_total_max));
170*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
171*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
174c4762a1bSJed Brown      Solve the problem
175c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   /*
178c4762a1bSJed Brown      Evaluate initial conditions
179c4762a1bSJed Brown   */
180*5f80ce2aSJacob Faibussowitsch   CHKERRQ(InitialConditions(u,&appctx));
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   /*
183c4762a1bSJed Brown      Run the timestepping solver
184c4762a1bSJed Brown   */
185*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,u));
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
189c4762a1bSJed Brown      are no longer needed.
190c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
191c4762a1bSJed Brown 
192*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetOptions((PetscObject)ts,&optionscopy));
1933c633725SBarry Smith   PetscCheck(options == optionscopy,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"PetscObjectGetOptions() failed");
194c4762a1bSJed Brown 
195*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
196*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
197*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
198*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&appctx.da));
199*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&appctx.localwork));
200*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&appctx.solution));
201*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&appctx.u_local));
202*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsDestroy(&options));
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   /*
205c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
206c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
207c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
208c4762a1bSJed Brown          options are chosen (e.g., -log_view).
209c4762a1bSJed Brown   */
210c4762a1bSJed Brown   ierr = PetscFinalize();
211c4762a1bSJed Brown   return ierr;
212c4762a1bSJed Brown }
213c4762a1bSJed Brown /* --------------------------------------------------------------------- */
214c4762a1bSJed Brown /*
215c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
216c4762a1bSJed Brown 
217c4762a1bSJed Brown    Input Parameters:
218c4762a1bSJed Brown    u - uninitialized solution vector (global)
219c4762a1bSJed Brown    appctx - user-defined application context
220c4762a1bSJed Brown 
221c4762a1bSJed Brown    Output Parameter:
222c4762a1bSJed Brown    u - vector with solution at initial time (global)
223c4762a1bSJed Brown */
224c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
225c4762a1bSJed Brown {
226c4762a1bSJed Brown   PetscScalar    *u_localptr,h = appctx->h,x;
227c4762a1bSJed Brown   PetscInt       i,mybase,myend;
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   /*
230c4762a1bSJed Brown      Determine starting point of each processor's range of
231c4762a1bSJed Brown      grid values.
232c4762a1bSJed Brown   */
233*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(u,&mybase,&myend));
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   /*
236c4762a1bSJed Brown     Get a pointer to vector data.
237c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
238c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
239c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
240c4762a1bSJed Brown       the array.
241c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
242c4762a1bSJed Brown       C version.  See the users manual for details.
243c4762a1bSJed Brown   */
244*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(u,&u_localptr));
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   /*
247c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
248c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
249c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
250c4762a1bSJed Brown   */
251c4762a1bSJed Brown   for (i=mybase; i<myend; i++) {
252c4762a1bSJed Brown     x = h*(PetscReal)i; /* current location in global grid */
253c4762a1bSJed Brown     u_localptr[i-mybase] = 1.0 + x*x;
254c4762a1bSJed Brown   }
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   /*
257c4762a1bSJed Brown      Restore vector
258c4762a1bSJed Brown   */
259*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(u,&u_localptr));
260c4762a1bSJed Brown 
261c4762a1bSJed Brown   return 0;
262c4762a1bSJed Brown }
263c4762a1bSJed Brown /* --------------------------------------------------------------------- */
264c4762a1bSJed Brown /*
265c4762a1bSJed Brown    ExactSolution - Computes the exact solution at a given time.
266c4762a1bSJed Brown 
267c4762a1bSJed Brown    Input Parameters:
268c4762a1bSJed Brown    t - current time
269c4762a1bSJed Brown    solution - vector in which exact solution will be computed
270c4762a1bSJed Brown    appctx - user-defined application context
271c4762a1bSJed Brown 
272c4762a1bSJed Brown    Output Parameter:
273c4762a1bSJed Brown    solution - vector with the newly computed exact solution
274c4762a1bSJed Brown */
275c4762a1bSJed Brown PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
276c4762a1bSJed Brown {
277c4762a1bSJed Brown   PetscScalar    *s_localptr,h = appctx->h,x;
278c4762a1bSJed Brown   PetscInt       i,mybase,myend;
279c4762a1bSJed Brown 
280c4762a1bSJed Brown   /*
281c4762a1bSJed Brown      Determine starting and ending points of each processor's
282c4762a1bSJed Brown      range of grid values
283c4762a1bSJed Brown   */
284*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(solution,&mybase,&myend));
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   /*
287c4762a1bSJed Brown      Get a pointer to vector data.
288c4762a1bSJed Brown   */
289*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(solution,&s_localptr));
290c4762a1bSJed Brown 
291c4762a1bSJed Brown   /*
292c4762a1bSJed Brown      Simply write the solution directly into the array locations.
293c4762a1bSJed Brown      Alternatively, we could use VecSetValues() or VecSetValuesLocal().
294c4762a1bSJed Brown   */
295c4762a1bSJed Brown   for (i=mybase; i<myend; i++) {
296c4762a1bSJed Brown     x = h*(PetscReal)i;
297c4762a1bSJed Brown     s_localptr[i-mybase] = (t + 1.0)*(1.0 + x*x);
298c4762a1bSJed Brown   }
299c4762a1bSJed Brown 
300c4762a1bSJed Brown   /*
301c4762a1bSJed Brown      Restore vector
302c4762a1bSJed Brown   */
303*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(solution,&s_localptr));
304c4762a1bSJed Brown   return 0;
305c4762a1bSJed Brown }
306c4762a1bSJed Brown /* --------------------------------------------------------------------- */
307c4762a1bSJed Brown /*
308c4762a1bSJed Brown    RHSFunction - User-provided routine that evalues the right-hand-side
309c4762a1bSJed Brown    function of the ODE.  This routine is set in the main program by
310c4762a1bSJed Brown    calling TSSetRHSFunction().  We compute:
311c4762a1bSJed Brown           global_out = F(global_in)
312c4762a1bSJed Brown 
313c4762a1bSJed Brown    Input Parameters:
314c4762a1bSJed Brown    ts         - timesteping context
315c4762a1bSJed Brown    t          - current time
316c4762a1bSJed Brown    global_in  - vector containing the current iterate
317c4762a1bSJed Brown    ctx        - (optional) user-provided context for function evaluation.
318c4762a1bSJed Brown                 In this case we use the appctx defined above.
319c4762a1bSJed Brown 
320c4762a1bSJed Brown    Output Parameter:
321c4762a1bSJed Brown    global_out - vector containing the newly evaluated function
322c4762a1bSJed Brown */
323c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec global_in,Vec global_out,void *ctx)
324c4762a1bSJed Brown {
325c4762a1bSJed Brown   AppCtx            *appctx   = (AppCtx*) ctx;     /* user-defined application context */
326c4762a1bSJed Brown   DM                da        = appctx->da;        /* distributed array */
327c4762a1bSJed Brown   Vec               local_in  = appctx->u_local;   /* local ghosted input vector */
328c4762a1bSJed Brown   Vec               localwork = appctx->localwork; /* local ghosted work vector */
329c4762a1bSJed Brown   PetscInt          i,localsize;
330c4762a1bSJed Brown   PetscMPIInt       rank,size;
331c4762a1bSJed Brown   PetscScalar       *copyptr,sc;
332c4762a1bSJed Brown   const PetscScalar *localptr;
333c4762a1bSJed Brown 
334c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
335c4762a1bSJed Brown      Get ready for local function computations
336c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
337c4762a1bSJed Brown   /*
338c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
339c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
340c4762a1bSJed Brown      By placing code between these two statements, computations can be
341c4762a1bSJed Brown      done while messages are in transition.
342c4762a1bSJed Brown   */
343*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in));
344*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in));
345c4762a1bSJed Brown 
346c4762a1bSJed Brown   /*
347c4762a1bSJed Brown       Access directly the values in our local INPUT work array
348c4762a1bSJed Brown   */
349*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(local_in,&localptr));
350c4762a1bSJed Brown 
351c4762a1bSJed Brown   /*
352c4762a1bSJed Brown       Access directly the values in our local OUTPUT work array
353c4762a1bSJed Brown   */
354*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(localwork,&copyptr));
355c4762a1bSJed Brown 
356c4762a1bSJed Brown   sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t));
357c4762a1bSJed Brown 
358c4762a1bSJed Brown   /*
359c4762a1bSJed Brown       Evaluate our function on the nodes owned by this processor
360c4762a1bSJed Brown   */
361*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(local_in,&localsize));
362c4762a1bSJed Brown 
363c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
364c4762a1bSJed Brown      Compute entries for the locally owned part
365c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
366c4762a1bSJed Brown 
367c4762a1bSJed Brown   /*
368c4762a1bSJed Brown      Handle boundary conditions: This is done by using the boundary condition
369c4762a1bSJed Brown         u(t,boundary) = g(t,boundary)
370c4762a1bSJed Brown      for some function g. Now take the derivative with respect to t to obtain
371c4762a1bSJed Brown         u_{t}(t,boundary) = g_{t}(t,boundary)
372c4762a1bSJed Brown 
373c4762a1bSJed Brown      In our case, u(t,0) = t + 1, so that u_{t}(t,0) = 1
374c4762a1bSJed Brown              and  u(t,1) = 2t+ 2, so that u_{t}(t,1) = 2
375c4762a1bSJed Brown   */
376*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(appctx->comm,&rank));
377*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(appctx->comm,&size));
378dd400576SPatrick Sanan   if (rank == 0)          copyptr[0]           = 1.0;
379c4762a1bSJed Brown   if (rank == size-1) copyptr[localsize-1] = 2.0;
380c4762a1bSJed Brown 
381c4762a1bSJed Brown   /*
382c4762a1bSJed Brown      Handle the interior nodes where the PDE is replace by finite
383c4762a1bSJed Brown      difference operators.
384c4762a1bSJed Brown   */
385c4762a1bSJed Brown   for (i=1; i<localsize-1; i++) copyptr[i] =  localptr[i] * sc * (localptr[i+1] + localptr[i-1] - 2.0*localptr[i]);
386c4762a1bSJed Brown 
387c4762a1bSJed Brown   /*
388c4762a1bSJed Brown      Restore vectors
389c4762a1bSJed Brown   */
390*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(local_in,&localptr));
391*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(localwork,&copyptr));
392c4762a1bSJed Brown 
393c4762a1bSJed Brown   /*
394c4762a1bSJed Brown      Insert values from the local OUTPUT vector into the global
395c4762a1bSJed Brown      output vector
396c4762a1bSJed Brown   */
397*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalBegin(da,localwork,INSERT_VALUES,global_out));
398*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalEnd(da,localwork,INSERT_VALUES,global_out));
399c4762a1bSJed Brown 
400c4762a1bSJed Brown   return 0;
401c4762a1bSJed Brown }
402c4762a1bSJed Brown /* --------------------------------------------------------------------- */
403c4762a1bSJed Brown /*
404c4762a1bSJed Brown    RHSJacobian - User-provided routine to compute the Jacobian of
405c4762a1bSJed Brown    the nonlinear right-hand-side function of the ODE.
406c4762a1bSJed Brown 
407c4762a1bSJed Brown    Input Parameters:
408c4762a1bSJed Brown    ts - the TS context
409c4762a1bSJed Brown    t - current time
410c4762a1bSJed Brown    global_in - global input vector
411c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
412c4762a1bSJed Brown 
413c4762a1bSJed Brown    Output Parameters:
414c4762a1bSJed Brown    AA - Jacobian matrix
415c4762a1bSJed Brown    BB - optionally different preconditioning matrix
416c4762a1bSJed Brown    str - flag indicating matrix structure
417c4762a1bSJed Brown 
418c4762a1bSJed Brown   Notes:
419c4762a1bSJed Brown   RHSJacobian computes entries for the locally owned part of the Jacobian.
420c4762a1bSJed Brown    - Currently, all PETSc parallel matrix formats are partitioned by
421c4762a1bSJed Brown      contiguous chunks of rows across the processors.
422c4762a1bSJed Brown    - Each processor needs to insert only elements that it owns
423c4762a1bSJed Brown      locally (but any non-local elements will be sent to the
424c4762a1bSJed Brown      appropriate processor during matrix assembly).
425c4762a1bSJed Brown    - Always specify global row and columns of matrix entries when
426c4762a1bSJed Brown      using MatSetValues().
427c4762a1bSJed Brown    - Here, we set all entries for a particular row at once.
428c4762a1bSJed Brown    - Note that MatSetValues() uses 0-based row and column numbers
429c4762a1bSJed Brown      in Fortran as well as in C.
430c4762a1bSJed Brown */
431c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec global_in,Mat AA,Mat BB,void *ctx)
432c4762a1bSJed Brown {
433c4762a1bSJed Brown   AppCtx            *appctx  = (AppCtx*)ctx;    /* user-defined application context */
434c4762a1bSJed Brown   Vec               local_in = appctx->u_local;   /* local ghosted input vector */
435c4762a1bSJed Brown   DM                da       = appctx->da;        /* distributed array */
436c4762a1bSJed Brown   PetscScalar       v[3],sc;
437c4762a1bSJed Brown   const PetscScalar *localptr;
438c4762a1bSJed Brown   PetscInt          i,mstart,mend,mstarts,mends,idx[3],is;
439c4762a1bSJed Brown 
440c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
441c4762a1bSJed Brown      Get ready for local Jacobian computations
442c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
443c4762a1bSJed Brown   /*
444c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
445c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
446c4762a1bSJed Brown      By placing code between these two statements, computations can be
447c4762a1bSJed Brown      done while messages are in transition.
448c4762a1bSJed Brown   */
449*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in));
450*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in));
451c4762a1bSJed Brown 
452c4762a1bSJed Brown   /*
453c4762a1bSJed Brown      Get pointer to vector data
454c4762a1bSJed Brown   */
455*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(local_in,&localptr));
456c4762a1bSJed Brown 
457c4762a1bSJed Brown   /*
458c4762a1bSJed Brown      Get starting and ending locally owned rows of the matrix
459c4762a1bSJed Brown   */
460*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(BB,&mstarts,&mends));
461c4762a1bSJed Brown   mstart = mstarts; mend = mends;
462c4762a1bSJed Brown 
463c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
464c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
465c4762a1bSJed Brown       - Currently, all PETSc parallel matrix formats are partitioned by
466c4762a1bSJed Brown         contiguous chunks of rows across the processors.
467c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
468c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
469c4762a1bSJed Brown         appropriate processor during matrix assembly).
470c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
471c4762a1bSJed Brown       - We can set matrix entries either using either
472c4762a1bSJed Brown         MatSetValuesLocal() or MatSetValues().
473c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
474c4762a1bSJed Brown 
475c4762a1bSJed Brown   /*
476c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
477c4762a1bSJed Brown   */
478c4762a1bSJed Brown   if (mstart == 0) {
479c4762a1bSJed Brown     v[0] = 0.0;
480*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(BB,1,&mstart,1,&mstart,v,INSERT_VALUES));
481c4762a1bSJed Brown     mstart++;
482c4762a1bSJed Brown   }
483c4762a1bSJed Brown   if (mend == appctx->m) {
484c4762a1bSJed Brown     mend--;
485c4762a1bSJed Brown     v[0] = 0.0;
486*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(BB,1,&mend,1,&mend,v,INSERT_VALUES));
487c4762a1bSJed Brown   }
488c4762a1bSJed Brown 
489c4762a1bSJed Brown   /*
490c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
491c4762a1bSJed Brown      matrix one row at a time.
492c4762a1bSJed Brown   */
493c4762a1bSJed Brown   sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t));
494c4762a1bSJed Brown   for (i=mstart; i<mend; i++) {
495c4762a1bSJed Brown     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
496c4762a1bSJed Brown     is     = i - mstart + 1;
497c4762a1bSJed Brown     v[0]   = sc*localptr[is];
498c4762a1bSJed Brown     v[1]   = sc*(localptr[is+1] + localptr[is-1] - 4.0*localptr[is]);
499c4762a1bSJed Brown     v[2]   = sc*localptr[is];
500*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES));
501c4762a1bSJed Brown   }
502c4762a1bSJed Brown 
503c4762a1bSJed Brown   /*
504c4762a1bSJed Brown      Restore vector
505c4762a1bSJed Brown   */
506*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(local_in,&localptr));
507c4762a1bSJed Brown 
508c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
509c4762a1bSJed Brown      Complete the matrix assembly process and set some options
510c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
511c4762a1bSJed Brown   /*
512c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
513c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
514c4762a1bSJed Brown      Computations can be done while messages are in transition
515c4762a1bSJed Brown      by placing code between these two statements.
516c4762a1bSJed Brown   */
517*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(BB,MAT_FINAL_ASSEMBLY));
518*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(BB,MAT_FINAL_ASSEMBLY));
519c4762a1bSJed Brown   if (BB != AA) {
520*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(AA,MAT_FINAL_ASSEMBLY));
521*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(AA,MAT_FINAL_ASSEMBLY));
522c4762a1bSJed Brown   }
523c4762a1bSJed Brown 
524c4762a1bSJed Brown   /*
525c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
526c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
527c4762a1bSJed Brown   */
528*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(BB,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
529c4762a1bSJed Brown 
530c4762a1bSJed Brown   return 0;
531c4762a1bSJed Brown }
532c4762a1bSJed Brown 
533c4762a1bSJed Brown /*TEST
534c4762a1bSJed Brown 
535c4762a1bSJed Brown     test:
536c4762a1bSJed Brown       requires: !single
537c4762a1bSJed Brown 
538c4762a1bSJed Brown TEST*/
539