1c4762a1bSJed Brown static char help[] = "Bratu nonlinear PDE in 3d.\n\
2c4762a1bSJed Brown We solve the Bratu (SFI - solid fuel ignition) problem in a 3D rectangular\n\
3c4762a1bSJed Brown domain, using distributed arrays (DMDAs) to partition the parallel grid.\n\
4c4762a1bSJed Brown The command line options include:\n\
5c4762a1bSJed Brown -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
6c4762a1bSJed Brown problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)\n\n";
7c4762a1bSJed Brown
8c4762a1bSJed Brown /* ------------------------------------------------------------------------
9c4762a1bSJed Brown
10c4762a1bSJed Brown Solid Fuel Ignition (SFI) problem. This problem is modeled by
11c4762a1bSJed Brown the partial differential equation
12c4762a1bSJed Brown
13c4762a1bSJed Brown -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1,
14c4762a1bSJed Brown
15c4762a1bSJed Brown with boundary conditions
16c4762a1bSJed Brown
17c4762a1bSJed Brown u = 0 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1
18c4762a1bSJed Brown
19c4762a1bSJed Brown A finite difference approximation with the usual 7-point stencil
20c4762a1bSJed Brown is used to discretize the boundary value problem to obtain a nonlinear
21c4762a1bSJed Brown system of equations.
22c4762a1bSJed Brown
23c4762a1bSJed Brown ------------------------------------------------------------------------- */
24c4762a1bSJed Brown
25c4762a1bSJed Brown /*
26c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
27c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this
28c4762a1bSJed Brown file automatically includes:
29c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors
30c4762a1bSJed Brown petscmat.h - matrices
31c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods
32c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners
33c4762a1bSJed Brown petscksp.h - linear solvers
34c4762a1bSJed Brown */
35c4762a1bSJed Brown #include <petscdm.h>
36c4762a1bSJed Brown #include <petscdmda.h>
37c4762a1bSJed Brown #include <petscsnes.h>
38c4762a1bSJed Brown
39c4762a1bSJed Brown /*
40c4762a1bSJed Brown User-defined application context - contains data needed by the
41c4762a1bSJed Brown application-provided call-back routines, FormJacobian() and
42c4762a1bSJed Brown FormFunction().
43c4762a1bSJed Brown */
44c4762a1bSJed Brown typedef struct {
45c4762a1bSJed Brown PetscReal param; /* test problem parameter */
46c4762a1bSJed Brown DM da; /* distributed array data structure */
47c4762a1bSJed Brown } AppCtx;
48c4762a1bSJed Brown
49c4762a1bSJed Brown /*
50c4762a1bSJed Brown User-defined routines
51c4762a1bSJed Brown */
52c4762a1bSJed Brown extern PetscErrorCode FormFunctionLocal(SNES, Vec, Vec, void *);
53c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *);
54c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(AppCtx *, Vec);
55c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *);
56c4762a1bSJed Brown
main(int argc,char ** argv)57d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
58d71ae5a4SJacob Faibussowitsch {
59c4762a1bSJed Brown SNES snes; /* nonlinear solver */
60c4762a1bSJed Brown Vec x, r; /* solution, residual vectors */
61c4762a1bSJed Brown Mat J = NULL; /* Jacobian matrix */
62c4762a1bSJed Brown AppCtx user; /* user-defined work context */
63c4762a1bSJed Brown PetscInt its; /* iterations for convergence */
64c4762a1bSJed Brown MatFDColoring matfdcoloring = NULL;
65c4762a1bSJed Brown PetscBool matrix_free = PETSC_FALSE, coloring = PETSC_FALSE, coloring_ds = PETSC_FALSE, local_coloring = PETSC_FALSE;
66c4762a1bSJed Brown PetscReal bratu_lambda_max = 6.81, bratu_lambda_min = 0., fnorm;
67c4762a1bSJed Brown
68c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69c4762a1bSJed Brown Initialize program
70c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71c4762a1bSJed Brown
72327415f7SBarry Smith PetscFunctionBeginUser;
73c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
74c4762a1bSJed Brown
75c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76c4762a1bSJed Brown Initialize problem parameters
77c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
78c4762a1bSJed Brown user.param = 6.0;
799566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, NULL));
80e00437b9SBarry Smith PetscCheck(user.param < bratu_lambda_max && user.param > bratu_lambda_min, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Lambda is out of range");
81c4762a1bSJed Brown
82c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83c4762a1bSJed Brown Create nonlinear solver context
84c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
859566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
86c4762a1bSJed Brown
87c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors
89c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
909566063dSJacob Faibussowitsch PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 4, 4, 4, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, NULL, &user.da));
919566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.da));
929566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.da));
93c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining
95c4762a1bSJed Brown vectors that are the same types
96c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
979566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.da, &x));
989566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r));
99c4762a1bSJed Brown
100c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101c4762a1bSJed Brown Set function evaluation routine and vector
102c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1039566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)&user));
104c4762a1bSJed Brown
105c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine
107c4762a1bSJed Brown
108c4762a1bSJed Brown Set Jacobian matrix data structure and default Jacobian evaluation
109c4762a1bSJed Brown routine. User can override with:
110c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning
111c4762a1bSJed Brown (unless user explicitly sets preconditioner)
112*7addb90fSBarry Smith -snes_mf_operator : form matrix used to construct the preconditioner as set by the user,
113c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector
114c4762a1bSJed Brown products within Newton-Krylov method
115c4762a1bSJed Brown -fdcoloring : using finite differences with coloring to compute the Jacobian
116c4762a1bSJed Brown
117c4762a1bSJed Brown Note one can use -matfd_coloring wp or ds the only reason for the -fdcoloring_ds option
118c4762a1bSJed Brown below is to test the call to MatFDColoringSetType().
119c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf", &matrix_free, NULL));
1219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-fdcoloring", &coloring, NULL));
1229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-fdcoloring_ds", &coloring_ds, NULL));
1239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-fdcoloring_local", &local_coloring, NULL));
124c4762a1bSJed Brown if (!matrix_free) {
1259566063dSJacob Faibussowitsch PetscCall(DMSetMatType(user.da, MATAIJ));
1269566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(user.da, &J));
127c4762a1bSJed Brown if (coloring) {
128c4762a1bSJed Brown ISColoring iscoloring;
129c4762a1bSJed Brown if (!local_coloring) {
1309566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(user.da, IS_COLORING_GLOBAL, &iscoloring));
1319566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(J, iscoloring, &matfdcoloring));
1322ba42892SBarry Smith PetscCall(MatFDColoringSetFunction(matfdcoloring, (MatFDColoringFn *)FormFunction, &user));
133c4762a1bSJed Brown } else {
1349566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(user.da, IS_COLORING_LOCAL, &iscoloring));
1359566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(J, iscoloring, &matfdcoloring));
1369566063dSJacob Faibussowitsch PetscCall(MatFDColoringUseDM(J, matfdcoloring));
1372ba42892SBarry Smith PetscCall(MatFDColoringSetFunction(matfdcoloring, (MatFDColoringFn *)FormFunctionLocal, &user));
138c4762a1bSJed Brown }
1391baa6e33SBarry Smith if (coloring_ds) PetscCall(MatFDColoringSetType(matfdcoloring, MATMFFD_DS));
1409566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFromOptions(matfdcoloring));
1419566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetUp(J, iscoloring, matfdcoloring));
1429566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, matfdcoloring));
1439566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring));
144c4762a1bSJed Brown } else {
1459566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, &user));
146c4762a1bSJed Brown }
147c4762a1bSJed Brown }
148c4762a1bSJed Brown
149c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150c4762a1bSJed Brown Customize nonlinear solver; set runtime options
151c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1529566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, user.da));
1539566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes));
154c4762a1bSJed Brown
155c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156c4762a1bSJed Brown Evaluate initial guess
157c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess
158c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular,
159c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set
160c4762a1bSJed Brown this vector to zero by calling VecSet().
161c4762a1bSJed Brown */
1629566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(&user, x));
163c4762a1bSJed Brown
164c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165c4762a1bSJed Brown Solve nonlinear system
166c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1679566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x));
1689566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its));
169c4762a1bSJed Brown
170c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171c4762a1bSJed Brown Explicitly check norm of the residual of the solution
172c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1739566063dSJacob Faibussowitsch PetscCall(FormFunction(snes, x, r, (void *)&user));
1749566063dSJacob Faibussowitsch PetscCall(VecNorm(r, NORM_2, &fnorm));
17563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT " fnorm %g\n", its, (double)fnorm));
176c4762a1bSJed Brown
177c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they
179c4762a1bSJed Brown are no longer needed.
180c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
181c4762a1bSJed Brown
1829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J));
1839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r));
1859566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes));
1869566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.da));
1879566063dSJacob Faibussowitsch PetscCall(MatFDColoringDestroy(&matfdcoloring));
1889566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
189b122ec5aSJacob Faibussowitsch return 0;
190c4762a1bSJed Brown }
191c4762a1bSJed Brown /* ------------------------------------------------------------------- */
192c4762a1bSJed Brown /*
193c4762a1bSJed Brown FormInitialGuess - Forms initial approximation.
194c4762a1bSJed Brown
195c4762a1bSJed Brown Input Parameters:
196c4762a1bSJed Brown user - user-defined application context
197c4762a1bSJed Brown X - vector
198c4762a1bSJed Brown
199c4762a1bSJed Brown Output Parameter:
200c4762a1bSJed Brown X - vector
201c4762a1bSJed Brown */
FormInitialGuess(AppCtx * user,Vec X)202d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(AppCtx *user, Vec X)
203d71ae5a4SJacob Faibussowitsch {
204c4762a1bSJed Brown PetscInt i, j, k, Mx, My, Mz, xs, ys, zs, xm, ym, zm;
205c4762a1bSJed Brown PetscReal lambda, temp1, hx, hy, hz, tempk, tempj;
206c4762a1bSJed Brown PetscScalar ***x;
207c4762a1bSJed Brown
208c4762a1bSJed Brown PetscFunctionBeginUser;
2099566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da, PETSC_IGNORE, &Mx, &My, &Mz, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
210c4762a1bSJed Brown
211c4762a1bSJed Brown lambda = user->param;
212c4762a1bSJed Brown hx = 1.0 / (PetscReal)(Mx - 1);
213c4762a1bSJed Brown hy = 1.0 / (PetscReal)(My - 1);
214c4762a1bSJed Brown hz = 1.0 / (PetscReal)(Mz - 1);
215c4762a1bSJed Brown temp1 = lambda / (lambda + 1.0);
216c4762a1bSJed Brown
217c4762a1bSJed Brown /*
218c4762a1bSJed Brown Get a pointer to vector data.
219c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to
220c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent.
221c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to
222c4762a1bSJed Brown the array.
223c4762a1bSJed Brown */
2249566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->da, X, &x));
225c4762a1bSJed Brown
226c4762a1bSJed Brown /*
227c4762a1bSJed Brown Get local grid boundaries (for 3-dimensional DMDA):
228c4762a1bSJed Brown xs, ys, zs - starting grid indices (no ghost points)
229c4762a1bSJed Brown xm, ym, zm - widths of local grid (no ghost points)
230c4762a1bSJed Brown
231c4762a1bSJed Brown */
2329566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->da, &xs, &ys, &zs, &xm, &ym, &zm));
233c4762a1bSJed Brown
234c4762a1bSJed Brown /*
235c4762a1bSJed Brown Compute initial guess over the locally owned part of the grid
236c4762a1bSJed Brown */
237c4762a1bSJed Brown for (k = zs; k < zs + zm; k++) {
238c4762a1bSJed Brown tempk = (PetscReal)(PetscMin(k, Mz - k - 1)) * hz;
239c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
240c4762a1bSJed Brown tempj = PetscMin((PetscReal)(PetscMin(j, My - j - 1)) * hy, tempk);
241c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
242c4762a1bSJed Brown if (i == 0 || j == 0 || k == 0 || i == Mx - 1 || j == My - 1 || k == Mz - 1) {
243c4762a1bSJed Brown /* boundary conditions are all zero Dirichlet */
244c4762a1bSJed Brown x[k][j][i] = 0.0;
245c4762a1bSJed Brown } else {
246c4762a1bSJed Brown x[k][j][i] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, Mx - i - 1)) * hx, tempj));
247c4762a1bSJed Brown }
248c4762a1bSJed Brown }
249c4762a1bSJed Brown }
250c4762a1bSJed Brown }
251c4762a1bSJed Brown
252c4762a1bSJed Brown /*
253c4762a1bSJed Brown Restore vector
254c4762a1bSJed Brown */
2559566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->da, X, &x));
2563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
257c4762a1bSJed Brown }
258c4762a1bSJed Brown /* ------------------------------------------------------------------- */
259c4762a1bSJed Brown /*
260c4762a1bSJed Brown FormFunctionLocal - Evaluates nonlinear function, F(x) on a ghosted patch
261c4762a1bSJed Brown
262c4762a1bSJed Brown Input Parameters:
263c4762a1bSJed Brown . snes - the SNES context
264c4762a1bSJed Brown . localX - input vector, this contains the ghosted region needed
265c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction()
266c4762a1bSJed Brown
267c4762a1bSJed Brown Output Parameter:
268c4762a1bSJed Brown . F - function vector, this does not contain a ghosted region
269c4762a1bSJed Brown */
FormFunctionLocal(SNES snes,Vec localX,Vec F,void * ptr)270d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionLocal(SNES snes, Vec localX, Vec F, void *ptr)
271d71ae5a4SJacob Faibussowitsch {
272c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
273c4762a1bSJed Brown PetscInt i, j, k, Mx, My, Mz, xs, ys, zs, xm, ym, zm;
274c4762a1bSJed Brown PetscReal two = 2.0, lambda, hx, hy, hz, hxhzdhy, hyhzdhx, hxhydhz, sc;
275c4762a1bSJed Brown PetscScalar u_north, u_south, u_east, u_west, u_up, u_down, u, u_xx, u_yy, u_zz, ***x, ***f;
276c4762a1bSJed Brown DM da;
277c4762a1bSJed Brown
278c4762a1bSJed Brown PetscFunctionBeginUser;
2799566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &da));
2809566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, &Mz, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
281c4762a1bSJed Brown
282c4762a1bSJed Brown lambda = user->param;
283c4762a1bSJed Brown hx = 1.0 / (PetscReal)(Mx - 1);
284c4762a1bSJed Brown hy = 1.0 / (PetscReal)(My - 1);
285c4762a1bSJed Brown hz = 1.0 / (PetscReal)(Mz - 1);
286c4762a1bSJed Brown sc = hx * hy * hz * lambda;
287c4762a1bSJed Brown hxhzdhy = hx * hz / hy;
288c4762a1bSJed Brown hyhzdhx = hy * hz / hx;
289c4762a1bSJed Brown hxhydhz = hx * hy / hz;
290c4762a1bSJed Brown
291c4762a1bSJed Brown /*
292c4762a1bSJed Brown Get pointers to vector data
293c4762a1bSJed Brown */
2949566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localX, &x));
2959566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f));
296c4762a1bSJed Brown
297c4762a1bSJed Brown /*
298c4762a1bSJed Brown Get local grid boundaries
299c4762a1bSJed Brown */
3009566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
301c4762a1bSJed Brown
302c4762a1bSJed Brown /*
303c4762a1bSJed Brown Compute function over the locally owned part of the grid
304c4762a1bSJed Brown */
305c4762a1bSJed Brown for (k = zs; k < zs + zm; k++) {
306c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
307c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
308c4762a1bSJed Brown if (i == 0 || j == 0 || k == 0 || i == Mx - 1 || j == My - 1 || k == Mz - 1) {
309c4762a1bSJed Brown f[k][j][i] = x[k][j][i];
310c4762a1bSJed Brown } else {
311c4762a1bSJed Brown u = x[k][j][i];
312c4762a1bSJed Brown u_east = x[k][j][i + 1];
313c4762a1bSJed Brown u_west = x[k][j][i - 1];
314c4762a1bSJed Brown u_north = x[k][j + 1][i];
315c4762a1bSJed Brown u_south = x[k][j - 1][i];
316c4762a1bSJed Brown u_up = x[k + 1][j][i];
317c4762a1bSJed Brown u_down = x[k - 1][j][i];
318c4762a1bSJed Brown u_xx = (-u_east + two * u - u_west) * hyhzdhx;
319c4762a1bSJed Brown u_yy = (-u_north + two * u - u_south) * hxhzdhy;
320c4762a1bSJed Brown u_zz = (-u_up + two * u - u_down) * hxhydhz;
321c4762a1bSJed Brown f[k][j][i] = u_xx + u_yy + u_zz - sc * PetscExpScalar(u);
322c4762a1bSJed Brown }
323c4762a1bSJed Brown }
324c4762a1bSJed Brown }
325c4762a1bSJed Brown }
326c4762a1bSJed Brown
327c4762a1bSJed Brown /*
328c4762a1bSJed Brown Restore vectors
329c4762a1bSJed Brown */
3309566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localX, &x));
3319566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f));
3329566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(11.0 * ym * xm));
3333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
334c4762a1bSJed Brown }
335c4762a1bSJed Brown /* ------------------------------------------------------------------- */
336c4762a1bSJed Brown /*
337c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x) on the entire domain
338c4762a1bSJed Brown
339c4762a1bSJed Brown Input Parameters:
340c4762a1bSJed Brown . snes - the SNES context
341c4762a1bSJed Brown . X - input vector
342c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction()
343c4762a1bSJed Brown
344c4762a1bSJed Brown Output Parameter:
345c4762a1bSJed Brown . F - function vector
346c4762a1bSJed Brown */
FormFunction(SNES snes,Vec X,Vec F,void * ptr)347d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *ptr)
348d71ae5a4SJacob Faibussowitsch {
349c4762a1bSJed Brown Vec localX;
350c4762a1bSJed Brown DM da;
351c4762a1bSJed Brown
352c4762a1bSJed Brown PetscFunctionBeginUser;
3539566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &da));
3549566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX));
355c4762a1bSJed Brown
356c4762a1bSJed Brown /*
357c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process
358c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
359c4762a1bSJed Brown By placing code between these two statements, computations can be
360c4762a1bSJed Brown done while messages are in transition.
361c4762a1bSJed Brown */
3629566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
3639566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
364c4762a1bSJed Brown
3659566063dSJacob Faibussowitsch PetscCall(FormFunctionLocal(snes, localX, F, ptr));
3669566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX));
3673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
368c4762a1bSJed Brown }
369c4762a1bSJed Brown /* ------------------------------------------------------------------- */
370c4762a1bSJed Brown /*
371c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix.
372c4762a1bSJed Brown
373c4762a1bSJed Brown Input Parameters:
374c4762a1bSJed Brown . snes - the SNES context
375c4762a1bSJed Brown . x - input vector
376c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetJacobian()
377c4762a1bSJed Brown
378c4762a1bSJed Brown Output Parameters:
379c4762a1bSJed Brown . A - Jacobian matrix
380*7addb90fSBarry Smith . B - optionally different matrix used to construct the preconditioner
381c4762a1bSJed Brown
382c4762a1bSJed Brown */
FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void * ptr)383d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec X, Mat J, Mat jac, void *ptr)
384d71ae5a4SJacob Faibussowitsch {
385c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; /* user-defined application context */
386c4762a1bSJed Brown Vec localX;
387c4762a1bSJed Brown PetscInt i, j, k, Mx, My, Mz;
388c4762a1bSJed Brown MatStencil col[7], row;
389c4762a1bSJed Brown PetscInt xs, ys, zs, xm, ym, zm;
390c4762a1bSJed Brown PetscScalar lambda, v[7], hx, hy, hz, hxhzdhy, hyhzdhx, hxhydhz, sc, ***x;
391c4762a1bSJed Brown DM da;
392c4762a1bSJed Brown
393c4762a1bSJed Brown PetscFunctionBeginUser;
3949566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &da));
3959566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX));
3969566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, &Mz, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
397c4762a1bSJed Brown
398c4762a1bSJed Brown lambda = user->param;
399c4762a1bSJed Brown hx = 1.0 / (PetscReal)(Mx - 1);
400c4762a1bSJed Brown hy = 1.0 / (PetscReal)(My - 1);
401c4762a1bSJed Brown hz = 1.0 / (PetscReal)(Mz - 1);
402c4762a1bSJed Brown sc = hx * hy * hz * lambda;
403c4762a1bSJed Brown hxhzdhy = hx * hz / hy;
404c4762a1bSJed Brown hyhzdhx = hy * hz / hx;
405c4762a1bSJed Brown hxhydhz = hx * hy / hz;
406c4762a1bSJed Brown
407c4762a1bSJed Brown /*
408c4762a1bSJed Brown Scatter ghost points to local vector, using the 2-step process
409c4762a1bSJed Brown DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
410c4762a1bSJed Brown By placing code between these two statements, computations can be
411c4762a1bSJed Brown done while messages are in transition.
412c4762a1bSJed Brown */
4139566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
4149566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
415c4762a1bSJed Brown
416c4762a1bSJed Brown /*
417c4762a1bSJed Brown Get pointer to vector data
418c4762a1bSJed Brown */
4199566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localX, &x));
420c4762a1bSJed Brown
421c4762a1bSJed Brown /*
422c4762a1bSJed Brown Get local grid boundaries
423c4762a1bSJed Brown */
4249566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
425c4762a1bSJed Brown
426c4762a1bSJed Brown /*
427c4762a1bSJed Brown Compute entries for the locally owned part of the Jacobian.
428c4762a1bSJed Brown - Currently, all PETSc parallel matrix formats are partitioned by
429c4762a1bSJed Brown contiguous chunks of rows across the processors.
430c4762a1bSJed Brown - Each processor needs to insert only elements that it owns
431c4762a1bSJed Brown locally (but any non-local elements will be sent to the
432c4762a1bSJed Brown appropriate processor during matrix assembly).
433c4762a1bSJed Brown - Here, we set all entries for a particular row at once.
434c4762a1bSJed Brown - We can set matrix entries either using either
435c4762a1bSJed Brown MatSetValuesLocal() or MatSetValues(), as discussed above.
436c4762a1bSJed Brown */
437c4762a1bSJed Brown for (k = zs; k < zs + zm; k++) {
438c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
439c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
4409371c9d4SSatish Balay row.k = k;
4419371c9d4SSatish Balay row.j = j;
4429371c9d4SSatish Balay row.i = i;
443c4762a1bSJed Brown /* boundary points */
444c4762a1bSJed Brown if (i == 0 || j == 0 || k == 0 || i == Mx - 1 || j == My - 1 || k == Mz - 1) {
445c4762a1bSJed Brown v[0] = 1.0;
4469566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 1, &row, v, INSERT_VALUES));
447c4762a1bSJed Brown } else {
448c4762a1bSJed Brown /* interior grid points */
4499371c9d4SSatish Balay v[0] = -hxhydhz;
4509371c9d4SSatish Balay col[0].k = k - 1;
4519371c9d4SSatish Balay col[0].j = j;
4529371c9d4SSatish Balay col[0].i = i;
4539371c9d4SSatish Balay v[1] = -hxhzdhy;
4549371c9d4SSatish Balay col[1].k = k;
4559371c9d4SSatish Balay col[1].j = j - 1;
4569371c9d4SSatish Balay col[1].i = i;
4579371c9d4SSatish Balay v[2] = -hyhzdhx;
4589371c9d4SSatish Balay col[2].k = k;
4599371c9d4SSatish Balay col[2].j = j;
4609371c9d4SSatish Balay col[2].i = i - 1;
4619371c9d4SSatish Balay v[3] = 2.0 * (hyhzdhx + hxhzdhy + hxhydhz) - sc * PetscExpScalar(x[k][j][i]);
4629371c9d4SSatish Balay col[3].k = row.k;
4639371c9d4SSatish Balay col[3].j = row.j;
4649371c9d4SSatish Balay col[3].i = row.i;
4659371c9d4SSatish Balay v[4] = -hyhzdhx;
4669371c9d4SSatish Balay col[4].k = k;
4679371c9d4SSatish Balay col[4].j = j;
4689371c9d4SSatish Balay col[4].i = i + 1;
4699371c9d4SSatish Balay v[5] = -hxhzdhy;
4709371c9d4SSatish Balay col[5].k = k;
4719371c9d4SSatish Balay col[5].j = j + 1;
4729371c9d4SSatish Balay col[5].i = i;
4739371c9d4SSatish Balay v[6] = -hxhydhz;
4749371c9d4SSatish Balay col[6].k = k + 1;
4759371c9d4SSatish Balay col[6].j = j;
4769371c9d4SSatish Balay col[6].i = i;
4779566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 7, col, v, INSERT_VALUES));
478c4762a1bSJed Brown }
479c4762a1bSJed Brown }
480c4762a1bSJed Brown }
481c4762a1bSJed Brown }
4829566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localX, &x));
4839566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX));
484c4762a1bSJed Brown
485c4762a1bSJed Brown /*
486c4762a1bSJed Brown Assemble matrix, using the 2-step process:
487c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd().
488c4762a1bSJed Brown */
4899566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
4909566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
491c4762a1bSJed Brown
492c4762a1bSJed Brown /*
493c4762a1bSJed Brown Normally since the matrix has already been assembled above; this
49401c1178eSBarry Smith would do nothing. But in the matrix-free mode -snes_mf_operator
495c4762a1bSJed Brown this tells the "matrix-free" matrix that a new linear system solve
496c4762a1bSJed Brown is about to be done.
497c4762a1bSJed Brown */
498c4762a1bSJed Brown
4999566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
5009566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
501c4762a1bSJed Brown
502c4762a1bSJed Brown /*
503c4762a1bSJed Brown Tell the matrix we will never add a new nonzero location to the
504c4762a1bSJed Brown matrix. If we do, it will generate an error.
505c4762a1bSJed Brown */
5069566063dSJacob Faibussowitsch PetscCall(MatSetOption(jac, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
5073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
508c4762a1bSJed Brown }
509c4762a1bSJed Brown
510c4762a1bSJed Brown /*TEST
511c4762a1bSJed Brown
512c4762a1bSJed Brown test:
513c4762a1bSJed Brown nsize: 4
514c4762a1bSJed Brown args: -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
515c4762a1bSJed Brown
516c4762a1bSJed Brown test:
517c4762a1bSJed Brown suffix: 2
518c4762a1bSJed Brown nsize: 4
519c4762a1bSJed Brown args: -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
520c4762a1bSJed Brown
521c4762a1bSJed Brown test:
522c4762a1bSJed Brown suffix: 3
523c4762a1bSJed Brown nsize: 4
524c4762a1bSJed Brown args: -fdcoloring -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
525c4762a1bSJed Brown
526c4762a1bSJed Brown test:
527c4762a1bSJed Brown suffix: 3_ds
528c4762a1bSJed Brown nsize: 4
529c4762a1bSJed Brown args: -fdcoloring -fdcoloring_ds -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
530c4762a1bSJed Brown
531c4762a1bSJed Brown test:
532c4762a1bSJed Brown suffix: 4
533c4762a1bSJed Brown nsize: 4
534c4762a1bSJed Brown args: -fdcoloring_local -fdcoloring -ksp_monitor_short -da_refine 1
535c4762a1bSJed Brown requires: !single
536c4762a1bSJed Brown
53741ba4c6cSHeeho Park test:
53841ba4c6cSHeeho Park suffix: 5
53941ba4c6cSHeeho Park nsize: 4
54041ba4c6cSHeeho Park args: -fdcoloring_local -fdcoloring -ksp_monitor_short -da_refine 1 -snes_type newtontrdc
54141ba4c6cSHeeho Park requires: !single
54241ba4c6cSHeeho Park
5434a221d59SStefano Zampini test:
5444a221d59SStefano Zampini suffix: 6
5454a221d59SStefano Zampini nsize: 4
5464a221d59SStefano Zampini args: -fdcoloring_local -fdcoloring -da_refine 1 -snes_type newtontr -snes_tr_fallback_type dogleg
5474a221d59SStefano Zampini requires: !single
5484a221d59SStefano Zampini
549c4762a1bSJed Brown TEST*/
550