xref: /petsc/src/snes/tutorials/ex14.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Bratu nonlinear PDE in 3d.\n\
3c4762a1bSJed Brown We solve the  Bratu (SFI - solid fuel ignition) problem in a 3D rectangular\n\
4c4762a1bSJed Brown domain, using distributed arrays (DMDAs) to partition the parallel grid.\n\
5c4762a1bSJed Brown The command line options include:\n\
6c4762a1bSJed Brown   -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
7c4762a1bSJed Brown      problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)\n\n";
8c4762a1bSJed Brown 
9c4762a1bSJed Brown /*T
10c4762a1bSJed Brown    Concepts: SNES^parallel Bratu example
11c4762a1bSJed Brown    Concepts: DMDA^using distributed arrays;
12c4762a1bSJed Brown    Processors: n
13c4762a1bSJed Brown T*/
14c4762a1bSJed Brown 
15c4762a1bSJed Brown /* ------------------------------------------------------------------------
16c4762a1bSJed Brown 
17c4762a1bSJed Brown     Solid Fuel Ignition (SFI) problem.  This problem is modeled by
18c4762a1bSJed Brown     the partial differential equation
19c4762a1bSJed Brown 
20c4762a1bSJed Brown             -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
21c4762a1bSJed Brown 
22c4762a1bSJed Brown     with boundary conditions
23c4762a1bSJed Brown 
24c4762a1bSJed Brown              u = 0  for  x = 0, x = 1, y = 0, y = 1, z = 0, z = 1
25c4762a1bSJed Brown 
26c4762a1bSJed Brown     A finite difference approximation with the usual 7-point stencil
27c4762a1bSJed Brown     is used to discretize the boundary value problem to obtain a nonlinear
28c4762a1bSJed Brown     system of equations.
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   ------------------------------------------------------------------------- */
31c4762a1bSJed Brown 
32c4762a1bSJed Brown /*
33c4762a1bSJed Brown    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
34c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
35c4762a1bSJed Brown    file automatically includes:
36c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
37c4762a1bSJed Brown      petscmat.h - matrices
38c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
39c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
40c4762a1bSJed Brown      petscksp.h   - linear solvers
41c4762a1bSJed Brown */
42c4762a1bSJed Brown #include <petscdm.h>
43c4762a1bSJed Brown #include <petscdmda.h>
44c4762a1bSJed Brown #include <petscsnes.h>
45c4762a1bSJed Brown 
46c4762a1bSJed Brown /*
47c4762a1bSJed Brown    User-defined application context - contains data needed by the
48c4762a1bSJed Brown    application-provided call-back routines, FormJacobian() and
49c4762a1bSJed Brown    FormFunction().
50c4762a1bSJed Brown */
51c4762a1bSJed Brown typedef struct {
52c4762a1bSJed Brown   PetscReal param;             /* test problem parameter */
53c4762a1bSJed Brown   DM        da;                /* distributed array data structure */
54c4762a1bSJed Brown } AppCtx;
55c4762a1bSJed Brown 
56c4762a1bSJed Brown /*
57c4762a1bSJed Brown    User-defined routines
58c4762a1bSJed Brown */
59c4762a1bSJed Brown extern PetscErrorCode FormFunctionLocal(SNES,Vec,Vec,void*);
60c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
61c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(AppCtx*,Vec);
62c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
63c4762a1bSJed Brown 
64c4762a1bSJed Brown int main(int argc,char **argv)
65c4762a1bSJed Brown {
66c4762a1bSJed Brown   SNES           snes;                         /* nonlinear solver */
67c4762a1bSJed Brown   Vec            x,r;                          /* solution, residual vectors */
68c4762a1bSJed Brown   Mat            J = NULL;                            /* Jacobian matrix */
69c4762a1bSJed Brown   AppCtx         user;                         /* user-defined work context */
70c4762a1bSJed Brown   PetscInt       its;                          /* iterations for convergence */
71c4762a1bSJed Brown   MatFDColoring  matfdcoloring = NULL;
72c4762a1bSJed Brown   PetscBool      matrix_free = PETSC_FALSE,coloring = PETSC_FALSE, coloring_ds = PETSC_FALSE,local_coloring = PETSC_FALSE;
73c4762a1bSJed Brown   PetscErrorCode ierr;
74c4762a1bSJed Brown   PetscReal      bratu_lambda_max = 6.81,bratu_lambda_min = 0.,fnorm;
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77c4762a1bSJed Brown      Initialize program
78c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83c4762a1bSJed Brown      Initialize problem parameters
84c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85c4762a1bSJed Brown   user.param = 6.0;
86c4762a1bSJed Brown   ierr       = PetscOptionsGetReal(NULL,NULL,"-par",&user.param,NULL);CHKERRQ(ierr);
87*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(user.param >= bratu_lambda_max || user.param <= bratu_lambda_min,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lambda is out of range");
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
90c4762a1bSJed Brown      Create nonlinear solver context
91c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
92c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
96c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97c4762a1bSJed Brown   ierr = 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);CHKERRQ(ierr);
98c4762a1bSJed Brown   ierr = DMSetFromOptions(user.da);CHKERRQ(ierr);
99c4762a1bSJed Brown   ierr = DMSetUp(user.da);CHKERRQ(ierr);
100c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
102c4762a1bSJed Brown      vectors that are the same types
103c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104c4762a1bSJed Brown   ierr = DMCreateGlobalVector(user.da,&x);CHKERRQ(ierr);
105c4762a1bSJed Brown   ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108c4762a1bSJed Brown      Set function evaluation routine and vector
109c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110c4762a1bSJed Brown   ierr = SNESSetFunction(snes,r,FormFunction,(void*)&user);CHKERRQ(ierr);
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
114c4762a1bSJed Brown 
115c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
116c4762a1bSJed Brown      routine. User can override with:
117c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
118c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
119c4762a1bSJed Brown      -snes_mf_operator : form preconditioning matrix as set by the user,
120c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
121c4762a1bSJed Brown                          products within Newton-Krylov method
122c4762a1bSJed Brown      -fdcoloring : using finite differences with coloring to compute the Jacobian
123c4762a1bSJed Brown 
124c4762a1bSJed Brown      Note one can use -matfd_coloring wp or ds the only reason for the -fdcoloring_ds option
125c4762a1bSJed Brown      below is to test the call to MatFDColoringSetType().
126c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-snes_mf",&matrix_free,NULL);CHKERRQ(ierr);
128c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-fdcoloring",&coloring,NULL);CHKERRQ(ierr);
129c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-fdcoloring_ds",&coloring_ds,NULL);CHKERRQ(ierr);
130c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-fdcoloring_local",&local_coloring,NULL);CHKERRQ(ierr);
131c4762a1bSJed Brown   if (!matrix_free) {
132c4762a1bSJed Brown     ierr = DMSetMatType(user.da,MATAIJ);CHKERRQ(ierr);
133c4762a1bSJed Brown     ierr = DMCreateMatrix(user.da,&J);CHKERRQ(ierr);
134c4762a1bSJed Brown     if (coloring) {
135c4762a1bSJed Brown       ISColoring iscoloring;
136c4762a1bSJed Brown       if (!local_coloring) {
137c4762a1bSJed Brown         ierr = DMCreateColoring(user.da,IS_COLORING_GLOBAL,&iscoloring);CHKERRQ(ierr);
138c4762a1bSJed Brown         ierr = MatFDColoringCreate(J,iscoloring,&matfdcoloring);CHKERRQ(ierr);
139c4762a1bSJed Brown         ierr = MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))FormFunction,&user);CHKERRQ(ierr);
140c4762a1bSJed Brown       } else {
141c4762a1bSJed Brown         ierr = DMCreateColoring(user.da,IS_COLORING_LOCAL,&iscoloring);CHKERRQ(ierr);
142c4762a1bSJed Brown         ierr = MatFDColoringCreate(J,iscoloring,&matfdcoloring);CHKERRQ(ierr);
143c4762a1bSJed Brown         ierr = MatFDColoringUseDM(J,matfdcoloring);CHKERRQ(ierr);
144c4762a1bSJed Brown         ierr = MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))FormFunctionLocal,&user);CHKERRQ(ierr);
145c4762a1bSJed Brown       }
146c4762a1bSJed Brown       if (coloring_ds) {
147c4762a1bSJed Brown         ierr = MatFDColoringSetType(matfdcoloring,MATMFFD_DS);CHKERRQ(ierr);
148c4762a1bSJed Brown       }
149c4762a1bSJed Brown       ierr = MatFDColoringSetFromOptions(matfdcoloring);CHKERRQ(ierr);
150c4762a1bSJed Brown       ierr = MatFDColoringSetUp(J,iscoloring,matfdcoloring);CHKERRQ(ierr);
151c4762a1bSJed Brown       ierr = SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);CHKERRQ(ierr);
152c4762a1bSJed Brown       ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
153c4762a1bSJed Brown     } else {
154c4762a1bSJed Brown       ierr = SNESSetJacobian(snes,J,J,FormJacobian,&user);CHKERRQ(ierr);
155c4762a1bSJed Brown     }
156c4762a1bSJed Brown   }
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
160c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
161c4762a1bSJed Brown   ierr = SNESSetDM(snes,user.da);CHKERRQ(ierr);
162c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165c4762a1bSJed Brown      Evaluate initial guess
166c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
167c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
168c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
169c4762a1bSJed Brown      this vector to zero by calling VecSet().
170c4762a1bSJed Brown   */
171c4762a1bSJed Brown   ierr = FormInitialGuess(&user,x);CHKERRQ(ierr);
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
174c4762a1bSJed Brown      Solve nonlinear system
175c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
176c4762a1bSJed Brown   ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr);
177c4762a1bSJed Brown   ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180c4762a1bSJed Brown      Explicitly check norm of the residual of the solution
181c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
182c4762a1bSJed Brown   ierr = FormFunction(snes,x,r,(void*)&user);CHKERRQ(ierr);
183c4762a1bSJed Brown   ierr = VecNorm(r,NORM_2,&fnorm);CHKERRQ(ierr);
184c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D fnorm %g\n",its,(double)fnorm);CHKERRQ(ierr);
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
188c4762a1bSJed Brown      are no longer needed.
189c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
192c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
193c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
194c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
195c4762a1bSJed Brown   ierr = DMDestroy(&user.da);CHKERRQ(ierr);
196c4762a1bSJed Brown   ierr = MatFDColoringDestroy(&matfdcoloring);CHKERRQ(ierr);
197c4762a1bSJed Brown   ierr = PetscFinalize();
198c4762a1bSJed Brown   return ierr;
199c4762a1bSJed Brown }
200c4762a1bSJed Brown /* ------------------------------------------------------------------- */
201c4762a1bSJed Brown /*
202c4762a1bSJed Brown    FormInitialGuess - Forms initial approximation.
203c4762a1bSJed Brown 
204c4762a1bSJed Brown    Input Parameters:
205c4762a1bSJed Brown    user - user-defined application context
206c4762a1bSJed Brown    X - vector
207c4762a1bSJed Brown 
208c4762a1bSJed Brown    Output Parameter:
209c4762a1bSJed Brown    X - vector
210c4762a1bSJed Brown  */
211c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
212c4762a1bSJed Brown {
213c4762a1bSJed Brown   PetscInt       i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm;
214c4762a1bSJed Brown   PetscErrorCode ierr;
215c4762a1bSJed Brown   PetscReal      lambda,temp1,hx,hy,hz,tempk,tempj;
216c4762a1bSJed Brown   PetscScalar    ***x;
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   PetscFunctionBeginUser;
219c4762a1bSJed Brown   ierr = 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);CHKERRQ(ierr);
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   lambda = user->param;
222c4762a1bSJed Brown   hx     = 1.0/(PetscReal)(Mx-1);
223c4762a1bSJed Brown   hy     = 1.0/(PetscReal)(My-1);
224c4762a1bSJed Brown   hz     = 1.0/(PetscReal)(Mz-1);
225c4762a1bSJed Brown   temp1  = lambda/(lambda + 1.0);
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   /*
228c4762a1bSJed Brown      Get a pointer to vector data.
229c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
230c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
231c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
232c4762a1bSJed Brown          the array.
233c4762a1bSJed Brown   */
234c4762a1bSJed Brown   ierr = DMDAVecGetArray(user->da,X,&x);CHKERRQ(ierr);
235c4762a1bSJed Brown 
236c4762a1bSJed Brown   /*
237c4762a1bSJed Brown      Get local grid boundaries (for 3-dimensional DMDA):
238c4762a1bSJed Brown        xs, ys, zs   - starting grid indices (no ghost points)
239c4762a1bSJed Brown        xm, ym, zm   - widths of local grid (no ghost points)
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   */
242c4762a1bSJed Brown   ierr = DMDAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   /*
245c4762a1bSJed Brown      Compute initial guess over the locally owned part of the grid
246c4762a1bSJed Brown   */
247c4762a1bSJed Brown   for (k=zs; k<zs+zm; k++) {
248c4762a1bSJed Brown     tempk = (PetscReal)(PetscMin(k,Mz-k-1))*hz;
249c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
250c4762a1bSJed Brown       tempj = PetscMin((PetscReal)(PetscMin(j,My-j-1))*hy,tempk);
251c4762a1bSJed Brown       for (i=xs; i<xs+xm; i++) {
252c4762a1bSJed Brown         if (i == 0 || j == 0 || k == 0 || i == Mx-1 || j == My-1 || k == Mz-1) {
253c4762a1bSJed Brown           /* boundary conditions are all zero Dirichlet */
254c4762a1bSJed Brown           x[k][j][i] = 0.0;
255c4762a1bSJed Brown         } else {
256c4762a1bSJed Brown           x[k][j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,tempj));
257c4762a1bSJed Brown         }
258c4762a1bSJed Brown       }
259c4762a1bSJed Brown     }
260c4762a1bSJed Brown   }
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   /*
263c4762a1bSJed Brown      Restore vector
264c4762a1bSJed Brown   */
265c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(user->da,X,&x);CHKERRQ(ierr);
266c4762a1bSJed Brown   PetscFunctionReturn(0);
267c4762a1bSJed Brown }
268c4762a1bSJed Brown /* ------------------------------------------------------------------- */
269c4762a1bSJed Brown /*
270c4762a1bSJed Brown    FormFunctionLocal - Evaluates nonlinear function, F(x) on a ghosted patch
271c4762a1bSJed Brown 
272c4762a1bSJed Brown    Input Parameters:
273c4762a1bSJed Brown .  snes - the SNES context
274c4762a1bSJed Brown .  localX - input vector, this contains the ghosted region needed
275c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
276c4762a1bSJed Brown 
277c4762a1bSJed Brown    Output Parameter:
278c4762a1bSJed Brown .  F - function vector, this does not contain a ghosted region
279c4762a1bSJed Brown  */
280c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(SNES snes,Vec localX,Vec F,void *ptr)
281c4762a1bSJed Brown {
282c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
283c4762a1bSJed Brown   PetscErrorCode ierr;
284c4762a1bSJed Brown   PetscInt       i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm;
285c4762a1bSJed Brown   PetscReal      two = 2.0,lambda,hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc;
286c4762a1bSJed Brown   PetscScalar    u_north,u_south,u_east,u_west,u_up,u_down,u,u_xx,u_yy,u_zz,***x,***f;
287c4762a1bSJed Brown   DM             da;
288c4762a1bSJed Brown 
289c4762a1bSJed Brown   PetscFunctionBeginUser;
290c4762a1bSJed Brown   ierr = SNESGetDM(snes,&da);CHKERRQ(ierr);
291c4762a1bSJed Brown   ierr = 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);CHKERRQ(ierr);
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   lambda  = user->param;
294c4762a1bSJed Brown   hx      = 1.0/(PetscReal)(Mx-1);
295c4762a1bSJed Brown   hy      = 1.0/(PetscReal)(My-1);
296c4762a1bSJed Brown   hz      = 1.0/(PetscReal)(Mz-1);
297c4762a1bSJed Brown   sc      = hx*hy*hz*lambda;
298c4762a1bSJed Brown   hxhzdhy = hx*hz/hy;
299c4762a1bSJed Brown   hyhzdhx = hy*hz/hx;
300c4762a1bSJed Brown   hxhydhz = hx*hy/hz;
301c4762a1bSJed Brown 
302c4762a1bSJed Brown   /*
303c4762a1bSJed Brown      Get pointers to vector data
304c4762a1bSJed Brown   */
305c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localX,&x);CHKERRQ(ierr);
306c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
307c4762a1bSJed Brown 
308c4762a1bSJed Brown   /*
309c4762a1bSJed Brown      Get local grid boundaries
310c4762a1bSJed Brown   */
311c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
312c4762a1bSJed Brown 
313c4762a1bSJed Brown   /*
314c4762a1bSJed Brown      Compute function over the locally owned part of the grid
315c4762a1bSJed Brown   */
316c4762a1bSJed Brown   for (k=zs; k<zs+zm; k++) {
317c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
318c4762a1bSJed Brown       for (i=xs; i<xs+xm; i++) {
319c4762a1bSJed Brown         if (i == 0 || j == 0 || k == 0 || i == Mx-1 || j == My-1 || k == Mz-1) {
320c4762a1bSJed Brown           f[k][j][i] = x[k][j][i];
321c4762a1bSJed Brown         } else {
322c4762a1bSJed Brown           u          = x[k][j][i];
323c4762a1bSJed Brown           u_east     = x[k][j][i+1];
324c4762a1bSJed Brown           u_west     = x[k][j][i-1];
325c4762a1bSJed Brown           u_north    = x[k][j+1][i];
326c4762a1bSJed Brown           u_south    = x[k][j-1][i];
327c4762a1bSJed Brown           u_up       = x[k+1][j][i];
328c4762a1bSJed Brown           u_down     = x[k-1][j][i];
329c4762a1bSJed Brown           u_xx       = (-u_east + two*u - u_west)*hyhzdhx;
330c4762a1bSJed Brown           u_yy       = (-u_north + two*u - u_south)*hxhzdhy;
331c4762a1bSJed Brown           u_zz       = (-u_up + two*u - u_down)*hxhydhz;
332c4762a1bSJed Brown           f[k][j][i] = u_xx + u_yy + u_zz - sc*PetscExpScalar(u);
333c4762a1bSJed Brown         }
334c4762a1bSJed Brown       }
335c4762a1bSJed Brown     }
336c4762a1bSJed Brown   }
337c4762a1bSJed Brown 
338c4762a1bSJed Brown   /*
339c4762a1bSJed Brown      Restore vectors
340c4762a1bSJed Brown   */
341c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localX,&x);CHKERRQ(ierr);
342c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
343c4762a1bSJed Brown   ierr = PetscLogFlops(11.0*ym*xm);CHKERRQ(ierr);
344c4762a1bSJed Brown   PetscFunctionReturn(0);
345c4762a1bSJed Brown }
346c4762a1bSJed Brown /* ------------------------------------------------------------------- */
347c4762a1bSJed Brown /*
348c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x) on the entire domain
349c4762a1bSJed Brown 
350c4762a1bSJed Brown    Input Parameters:
351c4762a1bSJed Brown .  snes - the SNES context
352c4762a1bSJed Brown .  X - input vector
353c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
354c4762a1bSJed Brown 
355c4762a1bSJed Brown    Output Parameter:
356c4762a1bSJed Brown .  F - function vector
357c4762a1bSJed Brown  */
358c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr)
359c4762a1bSJed Brown {
360c4762a1bSJed Brown   PetscErrorCode ierr;
361c4762a1bSJed Brown   Vec            localX;
362c4762a1bSJed Brown   DM             da;
363c4762a1bSJed Brown 
364c4762a1bSJed Brown   PetscFunctionBeginUser;
365c4762a1bSJed Brown   ierr = SNESGetDM(snes,&da);CHKERRQ(ierr);
366c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
367c4762a1bSJed Brown 
368c4762a1bSJed Brown   /*
369c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
370c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
371c4762a1bSJed Brown      By placing code between these two statements, computations can be
372c4762a1bSJed Brown      done while messages are in transition.
373c4762a1bSJed Brown   */
374c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
375c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
376c4762a1bSJed Brown 
377c4762a1bSJed Brown   ierr = FormFunctionLocal(snes,localX,F,ptr);CHKERRQ(ierr);
378c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
379c4762a1bSJed Brown   PetscFunctionReturn(0);
380c4762a1bSJed Brown }
381c4762a1bSJed Brown /* ------------------------------------------------------------------- */
382c4762a1bSJed Brown /*
383c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
384c4762a1bSJed Brown 
385c4762a1bSJed Brown    Input Parameters:
386c4762a1bSJed Brown .  snes - the SNES context
387c4762a1bSJed Brown .  x - input vector
388c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetJacobian()
389c4762a1bSJed Brown 
390c4762a1bSJed Brown    Output Parameters:
391c4762a1bSJed Brown .  A - Jacobian matrix
392c4762a1bSJed Brown .  B - optionally different preconditioning matrix
393c4762a1bSJed Brown 
394c4762a1bSJed Brown */
395c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void *ptr)
396c4762a1bSJed Brown {
397c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;  /* user-defined application context */
398c4762a1bSJed Brown   Vec            localX;
399c4762a1bSJed Brown   PetscErrorCode ierr;
400c4762a1bSJed Brown   PetscInt       i,j,k,Mx,My,Mz;
401c4762a1bSJed Brown   MatStencil     col[7],row;
402c4762a1bSJed Brown   PetscInt       xs,ys,zs,xm,ym,zm;
403c4762a1bSJed Brown   PetscScalar    lambda,v[7],hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc,***x;
404c4762a1bSJed Brown   DM             da;
405c4762a1bSJed Brown 
406c4762a1bSJed Brown   PetscFunctionBeginUser;
407c4762a1bSJed Brown   ierr = SNESGetDM(snes,&da);CHKERRQ(ierr);
408c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
409c4762a1bSJed Brown   ierr = 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);CHKERRQ(ierr);
410c4762a1bSJed Brown 
411c4762a1bSJed Brown   lambda  = user->param;
412c4762a1bSJed Brown   hx      = 1.0/(PetscReal)(Mx-1);
413c4762a1bSJed Brown   hy      = 1.0/(PetscReal)(My-1);
414c4762a1bSJed Brown   hz      = 1.0/(PetscReal)(Mz-1);
415c4762a1bSJed Brown   sc      = hx*hy*hz*lambda;
416c4762a1bSJed Brown   hxhzdhy = hx*hz/hy;
417c4762a1bSJed Brown   hyhzdhx = hy*hz/hx;
418c4762a1bSJed Brown   hxhydhz = hx*hy/hz;
419c4762a1bSJed Brown 
420c4762a1bSJed Brown   /*
421c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
422c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
423c4762a1bSJed Brown      By placing code between these two statements, computations can be
424c4762a1bSJed Brown      done while messages are in transition.
425c4762a1bSJed Brown   */
426c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
427c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
428c4762a1bSJed Brown 
429c4762a1bSJed Brown   /*
430c4762a1bSJed Brown      Get pointer to vector data
431c4762a1bSJed Brown   */
432c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localX,&x);CHKERRQ(ierr);
433c4762a1bSJed Brown 
434c4762a1bSJed Brown   /*
435c4762a1bSJed Brown      Get local grid boundaries
436c4762a1bSJed Brown   */
437c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
438c4762a1bSJed Brown 
439c4762a1bSJed Brown   /*
440c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
441c4762a1bSJed Brown       - Currently, all PETSc parallel matrix formats are partitioned by
442c4762a1bSJed Brown         contiguous chunks of rows across the processors.
443c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
444c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
445c4762a1bSJed Brown         appropriate processor during matrix assembly).
446c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
447c4762a1bSJed Brown       - We can set matrix entries either using either
448c4762a1bSJed Brown         MatSetValuesLocal() or MatSetValues(), as discussed above.
449c4762a1bSJed Brown   */
450c4762a1bSJed Brown   for (k=zs; k<zs+zm; k++) {
451c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
452c4762a1bSJed Brown       for (i=xs; i<xs+xm; i++) {
453c4762a1bSJed Brown         row.k = k; row.j = j; row.i = i;
454c4762a1bSJed Brown         /* boundary points */
455c4762a1bSJed Brown         if (i == 0 || j == 0 || k == 0|| i == Mx-1 || j == My-1 || k == Mz-1) {
456c4762a1bSJed Brown           v[0] = 1.0;
457c4762a1bSJed Brown           ierr = MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr);
458c4762a1bSJed Brown         } else {
459c4762a1bSJed Brown           /* interior grid points */
460c4762a1bSJed Brown           v[0] = -hxhydhz; col[0].k=k-1;col[0].j=j;  col[0].i = i;
461c4762a1bSJed Brown           v[1] = -hxhzdhy; col[1].k=k;  col[1].j=j-1;col[1].i = i;
462c4762a1bSJed Brown           v[2] = -hyhzdhx; col[2].k=k;  col[2].j=j;  col[2].i = i-1;
463c4762a1bSJed Brown           v[3] = 2.0*(hyhzdhx+hxhzdhy+hxhydhz)-sc*PetscExpScalar(x[k][j][i]);col[3].k=row.k;col[3].j=row.j;col[3].i = row.i;
464c4762a1bSJed Brown           v[4] = -hyhzdhx; col[4].k=k;  col[4].j=j;  col[4].i = i+1;
465c4762a1bSJed Brown           v[5] = -hxhzdhy; col[5].k=k;  col[5].j=j+1;col[5].i = i;
466c4762a1bSJed Brown           v[6] = -hxhydhz; col[6].k=k+1;col[6].j=j;  col[6].i = i;
467c4762a1bSJed Brown           ierr = MatSetValuesStencil(jac,1,&row,7,col,v,INSERT_VALUES);CHKERRQ(ierr);
468c4762a1bSJed Brown         }
469c4762a1bSJed Brown       }
470c4762a1bSJed Brown     }
471c4762a1bSJed Brown   }
472c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localX,&x);CHKERRQ(ierr);
473c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
474c4762a1bSJed Brown 
475c4762a1bSJed Brown   /*
476c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
477c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd().
478c4762a1bSJed Brown   */
479c4762a1bSJed Brown   ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
480c4762a1bSJed Brown   ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
481c4762a1bSJed Brown 
482c4762a1bSJed Brown   /*
483c4762a1bSJed Brown      Normally since the matrix has already been assembled above; this
484c4762a1bSJed Brown      would do nothing. But in the matrix free mode -snes_mf_operator
485c4762a1bSJed Brown      this tells the "matrix-free" matrix that a new linear system solve
486c4762a1bSJed Brown      is about to be done.
487c4762a1bSJed Brown   */
488c4762a1bSJed Brown 
489c4762a1bSJed Brown   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
490c4762a1bSJed Brown   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
491c4762a1bSJed Brown 
492c4762a1bSJed Brown   /*
493c4762a1bSJed Brown      Tell the matrix we will never add a new nonzero location to the
494c4762a1bSJed Brown      matrix. If we do, it will generate an error.
495c4762a1bSJed Brown   */
496c4762a1bSJed Brown   ierr = MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
497c4762a1bSJed Brown   PetscFunctionReturn(0);
498c4762a1bSJed Brown }
499c4762a1bSJed Brown 
500c4762a1bSJed Brown /*TEST
501c4762a1bSJed Brown 
502c4762a1bSJed Brown    test:
503c4762a1bSJed Brown       nsize: 4
504c4762a1bSJed Brown       args: -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
505c4762a1bSJed Brown 
506c4762a1bSJed Brown    test:
507c4762a1bSJed Brown       suffix: 2
508c4762a1bSJed Brown       nsize: 4
509c4762a1bSJed Brown       args: -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
510c4762a1bSJed Brown 
511c4762a1bSJed Brown    test:
512c4762a1bSJed Brown       suffix: 3
513c4762a1bSJed Brown       nsize: 4
514c4762a1bSJed Brown       args: -fdcoloring -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
515c4762a1bSJed Brown 
516c4762a1bSJed Brown    test:
517c4762a1bSJed Brown       suffix: 3_ds
518c4762a1bSJed Brown       nsize: 4
519c4762a1bSJed Brown       args: -fdcoloring -fdcoloring_ds -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
520c4762a1bSJed Brown 
521c4762a1bSJed Brown    test:
522c4762a1bSJed Brown       suffix: 4
523c4762a1bSJed Brown       nsize: 4
524c4762a1bSJed Brown       args: -fdcoloring_local -fdcoloring -ksp_monitor_short -da_refine 1
525c4762a1bSJed Brown       requires: !single
526c4762a1bSJed Brown 
52741ba4c6cSHeeho Park    test:
52841ba4c6cSHeeho Park       suffix: 5
52941ba4c6cSHeeho Park       nsize: 4
53041ba4c6cSHeeho Park       args: -fdcoloring_local -fdcoloring -ksp_monitor_short -da_refine 1 -snes_type newtontrdc
53141ba4c6cSHeeho Park       requires: !single
53241ba4c6cSHeeho Park 
533c4762a1bSJed Brown TEST*/
534