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