Actual source code: eptorsion2.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  1: /* Program usage: mpirun -np <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    - sysem 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: TaoGetConvergedReason(); 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;


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


 79: int main(int argc, char **argv)
 80: {
 81:     PetscErrorCode     ierr;
 82:     Vec                x;
 83:     Mat                H;
 84:     PetscInt           Nx, Ny;
 85:     Tao                tao;
 86:     TaoConvergedReason reason;
 87:     PetscBool          flg;
 88:     KSP                ksp;
 89:     PC                 pc;
 90:     AppCtx             user;

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

 94:     /* Specify default dimension of the problem */
 95:     user.param = 5.0; user.mx = 10; user.my = 10;
 96:     Nx = Ny = PETSC_DECIDE;

 98:     /* Check for any command line arguments that override defaults */
 99:     PetscOptionsGetReal(NULL,"-par",&user.param,&flg);
100:     PetscOptionsGetInt(NULL,"-my",&user.my,&flg);
101:     PetscOptionsGetInt(NULL,"-mx",&user.mx,&flg);

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

106:     /* Set up distributed array */
107:     DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.mx,user.my,Nx,Ny,1,1,NULL,NULL,
108:                         &user.dm);

110:     /* Create vectors */
111:     DMCreateGlobalVector(user.dm,&x);

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

115:     /* Create Hessian */
116:     DMCreateMatrix(user.dm,&H);
117:     MatSetOption(H,MAT_SYMMETRIC,PETSC_TRUE);

119:     /* The TAO code begins here */

121:     /* Create TAO solver and set desired solution method */
122:     TaoCreate(PETSC_COMM_WORLD,&tao);
123:     TaoSetType(tao,TAOCG);

125:     /* Set initial solution guess */
126:     FormInitialGuess(&user,x);
127:     TaoSetInitialVector(tao,x);

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

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


135:     /* Check for any TAO command line options */
136:     TaoSetFromOptions(tao);

138:     TaoGetKSP(tao,&ksp);
139:     if (ksp) {
140:       KSPGetPC(ksp,&pc);
141:       PCSetType(pc,PCNONE);
142:     }

144:     /* SOLVE THE APPLICATION */
145:     TaoSolve(tao);

147:     /* Get information on termination */
148:     TaoGetConvergedReason(tao,&reason);
149:     if (reason <= 0){
150:         ierr=PetscPrintf(MPI_COMM_WORLD, "Try another method! \n");
151:     }

153:     /* Free TAO data structures */
154:     TaoDestroy(&tao);

156:     /* Free PETSc data structures */
157:     VecDestroy(&x);
158:     MatDestroy(&H);

160:     VecDestroy(&user.localX);
161:     DMDestroy(&user.dm);

163:     PetscFinalize();
164:     return 0;
165: }


168: /* ------------------------------------------------------------------- */
171: /*
172:     FormInitialGuess - Computes an initial approximation to the solution.

174:     Input Parameters:
175: .   user - user-defined application context
176: .   X    - vector

178:     Output Parameters:
179:     X    - vector
180: */
181: PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
182: {
184:   PetscInt       i, j, k, mx = user->mx, my = user->my;
185:   PetscInt       xs, ys, xm, ym, gxm, gym, gxs, gys, xe, ye;
186:   PetscReal      hx = 1.0/(mx+1), hy = 1.0/(my+1), temp, val;

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

193:   /* Compute initial guess over locally owned part of mesh */
194:   xe = xs+xm;
195:   ye = ys+ym;
196:   for (j=ys; j<ye; j++) {  /*  for (j=0; j<my; j++) */
197:     temp = PetscMin(j+1,my-j)*hy;
198:     for (i=xs; i<xe; i++) {  /*  for (i=0; i<mx; i++) */
199:       k   = (j-gys)*gxm + i-gxs;
200:       val = PetscMin((PetscMin(i+1,mx-i))*hx,temp);
201:       VecSetValuesLocal(X,1,&k,&val,ADD_VALUES);
202:     }
203:   }
204:   VecAssemblyBegin(X);
205:   VecAssemblyEnd(X);
206:   return(0);
207: }


210: /* ------------------------------------------------------------------ */
213: /*
214:    FormFunctionGradient - Evaluates the function and corresponding gradient.

216:    Input Parameters:
217:    tao - the Tao context
218:    X   - the input vector
219:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

221:    Output Parameters:
222:    f   - the newly evaluated function
223:    G   - the newly evaluated gradient
224: */
225: PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f,Vec G,void *ptr){

227:   AppCtx         *user = (AppCtx *)ptr;
229:   PetscInt       i,j,k,ind;
230:   PetscInt       xe,ye,xsm,ysm,xep,yep;
231:   PetscInt       xs, ys, xm, ym, gxm, gym, gxs, gys;
232:   PetscInt       mx = user->mx, my = user->my;
233:   PetscReal      three = 3.0, zero = 0.0, *x, floc, cdiv3 = user->param/three;
234:   PetscReal      p5 = 0.5, area, val, flin, fquad;
235:   PetscReal      v,vb,vl,vr,vt,dvdx,dvdy;
236:   PetscReal      hx = 1.0/(user->mx + 1);
237:   PetscReal      hy = 1.0/(user->my + 1);
238:   Vec            localX = user->localX;

241:   /* Initialize */
242:   flin = fquad = zero;

244:   VecSet(G, zero);
245:   /*
246:      Scatter ghost points to local vector,using the 2-step process
247:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
248:      By placing code between these two statements, computations can be
249:      done while messages are in transition.
250:   */
251:   DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);
252:   DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);

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

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

261:   /* Set local loop dimensions */
262:   xe = xs+xm;
263:   ye = ys+ym;
264:   if (xs == 0)  xsm = xs-1;
265:   else          xsm = xs;
266:   if (ys == 0)  ysm = ys-1;
267:   else          ysm = ys;
268:   if (xe == mx) xep = xe+1;
269:   else          xep = xe;
270:   if (ye == my) yep = ye+1;
271:   else          yep = ye;

273:   /* Compute local gradient contributions over the lower triangular elements */
274:   for (j=ysm; j<ye; j++) {  /*  for (j=-1; j<my; j++) */
275:     for (i=xsm; i<xe; i++) {  /*   for (i=-1; i<mx; i++) */
276:       k = (j-gys)*gxm + i-gxs;
277:       v = zero;
278:       vr = zero;
279:       vt = zero;
280:       if (i >= 0 && j >= 0) v = x[k];
281:       if (i < mx-1 && j > -1) vr = x[k+1];
282:       if (i > -1 && j < my-1) vt = x[k+gxm];
283:       dvdx = (vr-v)/hx;
284:       dvdy = (vt-v)/hy;
285:       if (i != -1 && j != -1) {
286:         ind = k; val = - dvdx/hx - dvdy/hy - cdiv3;
287:         VecSetValuesLocal(G,1,&k,&val,ADD_VALUES);
288:       }
289:       if (i != mx-1 && j != -1) {
290:         ind = k+1; val =  dvdx/hx - cdiv3;
291:         VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES);
292:       }
293:       if (i != -1 && j != my-1) {
294:         ind = k+gxm; val = dvdy/hy - cdiv3;
295:         VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES);
296:       }
297:       fquad += dvdx*dvdx + dvdy*dvdy;
298:       flin -= cdiv3 * (v + vr + vt);
299:     }
300:   }

302:   /* Compute local gradient contributions over the upper triangular elements */
303:   for (j=ys; j<yep; j++) { /*  for (j=0; j<=my; j++) */
304:     for (i=xs; i<xep; i++) {  /*   for (i=0; i<=mx; i++) */
305:       k = (j-gys)*gxm + i-gxs;
306:       vb = zero;
307:       vl = zero;
308:       v  = zero;
309:       if (i < mx && j > 0) vb = x[k-gxm];
310:       if (i > 0 && j < my) vl = x[k-1];
311:       if (i < mx && j < my) v = x[k];
312:       dvdx = (v-vl)/hx;
313:       dvdy = (v-vb)/hy;
314:       if (i != mx && j != 0) {
315:         ind = k-gxm; val = - dvdy/hy - cdiv3;
316:         VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES);
317:       }
318:       if (i != 0 && j != my) {
319:         ind = k-1; val =  - dvdx/hx - cdiv3;
320:         VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES);
321:       }
322:       if (i != mx && j != my) {
323:         ind = k; val =  dvdx/hx + dvdy/hy - cdiv3;
324:         VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES);
325:       }
326:       fquad += dvdx*dvdx + dvdy*dvdy;
327:       flin -= cdiv3 * (vb + vl + v);
328:     }
329:   }


332:   /* Restore vector */
333:   VecRestoreArray(localX,&x);

335:   /* Assemble gradient vector */
336:   VecAssemblyBegin(G);
337:   VecAssemblyEnd(G);

339:   /* Scale the gradient */
340:   area = p5*hx*hy;
341:   floc = area * (p5 * fquad + flin);
342:   VecScale(G, area);

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

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

349:   return(0);
350: }



356: PetscErrorCode FormHessian(Tao tao, Vec X, Mat A, Mat Hpre, void*ctx)
357: {
358:   AppCtx         *user= (AppCtx*) ctx;
360:   PetscInt       i,j,k;
361:   PetscInt       col[5],row;
362:   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym;
363:   PetscReal      v[5];
364:   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;

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

368:   /*
369:      Get local grid boundaries
370:   */

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

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

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

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

382:       k=0;
383:       if (j>gys){
384:         v[k]=-2*hyhy; col[k]=row - gxm; k++;
385:       }

387:       if (i>gxs){
388:         v[k]= -2*hxhx; col[k]=row - 1; k++;
389:       }

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

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

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

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

403:     }
404:   }
405:   /*
406:      Assemble matrix, using the 2-step process:
407:      MatAssemblyBegin(), MatAssemblyEnd().
408:      By placing code between these two statements, computations can be
409:      done while messages are in transition.
410:   */
411:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
412:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
413:   /*
414:     Tell the matrix we will never add a new nonzero location to the
415:     matrix. If we do it will generate an error.
416:   */
417:   MatScale(A,area);
418:   MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
419:   MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);
420:   PetscLogFlops(9*xm*ym+49*xm);
421:   return(0);
422: }