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