xref: /petsc/src/snes/tutorials/ex46.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Surface processes in geophysics.\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown /*T
4*c4762a1bSJed Brown    Concepts: SNES^parallel Surface process example
5*c4762a1bSJed Brown    Concepts: DMDA^using distributed arrays;
6*c4762a1bSJed Brown    Concepts: IS coloirng types;
7*c4762a1bSJed Brown    Processors: n
8*c4762a1bSJed Brown T*/
9*c4762a1bSJed Brown 
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown #include <petscsnes.h>
14*c4762a1bSJed Brown #include <petscdm.h>
15*c4762a1bSJed Brown #include <petscdmda.h>
16*c4762a1bSJed Brown 
17*c4762a1bSJed Brown /*
18*c4762a1bSJed Brown    User-defined application context - contains data needed by the
19*c4762a1bSJed Brown    application-provided call-back routines, FormJacobianLocal() and
20*c4762a1bSJed Brown    FormFunctionLocal().
21*c4762a1bSJed Brown */
22*c4762a1bSJed Brown typedef struct {
23*c4762a1bSJed Brown   PetscReal   D;  /* The diffusion coefficient */
24*c4762a1bSJed Brown   PetscReal   K;  /* The advection coefficient */
25*c4762a1bSJed Brown   PetscInt    m;  /* Exponent for A */
26*c4762a1bSJed Brown } AppCtx;
27*c4762a1bSJed Brown 
28*c4762a1bSJed Brown /*
29*c4762a1bSJed Brown    User-defined routines
30*c4762a1bSJed Brown */
31*c4762a1bSJed Brown extern PetscErrorCode FormFunctionLocal(DMDALocalInfo*,PetscScalar**,PetscScalar**,AppCtx*);
32*c4762a1bSJed Brown extern PetscErrorCode FormJacobianLocal(DMDALocalInfo*,PetscScalar**,Mat,AppCtx*);
33*c4762a1bSJed Brown 
34*c4762a1bSJed Brown int main(int argc,char **argv)
35*c4762a1bSJed Brown {
36*c4762a1bSJed Brown   SNES           snes;                         /* nonlinear solver */
37*c4762a1bSJed Brown   AppCtx         user;                         /* user-defined work context */
38*c4762a1bSJed Brown   PetscInt       its;                          /* iterations for convergence */
39*c4762a1bSJed Brown   PetscErrorCode ierr;
40*c4762a1bSJed Brown   DM             da;
41*c4762a1bSJed Brown 
42*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43*c4762a1bSJed Brown      Initialize program
44*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
45*c4762a1bSJed Brown 
46*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49*c4762a1bSJed Brown      Initialize problem parameters
50*c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51*c4762a1bSJed Brown   ierr   = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Surface Process Problem Options", "SNES");CHKERRQ(ierr);
52*c4762a1bSJed Brown   user.D = 1.0;
53*c4762a1bSJed Brown   ierr   = PetscOptionsReal("-D", "The diffusion coefficient D", __FILE__, user.D, &user.D, NULL);CHKERRQ(ierr);
54*c4762a1bSJed Brown   user.K = 1.0;
55*c4762a1bSJed Brown   ierr   = PetscOptionsReal("-K", "The advection coefficient K", __FILE__, user.K, &user.K, NULL);CHKERRQ(ierr);
56*c4762a1bSJed Brown   user.m = 1;
57*c4762a1bSJed Brown   ierr   = PetscOptionsInt("-m", "The exponent for A", __FILE__, user.m, &user.m, NULL);CHKERRQ(ierr);
58*c4762a1bSJed Brown   ierr   = PetscOptionsEnd();CHKERRQ(ierr);
59*c4762a1bSJed Brown 
60*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61*c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
62*c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63*c4762a1bSJed Brown   ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr);
64*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
65*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
66*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);CHKERRQ(ierr);
67*c4762a1bSJed Brown   ierr = DMSetApplicationContext(da,&user);CHKERRQ(ierr);
68*c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
69*c4762a1bSJed Brown   ierr = SNESSetDM(snes, da);CHKERRQ(ierr);
70*c4762a1bSJed Brown 
71*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72*c4762a1bSJed Brown      Set local function evaluation routine
73*c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74*c4762a1bSJed Brown   ierr = DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user);CHKERRQ(ierr);
75*c4762a1bSJed Brown 
76*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77*c4762a1bSJed Brown      Customize solver; set runtime options
78*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
79*c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
80*c4762a1bSJed Brown 
81*c4762a1bSJed Brown 
82*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83*c4762a1bSJed Brown      Solve nonlinear system
84*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85*c4762a1bSJed Brown   ierr = SNESSolve(snes,0,0);CHKERRQ(ierr);
86*c4762a1bSJed Brown   ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
87*c4762a1bSJed Brown 
88*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its);CHKERRQ(ierr);
91*c4762a1bSJed Brown 
92*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
94*c4762a1bSJed Brown      are no longer needed.
95*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96*c4762a1bSJed Brown 
97*c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
99*c4762a1bSJed Brown 
100*c4762a1bSJed Brown   ierr = PetscFinalize();
101*c4762a1bSJed Brown   return ierr;
102*c4762a1bSJed Brown }
103*c4762a1bSJed Brown 
104*c4762a1bSJed Brown PetscScalar funcU(DMDACoor2d *coords)
105*c4762a1bSJed Brown {
106*c4762a1bSJed Brown   return coords->x + coords->y;
107*c4762a1bSJed Brown }
108*c4762a1bSJed Brown 
109*c4762a1bSJed Brown PetscScalar funcA(PetscScalar z, AppCtx *user)
110*c4762a1bSJed Brown {
111*c4762a1bSJed Brown   PetscScalar v = 1.0;
112*c4762a1bSJed Brown   PetscInt    i;
113*c4762a1bSJed Brown 
114*c4762a1bSJed Brown   for (i = 0; i < user->m; ++i) v *= z;
115*c4762a1bSJed Brown   return v;
116*c4762a1bSJed Brown }
117*c4762a1bSJed Brown 
118*c4762a1bSJed Brown PetscScalar funcADer(PetscScalar z, AppCtx *user)
119*c4762a1bSJed Brown {
120*c4762a1bSJed Brown   PetscScalar v = 1.0;
121*c4762a1bSJed Brown   PetscInt    i;
122*c4762a1bSJed Brown 
123*c4762a1bSJed Brown   for (i = 0; i < user->m-1; ++i) v *= z;
124*c4762a1bSJed Brown   return (PetscScalar)user->m*v;
125*c4762a1bSJed Brown }
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown /*
128*c4762a1bSJed Brown    FormFunctionLocal - Evaluates nonlinear function, F(x).
129*c4762a1bSJed Brown */
130*c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
131*c4762a1bSJed Brown {
132*c4762a1bSJed Brown   DM             coordDA;
133*c4762a1bSJed Brown   Vec            coordinates;
134*c4762a1bSJed Brown   DMDACoor2d     **coords;
135*c4762a1bSJed Brown   PetscScalar    u, ux, uy, uxx, uyy;
136*c4762a1bSJed Brown   PetscReal      D, K, hx, hy, hxdhy, hydhx;
137*c4762a1bSJed Brown   PetscInt       i,j;
138*c4762a1bSJed Brown   PetscErrorCode ierr;
139*c4762a1bSJed Brown 
140*c4762a1bSJed Brown   PetscFunctionBeginUser;
141*c4762a1bSJed Brown   D     = user->D;
142*c4762a1bSJed Brown   K     = user->K;
143*c4762a1bSJed Brown   hx    = 1.0/(PetscReal)(info->mx-1);
144*c4762a1bSJed Brown   hy    = 1.0/(PetscReal)(info->my-1);
145*c4762a1bSJed Brown   hxdhy = hx/hy;
146*c4762a1bSJed Brown   hydhx = hy/hx;
147*c4762a1bSJed Brown   /*
148*c4762a1bSJed Brown      Compute function over the locally owned part of the grid
149*c4762a1bSJed Brown   */
150*c4762a1bSJed Brown   ierr = DMGetCoordinateDM(info->da, &coordDA);CHKERRQ(ierr);
151*c4762a1bSJed Brown   ierr = DMGetCoordinates(info->da, &coordinates);CHKERRQ(ierr);
152*c4762a1bSJed Brown   ierr = DMDAVecGetArray(coordDA, coordinates, &coords);CHKERRQ(ierr);
153*c4762a1bSJed Brown   for (j=info->ys; j<info->ys+info->ym; j++) {
154*c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
155*c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) f[j][i] = x[j][i];
156*c4762a1bSJed Brown       else {
157*c4762a1bSJed Brown         u       = x[j][i];
158*c4762a1bSJed Brown         ux      = (x[j][i+1] - x[j][i])/hx;
159*c4762a1bSJed Brown         uy      = (x[j+1][i] - x[j][i])/hy;
160*c4762a1bSJed Brown         uxx     = (2.0*u - x[j][i-1] - x[j][i+1])*hydhx;
161*c4762a1bSJed Brown         uyy     = (2.0*u - x[j-1][i] - x[j+1][i])*hxdhy;
162*c4762a1bSJed Brown         f[j][i] = D*(uxx + uyy) - (K*funcA(x[j][i], user)*PetscSqrtScalar(ux*ux + uy*uy) + funcU(&coords[j][i]))*hx*hy;
163*c4762a1bSJed Brown         if (PetscIsInfOrNanScalar(f[j][i])) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FP, "Invalid residual: %g", (double)PetscRealPart(f[j][i]));
164*c4762a1bSJed Brown       }
165*c4762a1bSJed Brown     }
166*c4762a1bSJed Brown   }
167*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(coordDA, coordinates, &coords);CHKERRQ(ierr);
168*c4762a1bSJed Brown   ierr = PetscLogFlops(11*info->ym*info->xm);CHKERRQ(ierr);
169*c4762a1bSJed Brown   PetscFunctionReturn(0);
170*c4762a1bSJed Brown }
171*c4762a1bSJed Brown 
172*c4762a1bSJed Brown /*
173*c4762a1bSJed Brown    FormJacobianLocal - Evaluates Jacobian matrix.
174*c4762a1bSJed Brown */
175*c4762a1bSJed Brown PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,PetscScalar **x,Mat jac,AppCtx *user)
176*c4762a1bSJed Brown {
177*c4762a1bSJed Brown   MatStencil     col[5], row;
178*c4762a1bSJed Brown   PetscScalar    D, K, A, v[5], hx, hy, hxdhy, hydhx, ux, uy;
179*c4762a1bSJed Brown   PetscReal      normGradZ;
180*c4762a1bSJed Brown   PetscInt       i, j,k;
181*c4762a1bSJed Brown   PetscErrorCode ierr;
182*c4762a1bSJed Brown 
183*c4762a1bSJed Brown   PetscFunctionBeginUser;
184*c4762a1bSJed Brown   D     = user->D;
185*c4762a1bSJed Brown   K     = user->K;
186*c4762a1bSJed Brown   hx    = 1.0/(PetscReal)(info->mx-1);
187*c4762a1bSJed Brown   hy    = 1.0/(PetscReal)(info->my-1);
188*c4762a1bSJed Brown   hxdhy = hx/hy;
189*c4762a1bSJed Brown   hydhx = hy/hx;
190*c4762a1bSJed Brown 
191*c4762a1bSJed Brown   /*
192*c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
193*c4762a1bSJed Brown       - Currently, all PETSc parallel matrix formats are partitioned by
194*c4762a1bSJed Brown         contiguous chunks of rows across the processors.
195*c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
196*c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
197*c4762a1bSJed Brown         appropriate processor during matrix assembly).
198*c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
199*c4762a1bSJed Brown       - We can set matrix entries either using either
200*c4762a1bSJed Brown         MatSetValuesLocal() or MatSetValues(), as discussed above.
201*c4762a1bSJed Brown   */
202*c4762a1bSJed Brown   for (j=info->ys; j<info->ys+info->ym; j++) {
203*c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
204*c4762a1bSJed Brown       row.j = j; row.i = i;
205*c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
206*c4762a1bSJed Brown         /* boundary points */
207*c4762a1bSJed Brown         v[0] = 1.0;
208*c4762a1bSJed Brown         ierr = MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr);
209*c4762a1bSJed Brown       } else {
210*c4762a1bSJed Brown         /* interior grid points */
211*c4762a1bSJed Brown         ux        = (x[j][i+1] - x[j][i])/hx;
212*c4762a1bSJed Brown         uy        = (x[j+1][i] - x[j][i])/hy;
213*c4762a1bSJed Brown         normGradZ = PetscRealPart(PetscSqrtScalar(ux*ux + uy*uy));
214*c4762a1bSJed Brown         if (normGradZ < 1.0e-8) normGradZ = 1.0e-8;
215*c4762a1bSJed Brown         A = funcA(x[j][i], user);
216*c4762a1bSJed Brown 
217*c4762a1bSJed Brown         v[0] = -D*hxdhy;                                                                          col[0].j = j - 1; col[0].i = i;
218*c4762a1bSJed Brown         v[1] = -D*hydhx;                                                                          col[1].j = j;     col[1].i = i-1;
219*c4762a1bSJed Brown         v[2] = D*2.0*(hydhx + hxdhy) + K*(funcADer(x[j][i], user)*normGradZ - A/normGradZ)*hx*hy; col[2].j = row.j; col[2].i = row.i;
220*c4762a1bSJed Brown         v[3] = -D*hydhx + K*A*hx*hy/(2.0*normGradZ);                                              col[3].j = j;     col[3].i = i+1;
221*c4762a1bSJed Brown         v[4] = -D*hxdhy + K*A*hx*hy/(2.0*normGradZ);                                              col[4].j = j + 1; col[4].i = i;
222*c4762a1bSJed Brown         for (k = 0; k < 5; ++k) {
223*c4762a1bSJed Brown           if (PetscIsInfOrNanScalar(v[k])) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FP, "Invalid residual: %g", (double)PetscRealPart(v[k]));
224*c4762a1bSJed Brown         }
225*c4762a1bSJed Brown         ierr = MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);CHKERRQ(ierr);
226*c4762a1bSJed Brown       }
227*c4762a1bSJed Brown     }
228*c4762a1bSJed Brown   }
229*c4762a1bSJed Brown 
230*c4762a1bSJed Brown   /*
231*c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
232*c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd().
233*c4762a1bSJed Brown   */
234*c4762a1bSJed Brown   ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
235*c4762a1bSJed Brown   ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
236*c4762a1bSJed Brown   /*
237*c4762a1bSJed Brown      Tell the matrix we will never add a new nonzero location to the
238*c4762a1bSJed Brown      matrix. If we do, it will generate an error.
239*c4762a1bSJed Brown   */
240*c4762a1bSJed Brown   ierr = MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
241*c4762a1bSJed Brown   PetscFunctionReturn(0);
242*c4762a1bSJed Brown }
243*c4762a1bSJed Brown 
244*c4762a1bSJed Brown 
245*c4762a1bSJed Brown /*TEST
246*c4762a1bSJed Brown 
247*c4762a1bSJed Brown    test:
248*c4762a1bSJed Brown       args: -snes_view -snes_monitor_short -da_refine 1 -pc_type mg -ksp_type fgmres -pc_mg_type full -mg_levels_ksp_chebyshev_esteig 0.5,1.1
249*c4762a1bSJed Brown 
250*c4762a1bSJed Brown    test:
251*c4762a1bSJed Brown       suffix: ew_1
252*c4762a1bSJed Brown       args: -snes_monitor_short -ksp_converged_reason -da_grid_x 20 -da_grid_y 20 -snes_ksp_ew -snes_ksp_ew_version 1
253*c4762a1bSJed Brown       requires: !single
254*c4762a1bSJed Brown 
255*c4762a1bSJed Brown    test:
256*c4762a1bSJed Brown       suffix: ew_2
257*c4762a1bSJed Brown       args: -snes_monitor_short -ksp_converged_reason -da_grid_x 20 -da_grid_y 20 -snes_ksp_ew -snes_ksp_ew_version 2
258*c4762a1bSJed Brown 
259*c4762a1bSJed Brown    test:
260*c4762a1bSJed Brown       suffix: ew_3
261*c4762a1bSJed Brown       args: -snes_monitor_short -ksp_converged_reason -da_grid_x 20 -da_grid_y 20 -snes_ksp_ew -snes_ksp_ew_version 3
262*c4762a1bSJed Brown       requires: !single
263*c4762a1bSJed Brown 
264*c4762a1bSJed Brown    test:
265*c4762a1bSJed Brown       suffix: fm_rise_2
266*c4762a1bSJed Brown       args: -K 3 -m 1 -D 0.2 -snes_monitor_short -snes_type ngmres -snes_npc_side right -npc_snes_type newtonls -npc_snes_linesearch_type basic -snes_ngmres_restart_it 1 -snes_ngmres_restart_fm_rise
267*c4762a1bSJed Brown       requires: !single
268*c4762a1bSJed Brown 
269*c4762a1bSJed Brown    test:
270*c4762a1bSJed Brown       suffix: fm_rise_4
271*c4762a1bSJed Brown       args: -K 3 -m 1 -D 0.2 -snes_monitor_short -snes_type ngmres -snes_npc_side right -npc_snes_type newtonls -npc_snes_linesearch_type basic -snes_ngmres_restart_it 2 -snes_ngmres_restart_fm_rise -snes_rtol 1.e-2 -snes_max_it 5
272*c4762a1bSJed Brown 
273*c4762a1bSJed Brown TEST*/
274