xref: /petsc/src/snes/tutorials/ex58.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown #include <petscsnes.h>
3c4762a1bSJed Brown #include <petscdm.h>
4c4762a1bSJed Brown #include <petscdmda.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown static const char help[] = "Parallel version of the minimum surface area problem in 2D using DMDA.\n\
7c4762a1bSJed Brown  It solves a system of nonlinear equations in mixed\n\
8c4762a1bSJed Brown complementarity form.This example is based on a\n\
9c4762a1bSJed Brown problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain and\n\
10c4762a1bSJed Brown boundary values along the edges of the domain, the objective is to find the\n\
11c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\
12c4762a1bSJed Brown This application solves this problem using complimentarity -- We are actually\n\
13c4762a1bSJed Brown solving the system  (grad f)_i >= 0, if x_i == l_i \n\
14c4762a1bSJed Brown                     (grad f)_i = 0, if l_i < x_i < u_i \n\
15c4762a1bSJed Brown                     (grad f)_i <= 0, if x_i == u_i  \n\
16c4762a1bSJed Brown where f is the function to be minimized. \n\
17c4762a1bSJed Brown \n\
18c4762a1bSJed Brown The command line options are:\n\
19c4762a1bSJed Brown   -da_grid_x <nx>, where <nx> = number of grid points in the 1st coordinate direction\n\
20c4762a1bSJed Brown   -da_grid_y <ny>, where <ny> = number of grid points in the 2nd coordinate direction\n\
21c4762a1bSJed Brown   -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise\n\
22c4762a1bSJed Brown   -lb <value>, lower bound on the variables\n\
23c4762a1bSJed Brown   -ub <value>, upper bound on the variables\n\n";
24c4762a1bSJed Brown 
25c4762a1bSJed Brown /*
26c4762a1bSJed Brown    User-defined application context - contains data needed by the
27c4762a1bSJed Brown    application-provided call-back routines, FormJacobian() and
28c4762a1bSJed Brown    FormFunction().
29c4762a1bSJed Brown */
30c4762a1bSJed Brown 
31c4762a1bSJed Brown /*
32c4762a1bSJed Brown      This is a new version of the ../tests/ex8.c code
33c4762a1bSJed Brown 
34c4762a1bSJed Brown      Run, for example, with the options ./ex58 -snes_vi_monitor -ksp_monitor -mg_levels_ksp_monitor -pc_type mg -pc_mg_levels 2 -pc_mg_galerkin pmat -ksp_type fgmres
35c4762a1bSJed Brown 
36c4762a1bSJed Brown      Or to run with grid sequencing on the nonlinear problem (note that you do not need to provide the number of
37c4762a1bSJed Brown          multigrid levels, it will be determined automatically based on the number of refinements done)
38c4762a1bSJed Brown 
39c4762a1bSJed Brown       ./ex58 -pc_type mg -ksp_monitor  -snes_view -pc_mg_galerkin pmat -snes_grid_sequence 3
40c4762a1bSJed Brown              -mg_levels_ksp_monitor -snes_vi_monitor -mg_levels_pc_type sor -pc_mg_type full
41c4762a1bSJed Brown 
42c4762a1bSJed Brown */
43c4762a1bSJed Brown 
44c4762a1bSJed Brown typedef struct {
45c4762a1bSJed Brown   PetscScalar *bottom, *top, *left, *right;
46c4762a1bSJed Brown   PetscScalar lb,ub;
47c4762a1bSJed Brown } AppCtx;
48c4762a1bSJed Brown 
49c4762a1bSJed Brown /* -------- User-defined Routines --------- */
50c4762a1bSJed Brown 
51c4762a1bSJed Brown extern PetscErrorCode FormBoundaryConditions(SNES,AppCtx**);
52c4762a1bSJed Brown extern PetscErrorCode DestroyBoundaryConditions(AppCtx**);
53c4762a1bSJed Brown extern PetscErrorCode ComputeInitialGuess(SNES,Vec,void*);
54c4762a1bSJed Brown extern PetscErrorCode FormGradient(SNES,Vec,Vec,void*);
55c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
56c4762a1bSJed Brown extern PetscErrorCode FormBounds(SNES,Vec,Vec);
57c4762a1bSJed Brown 
58c4762a1bSJed Brown int main(int argc, char **argv)
59c4762a1bSJed Brown {
60c4762a1bSJed Brown   PetscErrorCode ierr;              /* used to check for functions returning nonzeros */
61c4762a1bSJed Brown   Vec            x,r;               /* solution and residual vectors */
62c4762a1bSJed Brown   SNES           snes;              /* nonlinear solver context */
63c4762a1bSJed Brown   Mat            J;                 /* Jacobian matrix */
64c4762a1bSJed Brown   DM             da;
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, (char*)0, help);if (ierr) return ierr;
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* Create distributed array to manage the 2d grid */
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da));
70*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
71*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* Extract global vectors from DMDA; */
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&x));
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x, &r));
76c4762a1bSJed Brown 
77*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetMatType(da,MATAIJ));
78*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMatrix(da,&J));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* Create nonlinear solver context */
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetDM(snes,da));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /*  Set function evaluation and Jacobian evaluation  routines */
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFunction(snes,r,FormGradient,NULL));
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian,NULL));
87c4762a1bSJed Brown 
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetComputeApplicationContext(snes,(PetscErrorCode (*)(SNES,void**))FormBoundaryConditions,(PetscErrorCode (*)(void**))DestroyBoundaryConditions));
89c4762a1bSJed Brown 
90*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetComputeInitialGuess(snes,ComputeInitialGuess,NULL));
91c4762a1bSJed Brown 
92*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESVISetComputeVariableBounds(snes,FormBounds));
93c4762a1bSJed Brown 
94*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* Solve the application */
97*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,NULL,x));
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /* Free memory */
100*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
101*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&r));
102*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));
103*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* Free user-created data structures */
106*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   ierr = PetscFinalize();
109c4762a1bSJed Brown   return ierr;
110c4762a1bSJed Brown }
111c4762a1bSJed Brown 
112c4762a1bSJed Brown /* -------------------------------------------------------------------- */
113c4762a1bSJed Brown 
114c4762a1bSJed Brown /*  FormBounds - sets the upper and lower bounds
115c4762a1bSJed Brown 
116c4762a1bSJed Brown     Input Parameters:
117c4762a1bSJed Brown .   snes  - the SNES context
118c4762a1bSJed Brown 
119c4762a1bSJed Brown     Output Parameters:
120c4762a1bSJed Brown .   xl - lower bounds
121c4762a1bSJed Brown .   xu - upper bounds
122c4762a1bSJed Brown */
123c4762a1bSJed Brown PetscErrorCode FormBounds(SNES snes, Vec xl, Vec xu)
124c4762a1bSJed Brown {
125c4762a1bSJed Brown   AppCtx         *ctx;
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   PetscFunctionBeginUser;
128*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetApplicationContext(snes,&ctx));
129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(xl,ctx->lb));
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(xu,ctx->ub));
131c4762a1bSJed Brown   PetscFunctionReturn(0);
132c4762a1bSJed Brown }
133c4762a1bSJed Brown 
134c4762a1bSJed Brown /* -------------------------------------------------------------------- */
135c4762a1bSJed Brown 
136c4762a1bSJed Brown /*  FormGradient - Evaluates gradient of f.
137c4762a1bSJed Brown 
138c4762a1bSJed Brown     Input Parameters:
139c4762a1bSJed Brown .   snes  - the SNES context
140c4762a1bSJed Brown .   X     - input vector
141c4762a1bSJed Brown .   ptr   - optional user-defined context, as set by SNESSetFunction()
142c4762a1bSJed Brown 
143c4762a1bSJed Brown     Output Parameters:
144c4762a1bSJed Brown .   G - vector containing the newly evaluated gradient
145c4762a1bSJed Brown */
146c4762a1bSJed Brown PetscErrorCode FormGradient(SNES snes, Vec X, Vec G, void *ptr)
147c4762a1bSJed Brown {
148c4762a1bSJed Brown   AppCtx      *user;
149c4762a1bSJed Brown   PetscInt    i,j;
150c4762a1bSJed Brown   PetscInt    mx, my;
151c4762a1bSJed Brown   PetscScalar hx,hy, hydhx, hxdhy;
152c4762a1bSJed Brown   PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
153c4762a1bSJed Brown   PetscScalar df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
154c4762a1bSJed Brown   PetscScalar **g, **x;
155c4762a1bSJed Brown   PetscInt    xs,xm,ys,ym;
156c4762a1bSJed Brown   Vec         localX;
157c4762a1bSJed Brown   DM          da;
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   PetscFunctionBeginUser;
160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetApplicationContext(snes,&user));
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
163c4762a1bSJed Brown   hx   = 1.0/(mx+1);hy=1.0/(my+1); hydhx=hy/hx; hxdhy=hx/hy;
164c4762a1bSJed Brown 
165*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(G,0.0));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* Get local vector */
168*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localX));
169c4762a1bSJed Brown   /* Get ghost points */
170*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
171*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
172c4762a1bSJed Brown   /* Get pointer to local vector data */
173*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,localX, &x));
174*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,G, &g));
175c4762a1bSJed Brown 
176*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
177c4762a1bSJed Brown   /* Compute function over the locally owned part of the mesh */
178c4762a1bSJed Brown   for (j=ys; j < ys+ym; j++) {
179c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++) {
180c4762a1bSJed Brown 
181c4762a1bSJed Brown       xc = x[j][i];
182c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
183c4762a1bSJed Brown 
184c4762a1bSJed Brown       if (i==0) { /* left side */
185c4762a1bSJed Brown         xl  = user->left[j+1];
186c4762a1bSJed Brown         xlt = user->left[j+2];
187c4762a1bSJed Brown       } else xl = x[j][i-1];
188c4762a1bSJed Brown 
189c4762a1bSJed Brown       if (j==0) { /* bottom side */
190c4762a1bSJed Brown         xb  = user->bottom[i+1];
191c4762a1bSJed Brown         xrb = user->bottom[i+2];
192c4762a1bSJed Brown       } else xb = x[j-1][i];
193c4762a1bSJed Brown 
194c4762a1bSJed Brown       if (i+1 == mx) { /* right side */
195c4762a1bSJed Brown         xr  = user->right[j+1];
196c4762a1bSJed Brown         xrb = user->right[j];
197c4762a1bSJed Brown       } else xr = x[j][i+1];
198c4762a1bSJed Brown 
199c4762a1bSJed Brown       if (j+1==0+my) { /* top side */
200c4762a1bSJed Brown         xt  = user->top[i+1];
201c4762a1bSJed Brown         xlt = user->top[i];
202c4762a1bSJed Brown       } else xt = x[j+1][i];
203c4762a1bSJed Brown 
204c4762a1bSJed Brown       if (i>0 && j+1<my) xlt = x[j+1][i-1]; /* left top side */
205c4762a1bSJed Brown       if (j>0 && i+1<mx) xrb = x[j-1][i+1]; /* right bottom */
206c4762a1bSJed Brown 
207c4762a1bSJed Brown       d1 = (xc-xl);
208c4762a1bSJed Brown       d2 = (xc-xr);
209c4762a1bSJed Brown       d3 = (xc-xt);
210c4762a1bSJed Brown       d4 = (xc-xb);
211c4762a1bSJed Brown       d5 = (xr-xrb);
212c4762a1bSJed Brown       d6 = (xrb-xb);
213c4762a1bSJed Brown       d7 = (xlt-xl);
214c4762a1bSJed Brown       d8 = (xt-xlt);
215c4762a1bSJed Brown 
216c4762a1bSJed Brown       df1dxc = d1*hydhx;
217c4762a1bSJed Brown       df2dxc = (d1*hydhx + d4*hxdhy);
218c4762a1bSJed Brown       df3dxc = d3*hxdhy;
219c4762a1bSJed Brown       df4dxc = (d2*hydhx + d3*hxdhy);
220c4762a1bSJed Brown       df5dxc = d2*hydhx;
221c4762a1bSJed Brown       df6dxc = d4*hxdhy;
222c4762a1bSJed Brown 
223c4762a1bSJed Brown       d1 /= hx;
224c4762a1bSJed Brown       d2 /= hx;
225c4762a1bSJed Brown       d3 /= hy;
226c4762a1bSJed Brown       d4 /= hy;
227c4762a1bSJed Brown       d5 /= hy;
228c4762a1bSJed Brown       d6 /= hx;
229c4762a1bSJed Brown       d7 /= hy;
230c4762a1bSJed Brown       d8 /= hx;
231c4762a1bSJed Brown 
232c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
233c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
234c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
235c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
236c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
237c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
238c4762a1bSJed Brown 
239c4762a1bSJed Brown       df1dxc /= f1;
240c4762a1bSJed Brown       df2dxc /= f2;
241c4762a1bSJed Brown       df3dxc /= f3;
242c4762a1bSJed Brown       df4dxc /= f4;
243c4762a1bSJed Brown       df5dxc /= f5;
244c4762a1bSJed Brown       df6dxc /= f6;
245c4762a1bSJed Brown 
246c4762a1bSJed Brown       g[j][i] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0;
247c4762a1bSJed Brown 
248c4762a1bSJed Brown     }
249c4762a1bSJed Brown   }
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   /* Restore vectors */
252*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,localX, &x));
253*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,G, &g));
254*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localX));
255*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(67.0*mx*my));
256c4762a1bSJed Brown   PetscFunctionReturn(0);
257c4762a1bSJed Brown }
258c4762a1bSJed Brown 
259c4762a1bSJed Brown /* ------------------------------------------------------------------- */
260c4762a1bSJed Brown /*
261c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
262c4762a1bSJed Brown 
263c4762a1bSJed Brown    Input Parameters:
264c4762a1bSJed Brown .  snes - SNES context
265c4762a1bSJed Brown .  X    - input vector
266c4762a1bSJed Brown .  ptr  - optional user-defined context, as set by SNESSetJacobian()
267c4762a1bSJed Brown 
268c4762a1bSJed Brown    Output Parameters:
269c4762a1bSJed Brown .  tH    - Jacobian matrix
270c4762a1bSJed Brown 
271c4762a1bSJed Brown */
272c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes, Vec X, Mat H, Mat tHPre, void *ptr)
273c4762a1bSJed Brown {
274c4762a1bSJed Brown   AppCtx         *user;
275c4762a1bSJed Brown   PetscInt       i,j,k;
276c4762a1bSJed Brown   PetscInt       mx, my;
277c4762a1bSJed Brown   MatStencil     row,col[7];
278c4762a1bSJed Brown   PetscScalar    hx, hy, hydhx, hxdhy;
279c4762a1bSJed Brown   PetscScalar    f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
280c4762a1bSJed Brown   PetscScalar    hl,hr,ht,hb,hc,htl,hbr;
281c4762a1bSJed Brown   PetscScalar    **x, v[7];
282c4762a1bSJed Brown   PetscBool      assembled;
283c4762a1bSJed Brown   PetscInt       xs,xm,ys,ym;
284c4762a1bSJed Brown   Vec            localX;
285c4762a1bSJed Brown   DM             da;
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   PetscFunctionBeginUser;
288*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
289*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetApplicationContext(snes,&user));
290*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
291c4762a1bSJed Brown   hx   = 1.0/(mx+1); hy=1.0/(my+1); hydhx=hy/hx; hxdhy=hx/hy;
292c4762a1bSJed Brown 
293c4762a1bSJed Brown /* Set various matrix options */
294*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssembled(H,&assembled));
295*5f80ce2aSJacob Faibussowitsch   if (assembled) CHKERRQ(MatZeroEntries(H));
296c4762a1bSJed Brown 
297c4762a1bSJed Brown   /* Get local vector */
298*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localX));
299c4762a1bSJed Brown   /* Get ghost points */
300*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
301*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
302c4762a1bSJed Brown 
303c4762a1bSJed Brown   /* Get pointers to vector data */
304*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,localX, &x));
305c4762a1bSJed Brown 
306*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
307c4762a1bSJed Brown   /* Compute Jacobian over the locally owned part of the mesh */
308c4762a1bSJed Brown   for (j=ys; j< ys+ym; j++) {
309c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++) {
310c4762a1bSJed Brown       xc = x[j][i];
311c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
312c4762a1bSJed Brown 
313c4762a1bSJed Brown       /* Left */
314c4762a1bSJed Brown       if (i==0) {
315c4762a1bSJed Brown         xl  = user->left[j+1];
316c4762a1bSJed Brown         xlt = user->left[j+2];
317c4762a1bSJed Brown       } else xl = x[j][i-1];
318c4762a1bSJed Brown 
319c4762a1bSJed Brown       /* Bottom */
320c4762a1bSJed Brown       if (j==0) {
321c4762a1bSJed Brown         xb  =user->bottom[i+1];
322c4762a1bSJed Brown         xrb = user->bottom[i+2];
323c4762a1bSJed Brown       } else xb = x[j-1][i];
324c4762a1bSJed Brown 
325c4762a1bSJed Brown       /* Right */
326c4762a1bSJed Brown       if (i+1 == mx) {
327c4762a1bSJed Brown         xr  =user->right[j+1];
328c4762a1bSJed Brown         xrb = user->right[j];
329c4762a1bSJed Brown       } else xr = x[j][i+1];
330c4762a1bSJed Brown 
331c4762a1bSJed Brown       /* Top */
332c4762a1bSJed Brown       if (j+1==my) {
333c4762a1bSJed Brown         xt  =user->top[i+1];
334c4762a1bSJed Brown         xlt = user->top[i];
335c4762a1bSJed Brown       } else xt = x[j+1][i];
336c4762a1bSJed Brown 
337c4762a1bSJed Brown       /* Top left */
338c4762a1bSJed Brown       if (i>0 && j+1<my) xlt = x[j+1][i-1];
339c4762a1bSJed Brown 
340c4762a1bSJed Brown       /* Bottom right */
341c4762a1bSJed Brown       if (j>0 && i+1<mx) xrb = x[j-1][i+1];
342c4762a1bSJed Brown 
343c4762a1bSJed Brown       d1 = (xc-xl)/hx;
344c4762a1bSJed Brown       d2 = (xc-xr)/hx;
345c4762a1bSJed Brown       d3 = (xc-xt)/hy;
346c4762a1bSJed Brown       d4 = (xc-xb)/hy;
347c4762a1bSJed Brown       d5 = (xrb-xr)/hy;
348c4762a1bSJed Brown       d6 = (xrb-xb)/hx;
349c4762a1bSJed Brown       d7 = (xlt-xl)/hy;
350c4762a1bSJed Brown       d8 = (xlt-xt)/hx;
351c4762a1bSJed Brown 
352c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
353c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
354c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
355c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
356c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
357c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
358c4762a1bSJed Brown 
359c4762a1bSJed Brown       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
360c4762a1bSJed Brown            (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
361c4762a1bSJed Brown       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
362c4762a1bSJed Brown            (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
363c4762a1bSJed Brown       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
364c4762a1bSJed Brown            (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
365c4762a1bSJed Brown       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
366c4762a1bSJed Brown            (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
367c4762a1bSJed Brown 
368c4762a1bSJed Brown       hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
369c4762a1bSJed Brown       htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
370c4762a1bSJed Brown 
371c4762a1bSJed Brown       hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
372c4762a1bSJed Brown            hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
373c4762a1bSJed Brown            (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2.0*d1*d4)/(f2*f2*f2) +
374c4762a1bSJed Brown            (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2.0*d2*d3)/(f4*f4*f4);
375c4762a1bSJed Brown 
376c4762a1bSJed Brown       hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0;  hc/=2.0;
377c4762a1bSJed Brown 
378c4762a1bSJed Brown       k     =0;
379c4762a1bSJed Brown       row.i = i;row.j= j;
380c4762a1bSJed Brown       /* Bottom */
381c4762a1bSJed Brown       if (j>0) {
382c4762a1bSJed Brown         v[k]     =hb;
383c4762a1bSJed Brown         col[k].i = i; col[k].j=j-1; k++;
384c4762a1bSJed Brown       }
385c4762a1bSJed Brown 
386c4762a1bSJed Brown       /* Bottom right */
387c4762a1bSJed Brown       if (j>0 && i < mx -1) {
388c4762a1bSJed Brown         v[k]     =hbr;
389c4762a1bSJed Brown         col[k].i = i+1; col[k].j = j-1; k++;
390c4762a1bSJed Brown       }
391c4762a1bSJed Brown 
392c4762a1bSJed Brown       /* left */
393c4762a1bSJed Brown       if (i>0) {
394c4762a1bSJed Brown         v[k]     = hl;
395c4762a1bSJed Brown         col[k].i = i-1; col[k].j = j; k++;
396c4762a1bSJed Brown       }
397c4762a1bSJed Brown 
398c4762a1bSJed Brown       /* Centre */
399c4762a1bSJed Brown       v[k]= hc; col[k].i= row.i; col[k].j = row.j; k++;
400c4762a1bSJed Brown 
401c4762a1bSJed Brown       /* Right */
402c4762a1bSJed Brown       if (i < mx-1) {
403c4762a1bSJed Brown         v[k]    = hr;
404c4762a1bSJed Brown         col[k].i= i+1; col[k].j = j;k++;
405c4762a1bSJed Brown       }
406c4762a1bSJed Brown 
407c4762a1bSJed Brown       /* Top left */
408c4762a1bSJed Brown       if (i>0 && j < my-1) {
409c4762a1bSJed Brown         v[k]     = htl;
410c4762a1bSJed Brown         col[k].i = i-1;col[k].j = j+1; k++;
411c4762a1bSJed Brown       }
412c4762a1bSJed Brown 
413c4762a1bSJed Brown       /* Top */
414c4762a1bSJed Brown       if (j < my-1) {
415c4762a1bSJed Brown         v[k]     = ht;
416c4762a1bSJed Brown         col[k].i = i; col[k].j = j+1; k++;
417c4762a1bSJed Brown       }
418c4762a1bSJed Brown 
419*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValuesStencil(H,1,&row,k,col,v,INSERT_VALUES));
420c4762a1bSJed Brown     }
421c4762a1bSJed Brown   }
422c4762a1bSJed Brown 
423c4762a1bSJed Brown   /* Assemble the matrix */
424*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY));
425*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,localX,&x));
426*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY));
427*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localX));
428c4762a1bSJed Brown 
429*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(199.0*mx*my));
430c4762a1bSJed Brown   PetscFunctionReturn(0);
431c4762a1bSJed Brown }
432c4762a1bSJed Brown 
433c4762a1bSJed Brown /* ------------------------------------------------------------------- */
434c4762a1bSJed Brown /*
435c4762a1bSJed Brown    FormBoundaryConditions -  Calculates the boundary conditions for
436c4762a1bSJed Brown    the region.
437c4762a1bSJed Brown 
438c4762a1bSJed Brown    Input Parameter:
439c4762a1bSJed Brown .  user - user-defined application context
440c4762a1bSJed Brown 
441c4762a1bSJed Brown    Output Parameter:
442c4762a1bSJed Brown .  user - user-defined application context
443c4762a1bSJed Brown */
444c4762a1bSJed Brown PetscErrorCode FormBoundaryConditions(SNES snes,AppCtx **ouser)
445c4762a1bSJed Brown {
446c4762a1bSJed Brown   PetscInt       i,j,k,limit=0,maxits=5;
447c4762a1bSJed Brown   PetscInt       mx,my;
448c4762a1bSJed Brown   PetscInt       bsize=0, lsize=0, tsize=0, rsize=0;
449c4762a1bSJed Brown   PetscScalar    one  =1.0, two=2.0, three=3.0;
450c4762a1bSJed Brown   PetscScalar    det,hx,hy,xt=0,yt=0;
451c4762a1bSJed Brown   PetscReal      fnorm, tol=1e-10;
452c4762a1bSJed Brown   PetscScalar    u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
453c4762a1bSJed Brown   PetscScalar    b=-0.5, t=0.5, l=-0.5, r=0.5;
454c4762a1bSJed Brown   PetscScalar    *boundary;
455c4762a1bSJed Brown   AppCtx         *user;
456c4762a1bSJed Brown   DM             da;
457c4762a1bSJed Brown 
458c4762a1bSJed Brown   PetscFunctionBeginUser;
459*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
460*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&user));
461c4762a1bSJed Brown   *ouser   = user;
462c4762a1bSJed Brown   user->lb = .05;
463c4762a1bSJed Brown   user->ub = PETSC_INFINITY;
464*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
465c4762a1bSJed Brown 
466c4762a1bSJed Brown   /* Check if lower and upper bounds are set */
467*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetScalar(NULL,NULL, "-lb", &user->lb, 0));
468*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetScalar(NULL,NULL, "-ub", &user->ub, 0));
469c4762a1bSJed Brown   bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2;
470c4762a1bSJed Brown 
471*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(bsize, &user->bottom));
472*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(tsize, &user->top));
473*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(lsize, &user->left));
474*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(rsize, &user->right));
475c4762a1bSJed Brown 
476c4762a1bSJed Brown   hx= (r-l)/(mx+1.0); hy=(t-b)/(my+1.0);
477c4762a1bSJed Brown 
478c4762a1bSJed Brown   for (j=0; j<4; j++) {
479c4762a1bSJed Brown     if (j==0) {
480c4762a1bSJed Brown       yt       = b;
481c4762a1bSJed Brown       xt       = l;
482c4762a1bSJed Brown       limit    = bsize;
483c4762a1bSJed Brown       boundary = user->bottom;
484c4762a1bSJed Brown     } else if (j==1) {
485c4762a1bSJed Brown       yt       = t;
486c4762a1bSJed Brown       xt       = l;
487c4762a1bSJed Brown       limit    = tsize;
488c4762a1bSJed Brown       boundary = user->top;
489c4762a1bSJed Brown     } else if (j==2) {
490c4762a1bSJed Brown       yt       = b;
491c4762a1bSJed Brown       xt       = l;
492c4762a1bSJed Brown       limit    = lsize;
493c4762a1bSJed Brown       boundary = user->left;
494c4762a1bSJed Brown     } else { /* if  (j==3) */
495c4762a1bSJed Brown       yt       = b;
496c4762a1bSJed Brown       xt       = r;
497c4762a1bSJed Brown       limit    = rsize;
498c4762a1bSJed Brown       boundary = user->right;
499c4762a1bSJed Brown     }
500c4762a1bSJed Brown 
501c4762a1bSJed Brown     for (i=0; i<limit; i++) {
502c4762a1bSJed Brown       u1=xt;
503c4762a1bSJed Brown       u2=-yt;
504c4762a1bSJed Brown       for (k=0; k<maxits; k++) {
505c4762a1bSJed Brown         nf1   = u1 + u1*u2*u2 - u1*u1*u1/three-xt;
506c4762a1bSJed Brown         nf2   = -u2 - u1*u1*u2 + u2*u2*u2/three-yt;
507c4762a1bSJed Brown         fnorm = PetscRealPart(PetscSqrtScalar(nf1*nf1+nf2*nf2));
508c4762a1bSJed Brown         if (fnorm <= tol) break;
509c4762a1bSJed Brown         njac11=one+u2*u2-u1*u1;
510c4762a1bSJed Brown         njac12=two*u1*u2;
511c4762a1bSJed Brown         njac21=-two*u1*u2;
512c4762a1bSJed Brown         njac22=-one - u1*u1 + u2*u2;
513c4762a1bSJed Brown         det   = njac11*njac22-njac21*njac12;
514c4762a1bSJed Brown         u1    = u1-(njac22*nf1-njac12*nf2)/det;
515c4762a1bSJed Brown         u2    = u2-(njac11*nf2-njac21*nf1)/det;
516c4762a1bSJed Brown       }
517c4762a1bSJed Brown 
518c4762a1bSJed Brown       boundary[i]=u1*u1-u2*u2;
519c4762a1bSJed Brown       if (j==0 || j==1) xt=xt+hx;
520c4762a1bSJed Brown       else yt=yt+hy; /* if (j==2 || j==3) */
521c4762a1bSJed Brown     }
522c4762a1bSJed Brown   }
523c4762a1bSJed Brown   PetscFunctionReturn(0);
524c4762a1bSJed Brown }
525c4762a1bSJed Brown 
526c4762a1bSJed Brown PetscErrorCode DestroyBoundaryConditions(AppCtx **ouser)
527c4762a1bSJed Brown {
528c4762a1bSJed Brown   AppCtx         *user = *ouser;
529c4762a1bSJed Brown 
530c4762a1bSJed Brown   PetscFunctionBeginUser;
531*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user->bottom));
532*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user->top));
533*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user->left));
534*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user->right));
535*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(*ouser));
536c4762a1bSJed Brown   PetscFunctionReturn(0);
537c4762a1bSJed Brown }
538c4762a1bSJed Brown 
539c4762a1bSJed Brown /* ------------------------------------------------------------------- */
540c4762a1bSJed Brown /*
541c4762a1bSJed Brown    ComputeInitialGuess - Calculates the initial guess
542c4762a1bSJed Brown 
543c4762a1bSJed Brown    Input Parameters:
544c4762a1bSJed Brown .  user - user-defined application context
545c4762a1bSJed Brown .  X - vector for initial guess
546c4762a1bSJed Brown 
547c4762a1bSJed Brown    Output Parameters:
548c4762a1bSJed Brown .  X - newly computed initial guess
549c4762a1bSJed Brown */
550c4762a1bSJed Brown PetscErrorCode ComputeInitialGuess(SNES snes, Vec X,void *dummy)
551c4762a1bSJed Brown {
552c4762a1bSJed Brown   PetscInt       i,j,mx,my;
553c4762a1bSJed Brown   DM             da;
554c4762a1bSJed Brown   AppCtx         *user;
555c4762a1bSJed Brown   PetscScalar    **x;
556c4762a1bSJed Brown   PetscInt       xs,xm,ys,ym;
557c4762a1bSJed Brown 
558c4762a1bSJed Brown   PetscFunctionBeginUser;
559*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
560*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetApplicationContext(snes,&user));
561c4762a1bSJed Brown 
562*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
563*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
564c4762a1bSJed Brown 
565c4762a1bSJed Brown   /* Get pointers to vector data */
566*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,X,&x));
567c4762a1bSJed Brown   /* Perform local computations */
568c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
569c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++) {
570c4762a1bSJed Brown       x[j][i] = (((j+1.0)*user->bottom[i+1]+(my-j+1.0)*user->top[i+1])/(my+2.0)+((i+1.0)*user->left[j+1]+(mx-i+1.0)*user->right[j+1])/(mx+2.0))/2.0;
571c4762a1bSJed Brown     }
572c4762a1bSJed Brown   }
573c4762a1bSJed Brown   /* Restore vectors */
574*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,X,&x));
575c4762a1bSJed Brown   PetscFunctionReturn(0);
576c4762a1bSJed Brown }
577c4762a1bSJed Brown 
578c4762a1bSJed Brown /*TEST
579c4762a1bSJed Brown 
580c4762a1bSJed Brown    test:
581c4762a1bSJed Brown       args: -snes_type vinewtonrsls -pc_type mg -ksp_monitor_short -pc_mg_galerkin pmat -da_refine 5 -snes_vi_monitor -pc_mg_type full -snes_max_it 100 -snes_converged_reason
582c4762a1bSJed Brown       requires: !single
583c4762a1bSJed Brown 
584c4762a1bSJed Brown    test:
585c4762a1bSJed Brown       suffix: 2
586c4762a1bSJed Brown       args: -snes_type vinewtonssls -pc_type mg -ksp_monitor_short -pc_mg_galerkin pmat -da_refine 5 -snes_vi_monitor -pc_mg_type full -snes_max_it 100 -snes_converged_reason
587c4762a1bSJed Brown       requires: !single
588c4762a1bSJed Brown 
589c4762a1bSJed Brown TEST*/
590