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