Actual source code: eptorsion2.c

  1: /* Program usage: mpiexec -n <proc> eptorsion2 [-help] [all TAO options] */

  3: /* ----------------------------------------------------------------------

  5:   Elastic-plastic torsion problem.

  7:   The elastic plastic torsion problem arises from the determination
  8:   of the stress field on an infinitely long cylindrical bar, which is
  9:   equivalent to the solution of the following problem:

 11:   min{ .5 * integral(||gradient(v(x))||^2 dx) - C * integral(v(x) dx)}

 13:   where C is the torsion angle per unit length.

 15:   The uniprocessor version of this code is eptorsion1.c; the Fortran
 16:   version of this code is eptorsion2f.F.

 18:   This application solves the problem without calculating hessians
 19: ---------------------------------------------------------------------- */

 21: /*
 22:   Include "petsctao.h" so that we can use TAO solvers.  Note that this
 23:   file automatically includes files for lower-level support, such as those
 24:   provided by the PETSc library:
 25:      petsc.h       - base PETSc routines   petscvec.h - vectors
 26:      petscsys.h    - system routines        petscmat.h - matrices
 27:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 28:      petscviewer.h - viewers               petscpc.h  - preconditioners
 29:   Include "petscdmda.h" so that we can use distributed arrays (DMs) for managing
 30:   the parallel mesh.
 31: */

 33: #include <petsctao.h>
 34: #include <petscdmda.h>

 36: static  char help[] =
 37: "Demonstrates use of the TAO package to solve \n\
 38: unconstrained minimization problems in parallel.  This example is based on \n\
 39: the Elastic-Plastic Torsion (dept) problem from the MINPACK-2 test suite.\n\
 40: The command line options are:\n\
 41:   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
 42:   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
 43:   -par <param>, where <param> = angle of twist per unit length\n\n";

 45: /*T
 46:    Concepts: TAO^Solving an unconstrained minimization problem
 47:    Routines: TaoCreate(); TaoSetType();
 48:    Routines: TaoSetInitialVector();
 49:    Routines: TaoSetObjectiveAndGradientRoutine();
 50:    Routines: TaoSetHessianRoutine(); TaoSetFromOptions();
 51:    Routines: TaoSolve();
 52:    Routines: TaoDestroy();
 53:    Processors: n
 54: T*/

 56: /*
 57:    User-defined application context - contains data needed by the
 58:    application-provided call-back routines, FormFunction() and
 59:    FormGradient().
 60: */
 61: typedef struct {
 62:   /* parameters */
 63:    PetscInt      mx, my;         /* global discretization in x- and y-directions */
 64:    PetscReal     param;          /* nonlinearity parameter */

 66:   /* work space */
 67:    Vec           localX;         /* local vectors */
 68:    DM            dm;             /* distributed array data structure */
 69: } AppCtx;

 71: PetscErrorCode FormInitialGuess(AppCtx*, Vec);
 72: PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
 73: PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);

 75: int main(int argc, char **argv)
 76: {
 77:     PetscErrorCode     ierr;
 78:     Vec                x;
 79:     Mat                H;
 80:     PetscInt           Nx, Ny;
 81:     Tao                tao;
 82:     PetscBool          flg;
 83:     KSP                ksp;
 84:     PC                 pc;
 85:     AppCtx             user;

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

 89:     /* Specify default dimension of the problem */
 90:     user.param = 5.0; user.mx = 10; user.my = 10;
 91:     Nx = Ny = PETSC_DECIDE;

 93:     /* Check for any command line arguments that override defaults */
 94:     PetscOptionsGetReal(NULL,NULL,"-par",&user.param,&flg);
 95:     PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg);
 96:     PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg);

 98:     PetscPrintf(PETSC_COMM_WORLD,"\n---- Elastic-Plastic Torsion Problem -----\n");
 99:     PetscPrintf(PETSC_COMM_WORLD,"mx: %D     my: %D   \n\n",user.mx,user.my);

101:     /* Set up distributed array */
102:     DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.mx,user.my,Nx,Ny,1,1,NULL,NULL,&user.dm);
103:     DMSetFromOptions(user.dm);
104:     DMSetUp(user.dm);

106:     /* Create vectors */
107:     DMCreateGlobalVector(user.dm,&x);

109:     DMCreateLocalVector(user.dm,&user.localX);

111:     /* Create Hessian */
112:     DMCreateMatrix(user.dm,&H);
113:     MatSetOption(H,MAT_SYMMETRIC,PETSC_TRUE);

115:     /* The TAO code begins here */

117:     /* Create TAO solver and set desired solution method */
118:     TaoCreate(PETSC_COMM_WORLD,&tao);
119:     TaoSetType(tao,TAOCG);

121:     /* Set initial solution guess */
122:     FormInitialGuess(&user,x);
123:     TaoSetInitialVector(tao,x);

125:     /* Set routine for function and gradient evaluation */
126:     TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void *)&user);

128:     TaoSetHessianRoutine(tao,H,H,FormHessian,(void*)&user);

130:     /* Check for any TAO command line options */
131:     TaoSetFromOptions(tao);

133:     TaoGetKSP(tao,&ksp);
134:     if (ksp) {
135:       KSPGetPC(ksp,&pc);
136:       PCSetType(pc,PCNONE);
137:     }

139:     /* SOLVE THE APPLICATION */
140:     TaoSolve(tao);

142:     /* Free TAO data structures */
143:     TaoDestroy(&tao);

145:     /* Free PETSc data structures */
146:     VecDestroy(&x);
147:     MatDestroy(&H);

149:     VecDestroy(&user.localX);
150:     DMDestroy(&user.dm);

152:     PetscFinalize();
153:     return 0;
154: }

156: /* ------------------------------------------------------------------- */
157: /*
158:     FormInitialGuess - Computes an initial approximation to the solution.

160:     Input Parameters:
161: .   user - user-defined application context
162: .   X    - vector

164:     Output Parameters:
165:     X    - vector
166: */
167: PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
168: {
170:   PetscInt       i, j, k, mx = user->mx, my = user->my;
171:   PetscInt       xs, ys, xm, ym, gxm, gym, gxs, gys, xe, ye;
172:   PetscReal      hx = 1.0/(mx+1), hy = 1.0/(my+1), temp, val;

175:   /* Get local mesh boundaries */
176:   DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
177:   DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);

179:   /* Compute initial guess over locally owned part of mesh */
180:   xe = xs+xm;
181:   ye = ys+ym;
182:   for (j=ys; j<ye; j++) {  /*  for (j=0; j<my; j++) */
183:     temp = PetscMin(j+1,my-j)*hy;
184:     for (i=xs; i<xe; i++) {  /*  for (i=0; i<mx; i++) */
185:       k   = (j-gys)*gxm + i-gxs;
186:       val = PetscMin((PetscMin(i+1,mx-i))*hx,temp);
187:       VecSetValuesLocal(X,1,&k,&val,ADD_VALUES);
188:     }
189:   }
190:   VecAssemblyBegin(X);
191:   VecAssemblyEnd(X);
192:   return(0);
193: }

195: /* ------------------------------------------------------------------ */
196: /*
197:    FormFunctionGradient - Evaluates the function and corresponding gradient.

199:    Input Parameters:
200:    tao - the Tao context
201:    X   - the input vector
202:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

204:    Output Parameters:
205:    f   - the newly evaluated function
206:    G   - the newly evaluated gradient
207: */
208: PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f,Vec G,void *ptr)
209: {
210:   AppCtx         *user = (AppCtx *)ptr;
212:   PetscInt       i,j,k,ind;
213:   PetscInt       xe,ye,xsm,ysm,xep,yep;
214:   PetscInt       xs, ys, xm, ym, gxm, gym, gxs, gys;
215:   PetscInt       mx = user->mx, my = user->my;
216:   PetscReal      three = 3.0, zero = 0.0, *x, floc, cdiv3 = user->param/three;
217:   PetscReal      p5 = 0.5, area, val, flin, fquad;
218:   PetscReal      v,vb,vl,vr,vt,dvdx,dvdy;
219:   PetscReal      hx = 1.0/(user->mx + 1);
220:   PetscReal      hy = 1.0/(user->my + 1);
221:   Vec            localX = user->localX;

224:   /* Initialize */
225:   flin = fquad = zero;

227:   VecSet(G, zero);
228:   /*
229:      Scatter ghost points to local vector,using the 2-step process
230:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
231:      By placing code between these two statements, computations can be
232:      done while messages are in transition.
233:   */
234:   DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);
235:   DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);

237:   /* Get pointer to vector data */
238:   VecGetArray(localX,&x);

240:   /* Get local mesh boundaries */
241:   DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
242:   DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);

244:   /* Set local loop dimensions */
245:   xe = xs+xm;
246:   ye = ys+ym;
247:   if (xs == 0)  xsm = xs-1;
248:   else          xsm = xs;
249:   if (ys == 0)  ysm = ys-1;
250:   else          ysm = ys;
251:   if (xe == mx) xep = xe+1;
252:   else          xep = xe;
253:   if (ye == my) yep = ye+1;
254:   else          yep = ye;

256:   /* Compute local gradient contributions over the lower triangular elements */
257:   for (j=ysm; j<ye; j++) {  /*  for (j=-1; j<my; j++) */
258:     for (i=xsm; i<xe; i++) {  /*   for (i=-1; i<mx; i++) */
259:       k = (j-gys)*gxm + i-gxs;
260:       v = zero;
261:       vr = zero;
262:       vt = zero;
263:       if (i >= 0 && j >= 0) v = x[k];
264:       if (i < mx-1 && j > -1) vr = x[k+1];
265:       if (i > -1 && j < my-1) vt = x[k+gxm];
266:       dvdx = (vr-v)/hx;
267:       dvdy = (vt-v)/hy;
268:       if (i != -1 && j != -1) {
269:         ind = k; val = - dvdx/hx - dvdy/hy - cdiv3;
270:         VecSetValuesLocal(G,1,&k,&val,ADD_VALUES);
271:       }
272:       if (i != mx-1 && j != -1) {
273:         ind = k+1; val =  dvdx/hx - cdiv3;
274:         VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES);
275:       }
276:       if (i != -1 && j != my-1) {
277:         ind = k+gxm; val = dvdy/hy - cdiv3;
278:         VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES);
279:       }
280:       fquad += dvdx*dvdx + dvdy*dvdy;
281:       flin -= cdiv3 * (v + vr + vt);
282:     }
283:   }

285:   /* Compute local gradient contributions over the upper triangular elements */
286:   for (j=ys; j<yep; j++) { /*  for (j=0; j<=my; j++) */
287:     for (i=xs; i<xep; i++) {  /*   for (i=0; i<=mx; i++) */
288:       k = (j-gys)*gxm + i-gxs;
289:       vb = zero;
290:       vl = zero;
291:       v  = zero;
292:       if (i < mx && j > 0) vb = x[k-gxm];
293:       if (i > 0 && j < my) vl = x[k-1];
294:       if (i < mx && j < my) v = x[k];
295:       dvdx = (v-vl)/hx;
296:       dvdy = (v-vb)/hy;
297:       if (i != mx && j != 0) {
298:         ind = k-gxm; val = - dvdy/hy - cdiv3;
299:         VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES);
300:       }
301:       if (i != 0 && j != my) {
302:         ind = k-1; val =  - dvdx/hx - cdiv3;
303:         VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES);
304:       }
305:       if (i != mx && j != my) {
306:         ind = k; val =  dvdx/hx + dvdy/hy - cdiv3;
307:         VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES);
308:       }
309:       fquad += dvdx*dvdx + dvdy*dvdy;
310:       flin -= cdiv3 * (vb + vl + v);
311:     }
312:   }

314:   /* Restore vector */
315:   VecRestoreArray(localX,&x);

317:   /* Assemble gradient vector */
318:   VecAssemblyBegin(G);
319:   VecAssemblyEnd(G);

321:   /* Scale the gradient */
322:   area = p5*hx*hy;
323:   floc = area * (p5 * fquad + flin);
324:   VecScale(G, area);

326:   /* Sum function contributions from all processes */
327:   (PetscErrorCode)MPI_Allreduce((void*)&floc,(void*)f,1,MPIU_REAL,MPIU_SUM,MPI_COMM_WORLD);

329:   PetscLogFlops((ye-ysm)*(xe-xsm)*20+(xep-xs)*(yep-ys)*16);

331:   return(0);
332: }

334: PetscErrorCode FormHessian(Tao tao, Vec X, Mat A, Mat Hpre, void*ctx)
335: {
336:   AppCtx         *user= (AppCtx*) ctx;
338:   PetscInt       i,j,k;
339:   PetscInt       col[5],row;
340:   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym;
341:   PetscReal      v[5];
342:   PetscReal      hx=1.0/(user->mx+1), hy=1.0/(user->my+1), hxhx=1.0/(hx*hx), hyhy=1.0/(hy*hy), area=0.5*hx*hy;

344:   /* Compute the quadratic term in the objective function */

346:   /*
347:      Get local grid boundaries
348:   */

351:   DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
352:   DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);

354:   for (j=ys; j<ys+ym; j++) {

356:     for (i=xs; i< xs+xm; i++) {

358:       row=(j-gys)*gxm + (i-gxs);

360:       k=0;
361:       if (j>gys) {
362:         v[k]=-2*hyhy; col[k]=row - gxm; k++;
363:       }

365:       if (i>gxs) {
366:         v[k]= -2*hxhx; col[k]=row - 1; k++;
367:       }

369:       v[k]= 4.0*(hxhx+hyhy); col[k]=row; k++;

371:       if (i+1 < gxs+gxm) {
372:         v[k]= -2.0*hxhx; col[k]=row+1; k++;
373:       }

375:       if (j+1 <gys+gym) {
376:         v[k]= -2*hyhy; col[k] = row+gxm; k++;
377:       }

379:       MatSetValuesLocal(A,1,&row,k,col,v,INSERT_VALUES);

381:     }
382:   }
383:   /*
384:      Assemble matrix, using the 2-step process:
385:      MatAssemblyBegin(), MatAssemblyEnd().
386:      By placing code between these two statements, computations can be
387:      done while messages are in transition.
388:   */
389:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
390:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
391:   /*
392:     Tell the matrix we will never add a new nonzero location to the
393:     matrix. If we do it will generate an error.
394:   */
395:   MatScale(A,area);
396:   MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
397:   MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);
398:   PetscLogFlops(9*xm*ym+49*xm);
399:   return(0);
400: }

402: /*TEST

404:    build:
405:       requires: !complex

407:    test:
408:       suffix: 1
409:       nsize: 2
410:       args: -tao_smonitor -tao_type nls -mx 16 -my 16 -tao_gatol 1.e-4

412: TEST*/