Actual source code: ex46.c

petsc-3.8.4 2018-03-24
Report Typos and Errors
  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: T*/


 11:  #include <petscsnes.h>
 12:  #include <petscdm.h>
 13:  #include <petscdmda.h>

 15: /*
 16:    User-defined application context - contains data needed by the
 17:    application-provided call-back routines, FormJacobianLocal() and
 18:    FormFunctionLocal().
 19: */
 20: typedef struct {
 21:   PetscReal   D;  /* The diffusion coefficient */
 22:   PetscReal   K;  /* The advection coefficient */
 23:   PetscInt    m;  /* Exponent for A */
 24: } AppCtx;

 26: /*
 27:    User-defined routines
 28: */
 29: extern PetscErrorCode FormFunctionLocal(DMDALocalInfo*,PetscScalar**,PetscScalar**,AppCtx*);
 30: extern PetscErrorCode FormJacobianLocal(DMDALocalInfo*,PetscScalar**,Mat,AppCtx*);

 32: int main(int argc,char **argv)
 33: {
 34:   SNES           snes;                         /* nonlinear solver */
 35:   AppCtx         user;                         /* user-defined work context */
 36:   PetscInt       its;                          /* iterations for convergence */
 38:   DM             da;

 40:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 41:      Initialize program
 42:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 44:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 46:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 47:      Initialize problem parameters
 48:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 49:   PetscOptionsBegin(PETSC_COMM_WORLD, "", "Surface Process Problem Options", "SNES");
 50:   user.D = 1.0;
 51:   PetscOptionsReal("-D", "The diffusion coefficient D", __FILE__, user.D, &user.D, NULL);
 52:   user.K = 1.0;
 53:   PetscOptionsReal("-K", "The advection coefficient K", __FILE__, user.K, &user.K, NULL);
 54:   user.m = 1;
 55:   PetscOptionsInt("-m", "The exponent for A", __FILE__, user.m, &user.m, NULL);
 56:   PetscOptionsEnd();

 58:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 59:      Create distributed array (DMDA) to manage parallel grid and vectors
 60:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 61:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);
 62:   DMSetFromOptions(da);
 63:   DMSetUp(da);
 64:   DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);
 65:   DMSetApplicationContext(da,&user);
 66:   SNESCreate(PETSC_COMM_WORLD, &snes);
 67:   SNESSetDM(snes, da);

 69:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 70:      Set local function evaluation routine
 71:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 72:   DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user);

 74:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 75:      Customize solver; set runtime options
 76:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 77:   SNESSetFromOptions(snes);


 80:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 81:      Solve nonlinear system
 82:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 83:   SNESSolve(snes,0,0);
 84:   SNESGetIterationNumber(snes,&its);

 86:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 87:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 88:   PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its);

 90:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 91:      Free work space.  All PETSc objects should be destroyed when they
 92:      are no longer needed.
 93:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 95:   SNESDestroy(&snes);
 96:   DMDestroy(&da);

 98:   PetscFinalize();
 99:   return ierr;
100: }

102: PetscScalar funcU(DMDACoor2d *coords)
103: {
104:   return coords->x + coords->y;
105: }

107: PetscScalar funcA(PetscScalar z, AppCtx *user)
108: {
109:   PetscScalar v = 1.0;
110:   PetscInt    i;

112:   for (i = 0; i < user->m; ++i) v *= z;
113:   return v;
114: }

116: PetscScalar funcADer(PetscScalar z, AppCtx *user)
117: {
118:   PetscScalar v = 1.0;
119:   PetscInt    i;

121:   for (i = 0; i < user->m-1; ++i) v *= z;
122:   return (PetscScalar)user->m*v;
123: }

125: /*
126:    FormFunctionLocal - Evaluates nonlinear function, F(x).
127: */
128: PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
129: {
130:   DM             coordDA;
131:   Vec            coordinates;
132:   DMDACoor2d     **coords;
133:   PetscScalar    u, ux, uy, uxx, uyy;
134:   PetscReal      D, K, hx, hy, hxdhy, hydhx;
135:   PetscInt       i,j;

139:   D     = user->D;
140:   K     = user->K;
141:   hx    = 1.0/(PetscReal)(info->mx-1);
142:   hy    = 1.0/(PetscReal)(info->my-1);
143:   hxdhy = hx/hy;
144:   hydhx = hy/hx;
145:   /*
146:      Compute function over the locally owned part of the grid
147:   */
148:   DMGetCoordinateDM(info->da, &coordDA);
149:   DMGetCoordinates(info->da, &coordinates);
150:   DMDAVecGetArray(coordDA, coordinates, &coords);
151:   for (j=info->ys; j<info->ys+info->ym; j++) {
152:     for (i=info->xs; i<info->xs+info->xm; i++) {
153:       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) f[j][i] = x[j][i];
154:       else {
155:         u       = x[j][i];
156:         ux      = (x[j][i+1] - x[j][i])/hx;
157:         uy      = (x[j+1][i] - x[j][i])/hy;
158:         uxx     = (2.0*u - x[j][i-1] - x[j][i+1])*hydhx;
159:         uyy     = (2.0*u - x[j-1][i] - x[j+1][i])*hxdhy;
160:         f[j][i] = D*(uxx + uyy) - (K*funcA(x[j][i], user)*PetscSqrtScalar(ux*ux + uy*uy) + funcU(&coords[j][i]))*hx*hy;
161:         if (PetscIsInfOrNanScalar(f[j][i])) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FP, "Invalid residual: %g", PetscRealPart(f[j][i]));
162:       }
163:     }
164:   }
165:   DMDAVecRestoreArray(coordDA, coordinates, &coords);
166:   PetscLogFlops(11*info->ym*info->xm);
167:   return(0);
168: }

170: /*
171:    FormJacobianLocal - Evaluates Jacobian matrix.
172: */
173: PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,PetscScalar **x,Mat jac,AppCtx *user)
174: {
175:   MatStencil     col[5], row;
176:   PetscScalar    D, K, A, v[5], hx, hy, hxdhy, hydhx, ux, uy;
177:   PetscReal      normGradZ;
178:   PetscInt       i, j,k;

182:   D     = user->D;
183:   K     = user->K;
184:   hx    = 1.0/(PetscReal)(info->mx-1);
185:   hy    = 1.0/(PetscReal)(info->my-1);
186:   hxdhy = hx/hy;
187:   hydhx = hy/hx;

189:   /*
190:      Compute entries for the locally owned part of the Jacobian.
191:       - Currently, all PETSc parallel matrix formats are partitioned by
192:         contiguous chunks of rows across the processors.
193:       - Each processor needs to insert only elements that it owns
194:         locally (but any non-local elements will be sent to the
195:         appropriate processor during matrix assembly).
196:       - Here, we set all entries for a particular row at once.
197:       - We can set matrix entries either using either
198:         MatSetValuesLocal() or MatSetValues(), as discussed above.
199:   */
200:   for (j=info->ys; j<info->ys+info->ym; j++) {
201:     for (i=info->xs; i<info->xs+info->xm; i++) {
202:       row.j = j; row.i = i;
203:       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
204:         /* boundary points */
205:         v[0] = 1.0;
206:         MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
207:       } else {
208:         /* interior grid points */
209:         ux        = (x[j][i+1] - x[j][i])/hx;
210:         uy        = (x[j+1][i] - x[j][i])/hy;
211:         normGradZ = PetscRealPart(PetscSqrtScalar(ux*ux + uy*uy));
212:         /* PetscPrintf(PETSC_COMM_SELF, "i: %d j: %d normGradZ: %g\n", i, j, normGradZ); */
213:         if (normGradZ < 1.0e-8) normGradZ = 1.0e-8;
214:         A = funcA(x[j][i], user);

216:         v[0] = -D*hxdhy;                                                                          col[0].j = j - 1; col[0].i = i;
217:         v[1] = -D*hydhx;                                                                          col[1].j = j;     col[1].i = i-1;
218:         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;
219:         v[3] = -D*hydhx + K*A*hx*hy/(2.0*normGradZ);                                              col[3].j = j;     col[3].i = i+1;
220:         v[4] = -D*hxdhy + K*A*hx*hy/(2.0*normGradZ);                                              col[4].j = j + 1; col[4].i = i;
221:         for (k = 0; k < 5; ++k) {
222:           if (PetscIsInfOrNanScalar(v[k])) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FP, "Invalid residual: %g", PetscRealPart(v[k]));
223:         }
224:         MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
225:       }
226:     }
227:   }

229:   /*
230:      Assemble matrix, using the 2-step process:
231:        MatAssemblyBegin(), MatAssemblyEnd().
232:   */
233:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
234:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
235:   /*
236:      Tell the matrix we will never add a new nonzero location to the
237:      matrix. If we do, it will generate an error.
238:   */
239:   MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
240:   return(0);
241: }