Actual source code: eptorsion2.c

petsc-3.7.7 2017-09-25
Report Typos and Errors
  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    - 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: 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:     PetscBool          flg;
 87:     KSP                ksp;
 88:     PC                 pc;
 89:     AppCtx             user;

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

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

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

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

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

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

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

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

118:     /* The TAO code begins here */

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

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

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

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


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

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

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

146:     /* Free TAO data structures */
147:     TaoDestroy(&tao);

149:     /* Free PETSc data structures */
150:     VecDestroy(&x);
151:     MatDestroy(&H);

153:     VecDestroy(&user.localX);
154:     DMDestroy(&user.dm);

156:     PetscFinalize();
157:     return 0;
158: }


161: /* ------------------------------------------------------------------- */
164: /*
165:     FormInitialGuess - Computes an initial approximation to the solution.

167:     Input Parameters:
168: .   user - user-defined application context
169: .   X    - vector

171:     Output Parameters:
172:     X    - vector
173: */
174: PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
175: {
177:   PetscInt       i, j, k, mx = user->mx, my = user->my;
178:   PetscInt       xs, ys, xm, ym, gxm, gym, gxs, gys, xe, ye;
179:   PetscReal      hx = 1.0/(mx+1), hy = 1.0/(my+1), temp, val;

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

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


203: /* ------------------------------------------------------------------ */
206: /*
207:    FormFunctionGradient - Evaluates the function and corresponding gradient.

209:    Input Parameters:
210:    tao - the Tao context
211:    X   - the input vector
212:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

214:    Output Parameters:
215:    f   - the newly evaluated function
216:    G   - the newly evaluated gradient
217: */
218: PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f,Vec G,void *ptr){

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

234:   /* Initialize */
235:   flin = fquad = zero;

237:   VecSet(G, zero);
238:   /*
239:      Scatter ghost points to local vector,using the 2-step process
240:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
241:      By placing code between these two statements, computations can be
242:      done while messages are in transition.
243:   */
244:   DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);
245:   DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);

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

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

254:   /* Set local loop dimensions */
255:   xe = xs+xm;
256:   ye = ys+ym;
257:   if (xs == 0)  xsm = xs-1;
258:   else          xsm = xs;
259:   if (ys == 0)  ysm = ys-1;
260:   else          ysm = ys;
261:   if (xe == mx) xep = xe+1;
262:   else          xep = xe;
263:   if (ye == my) yep = ye+1;
264:   else          yep = ye;

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

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


325:   /* Restore vector */
326:   VecRestoreArray(localX,&x);

328:   /* Assemble gradient vector */
329:   VecAssemblyBegin(G);
330:   VecAssemblyEnd(G);

332:   /* Scale the gradient */
333:   area = p5*hx*hy;
334:   floc = area * (p5 * fquad + flin);
335:   VecScale(G, area);

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

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

342:   return(0);
343: }



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

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

361:   /*
362:      Get local grid boundaries
363:   */

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

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

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

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

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

380:       if (i>gxs){
381:         v[k]= -2*hxhx; col[k]=row - 1; k++;
382:       }

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

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

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

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

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