1c4762a1bSJed Brown static char help[] = "Nonlinear, time-dependent PDE in 2d.\n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown Solves the equation
4c4762a1bSJed Brown
5c4762a1bSJed Brown u_tt - \Delta u = 0
6c4762a1bSJed Brown
7c4762a1bSJed Brown which we split into two first-order systems
8c4762a1bSJed Brown
9c4762a1bSJed Brown u_t - v = 0 so that F(u,v).u = v
10c4762a1bSJed Brown v_t - \Delta u = 0 so that F(u,v).v = Delta u
11c4762a1bSJed Brown
12c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
13c4762a1bSJed Brown Include "petscts.h" so that we can use SNES solvers. Note that this
14c4762a1bSJed Brown file automatically includes:
15c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors
16c4762a1bSJed Brown petscmat.h - matrices
17c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods
18c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners
19c4762a1bSJed Brown petscksp.h - linear solvers
20c4762a1bSJed Brown */
21c4762a1bSJed Brown #include <petscdm.h>
22c4762a1bSJed Brown #include <petscdmda.h>
23c4762a1bSJed Brown #include <petscts.h>
24c4762a1bSJed Brown
25c4762a1bSJed Brown /*
26c4762a1bSJed Brown User-defined routines
27c4762a1bSJed Brown */
28c4762a1bSJed Brown extern PetscErrorCode FormFunction(TS, PetscReal, Vec, Vec, void *), FormInitialSolution(DM, Vec);
29c4762a1bSJed Brown extern PetscErrorCode MyTSMonitor(TS, PetscInt, PetscReal, Vec, void *);
30c4762a1bSJed Brown extern PetscErrorCode MySNESMonitor(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
31c4762a1bSJed Brown
main(int argc,char ** argv)32d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
33d71ae5a4SJacob Faibussowitsch {
34c4762a1bSJed Brown TS ts; /* nonlinear solver */
35c4762a1bSJed Brown Vec x, r; /* solution, residual vectors */
36c4762a1bSJed Brown PetscInt steps; /* iterations for convergence */
37c4762a1bSJed Brown DM da;
38c4762a1bSJed Brown PetscReal ftime;
39c4762a1bSJed Brown SNES ts_snes;
40c4762a1bSJed Brown PetscBool usemonitor = PETSC_TRUE;
41c4762a1bSJed Brown PetscViewerAndFormat *vf;
42c4762a1bSJed Brown
43c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
44c4762a1bSJed Brown Initialize program
45c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
46327415f7SBarry Smith PetscFunctionBeginUser;
47c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
489566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-usemonitor", &usemonitor, NULL));
49c4762a1bSJed Brown
50c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors
52c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
539566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 8, 8, PETSC_DECIDE, PETSC_DECIDE, 2, 1, NULL, NULL, &da));
549566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da));
559566063dSJacob Faibussowitsch PetscCall(DMSetUp(da));
569566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "u"));
579566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 1, "v"));
58c4762a1bSJed Brown
59c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining
61c4762a1bSJed Brown vectors that are the same types
62c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
639566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x));
649566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r));
65c4762a1bSJed Brown
66c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67c4762a1bSJed Brown Create timestepping solver context
68c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
699566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
709566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da));
719566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
729566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, FormFunction, da));
73c4762a1bSJed Brown
749566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 1.0));
751baa6e33SBarry Smith if (usemonitor) PetscCall(TSMonitorSet(ts, MyTSMonitor, 0, 0));
76c4762a1bSJed Brown
77c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78c4762a1bSJed Brown Customize nonlinear solver
79c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
809566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER));
819566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &ts_snes));
82c4762a1bSJed Brown if (usemonitor) {
839566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, &vf));
8449abdd8aSBarry Smith PetscCall(SNESMonitorSet(ts_snes, (PetscErrorCode (*)(SNES, PetscInt, PetscReal, void *))MySNESMonitor, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy));
85c4762a1bSJed Brown }
86c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87c4762a1bSJed Brown Set initial conditions
88c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
899566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(da, x));
909566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .0001));
919566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
929566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x));
93c4762a1bSJed Brown
94c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95c4762a1bSJed Brown Set runtime options
96c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
979566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
98c4762a1bSJed Brown
99c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100c4762a1bSJed Brown Solve nonlinear system
101c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1029566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x));
1039566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime));
1049566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps));
1059566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(x, NULL, "-final_sol"));
106c4762a1bSJed Brown
107c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they
109c4762a1bSJed Brown are no longer needed.
110c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r));
1139566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
1149566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da));
115c4762a1bSJed Brown
1169566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
117b122ec5aSJacob Faibussowitsch return 0;
118c4762a1bSJed Brown }
119c4762a1bSJed Brown /* ------------------------------------------------------------------- */
120c4762a1bSJed Brown /*
121c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x).
122c4762a1bSJed Brown
123c4762a1bSJed Brown Input Parameters:
124c4762a1bSJed Brown . ts - the TS context
125c4762a1bSJed Brown . X - input vector
126c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction()
127c4762a1bSJed Brown
128c4762a1bSJed Brown Output Parameter:
129c4762a1bSJed Brown . F - function vector
130c4762a1bSJed Brown */
FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void * ptr)131d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(TS ts, PetscReal ftime, Vec X, Vec F, void *ptr)
132d71ae5a4SJacob Faibussowitsch {
133c4762a1bSJed Brown DM da = (DM)ptr;
134c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym;
135c4762a1bSJed Brown PetscReal hx, hy, /*hxdhy,hydhx,*/ sx, sy;
136c4762a1bSJed Brown PetscScalar u, uxx, uyy, v, ***x, ***f;
137c4762a1bSJed Brown Vec localX;
138c4762a1bSJed Brown
139c4762a1bSJed Brown PetscFunctionBeginUser;
1409566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX));
1419566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
142c4762a1bSJed Brown
1439371c9d4SSatish Balay hx = 1.0 / (PetscReal)(Mx - 1);
1449371c9d4SSatish Balay sx = 1.0 / (hx * hx);
1459371c9d4SSatish Balay hy = 1.0 / (PetscReal)(My - 1);
1469371c9d4SSatish Balay sy = 1.0 / (hy * hy);
147c4762a1bSJed Brown /*hxdhy = hx/hy;*/
148c4762a1bSJed Brown /*hydhx = hy/hx;*/
149c4762a1bSJed Brown
150c4762a1bSJed Brown /*
151c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process
152c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
153c4762a1bSJed Brown By placing code between these two statements, computations can be
154c4762a1bSJed Brown done while messages are in transition.
155c4762a1bSJed Brown */
1569566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
1579566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
158c4762a1bSJed Brown
159c4762a1bSJed Brown /*
160c4762a1bSJed Brown Get pointers to vector data
161c4762a1bSJed Brown */
1629566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOF(da, localX, &x));
1639566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOF(da, F, &f));
164c4762a1bSJed Brown
165c4762a1bSJed Brown /*
166c4762a1bSJed Brown Get local grid boundaries
167c4762a1bSJed Brown */
1689566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
169c4762a1bSJed Brown
170c4762a1bSJed Brown /*
171c4762a1bSJed Brown Compute function over the locally owned part of the grid
172c4762a1bSJed Brown */
173c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
174c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
175c4762a1bSJed Brown if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) {
176c4762a1bSJed Brown f[j][i][0] = x[j][i][0];
177c4762a1bSJed Brown f[j][i][1] = x[j][i][1];
178c4762a1bSJed Brown continue;
179c4762a1bSJed Brown }
180c4762a1bSJed Brown u = x[j][i][0];
181c4762a1bSJed Brown v = x[j][i][1];
182c4762a1bSJed Brown uxx = (-2.0 * u + x[j][i - 1][0] + x[j][i + 1][0]) * sx;
183c4762a1bSJed Brown uyy = (-2.0 * u + x[j - 1][i][0] + x[j + 1][i][0]) * sy;
184c4762a1bSJed Brown f[j][i][0] = v;
185c4762a1bSJed Brown f[j][i][1] = uxx + uyy;
186c4762a1bSJed Brown }
187c4762a1bSJed Brown }
188c4762a1bSJed Brown
189c4762a1bSJed Brown /*
190c4762a1bSJed Brown Restore vectors
191c4762a1bSJed Brown */
1929566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOF(da, localX, &x));
1939566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOF(da, F, &f));
1949566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX));
1959566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(11.0 * ym * xm));
1963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
197c4762a1bSJed Brown }
198c4762a1bSJed Brown
199c4762a1bSJed Brown /* ------------------------------------------------------------------- */
FormInitialSolution(DM da,Vec U)200d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialSolution(DM da, Vec U)
201d71ae5a4SJacob Faibussowitsch {
202c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My;
203c4762a1bSJed Brown PetscScalar ***u;
204c4762a1bSJed Brown PetscReal hx, hy, x, y, r;
205c4762a1bSJed Brown
206c4762a1bSJed Brown PetscFunctionBeginUser;
2079566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
208c4762a1bSJed Brown
209c4762a1bSJed Brown hx = 1.0 / (PetscReal)(Mx - 1);
210c4762a1bSJed Brown hy = 1.0 / (PetscReal)(My - 1);
211c4762a1bSJed Brown
212c4762a1bSJed Brown /*
213c4762a1bSJed Brown Get pointers to vector data
214c4762a1bSJed Brown */
2159566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOF(da, U, &u));
216c4762a1bSJed Brown
217c4762a1bSJed Brown /*
218c4762a1bSJed Brown Get local grid boundaries
219c4762a1bSJed Brown */
2209566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
221c4762a1bSJed Brown
222c4762a1bSJed Brown /*
223c4762a1bSJed Brown Compute function over the locally owned part of the grid
224c4762a1bSJed Brown */
225c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
226c4762a1bSJed Brown y = j * hy;
227c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
228c4762a1bSJed Brown x = i * hx;
229c4762a1bSJed Brown r = PetscSqrtReal((x - .5) * (x - .5) + (y - .5) * (y - .5));
230c4762a1bSJed Brown if (r < .125) {
231c4762a1bSJed Brown u[j][i][0] = PetscExpReal(-30.0 * r * r * r);
232c4762a1bSJed Brown u[j][i][1] = 0.0;
233c4762a1bSJed Brown } else {
234c4762a1bSJed Brown u[j][i][0] = 0.0;
235c4762a1bSJed Brown u[j][i][1] = 0.0;
236c4762a1bSJed Brown }
237c4762a1bSJed Brown }
238c4762a1bSJed Brown }
239c4762a1bSJed Brown
240c4762a1bSJed Brown /*
241c4762a1bSJed Brown Restore vectors
242c4762a1bSJed Brown */
2439566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOF(da, U, &u));
2443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
245c4762a1bSJed Brown }
246c4762a1bSJed Brown
MyTSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscCtx ctx)247*2a8381b2SBarry Smith PetscErrorCode MyTSMonitor(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscCtx ctx)
248d71ae5a4SJacob Faibussowitsch {
249c4762a1bSJed Brown PetscReal norm;
250c4762a1bSJed Brown MPI_Comm comm;
251c4762a1bSJed Brown
252c4762a1bSJed Brown PetscFunctionBeginUser;
2539566063dSJacob Faibussowitsch PetscCall(VecNorm(v, NORM_2, &norm));
2549566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
255c4762a1bSJed Brown if (step > -1) { /* -1 is used to indicate an interpolated value */
25663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "timestep %" PetscInt_FMT " time %g norm %g\n", step, (double)ptime, (double)norm));
257c4762a1bSJed Brown }
2583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
259c4762a1bSJed Brown }
260c4762a1bSJed Brown
261c4762a1bSJed Brown /*
262c4762a1bSJed Brown MySNESMonitor - illustrate how to set user-defined monitoring routine for SNES.
263c4762a1bSJed Brown Input Parameters:
264c4762a1bSJed Brown snes - the SNES context
265c4762a1bSJed Brown its - iteration number
266c4762a1bSJed Brown fnorm - 2-norm function value (may be estimated)
267c4762a1bSJed Brown ctx - optional user-defined context for private data for the
268c4762a1bSJed Brown monitor routine, as set by SNESMonitorSet()
269c4762a1bSJed Brown */
MySNESMonitor(SNES snes,PetscInt its,PetscReal fnorm,PetscViewerAndFormat * vf)270d71ae5a4SJacob Faibussowitsch PetscErrorCode MySNESMonitor(SNES snes, PetscInt its, PetscReal fnorm, PetscViewerAndFormat *vf)
271d71ae5a4SJacob Faibussowitsch {
272c4762a1bSJed Brown PetscFunctionBeginUser;
2739566063dSJacob Faibussowitsch PetscCall(SNESMonitorDefaultShort(snes, its, fnorm, vf));
2743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
275c4762a1bSJed Brown }
276c4762a1bSJed Brown /*TEST
277c4762a1bSJed Brown
278c4762a1bSJed Brown test:
279188af4bfSBarry Smith args: -da_grid_x 20 -ts_max_time 3 -ts_time_step 1e-1 -ts_theta_initial_guess_extrapolate 0 -ts_monitor -ksp_monitor_short
280c4762a1bSJed Brown requires: !single
281c4762a1bSJed Brown
282c4762a1bSJed Brown test:
283c4762a1bSJed Brown suffix: 2
284188af4bfSBarry Smith args: -da_grid_x 20 -ts_max_time 0.11 -ts_time_step 1e-1 -ts_type glle -ts_monitor -ksp_monitor_short
285c4762a1bSJed Brown requires: !single
286c4762a1bSJed Brown
287c4762a1bSJed Brown test:
288c4762a1bSJed Brown suffix: glvis_da_2d_vect
289188af4bfSBarry Smith args: -usemonitor 0 -da_grid_x 20 -ts_max_time 0.3 -ts_time_step 1e-1 -ts_type glle -final_sol glvis: -viewer_glvis_dm_da_bs 2,0
290c4762a1bSJed Brown requires: !single
291c4762a1bSJed Brown
292c4762a1bSJed Brown test:
293c4762a1bSJed Brown suffix: glvis_da_2d_vect_ll
294188af4bfSBarry Smith args: -usemonitor 0 -da_grid_x 20 -ts_max_time 0.3 -ts_time_step 1e-1 -ts_type glle -final_sol glvis: -viewer_glvis_dm_da_bs 2,0 -viewer_glvis_dm_da_ll
295c4762a1bSJed Brown requires: !single
296c4762a1bSJed Brown
297c4762a1bSJed Brown TEST*/
298