Actual source code: ex46.c

petsc-3.7.5 2017-01-01
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*);

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

 42:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 43:      Initialize program
 44:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 46:   PetscInitialize(&argc,&argv,(char*)0,help);

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

 60:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 61:      Create distributed array (DMDA) to manage parallel grid and vectors
 62:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 63:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&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(0);
100: }

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

111: PetscScalar funcA(PetscScalar z, AppCtx *user)
112: {
113:   PetscScalar v = 1.0;
114:   PetscInt    i;

116:   for (i = 0; i < user->m; ++i) v *= z;
117:   return v;
118: }

122: PetscScalar funcADer(PetscScalar z, AppCtx *user)
123: {
124:   PetscScalar v = 1.0;
125:   PetscInt    i;

127:   for (i = 0; i < user->m-1; ++i) v *= z;
128:   return (PetscScalar)user->m*v;
129: }

133: /*
134:    FormFunctionLocal - Evaluates nonlinear function, F(x).
135: */
136: PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
137: {
138:   DM             coordDA;
139:   Vec            coordinates;
140:   DMDACoor2d     **coords;
141:   PetscScalar    u, ux, uy, uxx, uyy;
142:   PetscReal      D, K, hx, hy, hxdhy, hydhx;
143:   PetscInt       i,j;

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

180: /*
181:    FormJacobianLocal - Evaluates Jacobian matrix.
182: */
183: PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,PetscScalar **x,Mat jac,AppCtx *user)
184: {
185:   MatStencil     col[5], row;
186:   PetscScalar    D, K, A, v[5], hx, hy, hxdhy, hydhx, ux, uy;
187:   PetscReal      normGradZ;
188:   PetscInt       i, j,k;

192:   D     = user->D;
193:   K     = user->K;
194:   hx    = 1.0/(PetscReal)(info->mx-1);
195:   hy    = 1.0/(PetscReal)(info->my-1);
196:   hxdhy = hx/hy;
197:   hydhx = hy/hx;

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

226:         v[0] = -D*hxdhy;                                                                          col[0].j = j - 1; col[0].i = i;
227:         v[1] = -D*hydhx;                                                                          col[1].j = j;     col[1].i = i-1;
228:         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;
229:         v[3] = -D*hydhx + K*A*hx*hy/(2.0*normGradZ);                                              col[3].j = j;     col[3].i = i+1;
230:         v[4] = -D*hxdhy + K*A*hx*hy/(2.0*normGradZ);                                              col[4].j = j + 1; col[4].i = i;
231:         for (k = 0; k < 5; ++k) {
232:           if (PetscIsInfOrNanScalar(v[k])) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FP, "Invalid residual: %g", PetscRealPart(v[k]));
233:         }
234:         MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
235:       }
236:     }
237:   }

239:   /*
240:      Assemble matrix, using the 2-step process:
241:        MatAssemblyBegin(), MatAssemblyEnd().
242:   */
243:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
244:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
245:   /*
246:      Tell the matrix we will never add a new nonzero location to the
247:      matrix. If we do, it will generate an error.
248:   */
249:   MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
250:   return(0);
251: }