Actual source code: eptorsion2f.F

  1: ! "$Id$";
  2: !
  3: !  Program usage: mpirun -np <proc> eptorsion2f [all TAO options] 
  4: !
  5: !  Description:  This example demonstrates use of the TAO package to solve
  6: !  unconstrained minimization problems in parallel.  This example is based 
  7: !  on the Elastic-Plastic Torsion (dept) problem from the MINPACK-2 test suite.
  8: !  The command line options are:
  9: !    -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction
 10: !    -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction
 11: !    -par <param>, where <param> = angle of twist per unit length
 12: !
 13: !/*T
 14: !   Concepts: TAO - Solving an unconstrained minimization problem
 15: !   Routines: TaoInitialize(); TaoFinalize(); 
 16: !   Routines: TaoCreate(); TaoDestroy();
 17: !   Routines: TaoApplicationCreate(); TaoAppDestroy();
 18: !   Routines: TaoAppSetObjectiveAndGradientRoutine();
 19: !   Routines: TaoAppSetHessianMat(); TaoAppSetHessianRoutine();
 20: !   Routines: TaoSetOptions();
 21: !   Routines: TaoAppSetInitialSolutionVec(); TaoSolveApplication(); TaoDestroy();
 22: !   Routines: TaoGetSolutionStatus();
 23: !   Processors: n
 24: !T*/
 25: !
 26: ! ---------------------------------------------------------------------- 
 27: !
 28: !  Elastic-plastic torsion problem.  
 29: !
 30: !  The elastic plastic torsion problem arises from the determination 
 31: !  of the stress field on an infinitely long cylindrical bar, which is
 32: !  equivalent to the solution of the following problem:
 33: !     min{ .5 * integral(||gradient(v(x))||^2 dx) - C * integral(v(x) dx)}
 34: !  where C is the torsion angle per unit length.
 35: !
 36: !  The C version of this code is eptorsion2.c
 37: !
 38: ! ---------------------------------------------------------------------- 

 40:       implicit none
 41: #include "eptorsion2f.h"

 43: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 44: !                   Variable declarations
 45: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 46: !
 47: !  See additional variable declarations in the file eptorsion2f.h
 48: !
 49:       integer          info           ! used to check for functions returning nonzeros
 50:       Vec              x              ! solution vector
 51:       Mat              H              ! hessian matrix
 52:       PetscInt         Nx, Ny         ! number of processes in x- and y- directions
 53:       TAO_SOLVER       tao            ! TAO_SOLVER solver context
 54:       TAO_APPLICATION  torsionapp     ! TAO application context (PETSc)
 55:       TaoTerminateReason reason
 56:       PetscTruth       flg
 57:       PetscInt         iter           ! iteration information
 58:       PetscScalar      ff,gnorm,cnorm,xdiff
 59:       PetscInt         i1
 60:       

 62: !  Note: Any user-defined Fortran routines (such as FormGradient)
 63: !  MUST be declared as external.

 65:       external FormInitialGuess,FormFunctionGradient,ComputeHessian

 67:       i1 = 1

 69: !     Initialize TAO, PETSc  contexts
 70:       call PetscInitialize(PETSC_NULL_CHARACTER,info)
 71:       call TaoInitialize(PETSC_NULL_CHARACTER,info)

 73: !     Specify default parameters 
 74:       param = 5.0d0
 75:       mx = 10
 76:       my = 10
 77:       Nx = PETSC_DECIDE
 78:       Ny = PETSC_DECIDE

 80: !     Check for any command line arguments that might override defaults
 81:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,"-mx",mx,flg,info)
 82:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,"-my",my,flg,info)
 83:       call PetscOptionsGetReal(PETSC_NULL_CHARACTER,"-par",              &
 84:      &                         param,flg,info)

 86:       
 87: !     Set up distributed array and vectors
 88:       call DACreate2d(MPI_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_BOX,      &
 89:      &     mx,my,Nx,Ny,i1,i1,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,      &
 90:      &     da,info)

 92: !     Create vectors
 93:       call DACreateGlobalVector(da,x,info)
 94:       call DACreateLocalVector(da,localX,info)

 96: !     Create Hessian
 97:       call DAGetMatrix(da,MATAIJ,H,info)
 98:       call MatSetOption(H,MAT_SYMMETRIC,PETSC_TRUE,info)

100: !     The TAO code begins here

102: !     Create TAO solver
103:       call TaoCreate(MPI_COMM_WORLD,'tao_cg',tao,info)
104:       call TaoApplicationCreate(MPI_COMM_WORLD,torsionapp,info)

106: !     Set routines for function and gradient evaluation

108: !     TaoAppSetObjectiveAndGradientRoutine is shortened to 31 chars to comply with some compilers
109:       call TaoAppSetObjectiveAndGradientRo(torsionapp,                       &
110:      &     FormFunctionGradient,TAO_NULL_OBJECT,info)
111:       call TaoAppSetHessianMat(torsionapp,H,H,info)
112:       call TaoAppSetHessianRoutine(torsionapp,ComputeHessian,                &
113:      &    TAO_NULL_OBJECT,info)

115: !     Set initial guess
116:       call FormInitialGuess(x,info)
117:       call TaoAppSetInitialSolutionVec(torsionapp,x,info)

119: !     Check for any TAO command line options
120:       call TaoSetOptions(torsionapp, tao,info)


123: !     SOLVE THE APPLICATION
124:       call TaoSolveApplication(torsionapp,tao,info)

126: !     Get information on termination
127:       call TaoGetSolutionStatus(tao,iter,ff,gnorm,cnorm,xdiff,           &
128:      &                          reason,info)
129:       if (reason .lt. 0) then
130:          print *,'TAO did not terminate successfully'
131:       endif

133:       
134: !     Free TAO data structures
135:       call TaoDestroy(tao,info)
136:       call TaoAppDestroy(torsionapp,info);  

138:   
139: !     Free PETSc data structures 
140:       call VecDestroy(x,info)
141:       call VecDestroy(localX,info)
142:       call MatDestroy(H,info)
143:       call DADestroy(da,info)


146: !     Finalize TAO and PETSc
147:       call PetscFinalize(info)
148:       call TaoFinalize(info)

150:       end


153: ! ---------------------------------------------------------------------
154: !
155: !   FormInitialGuess - Computes an initial approximation to the solution.
156: !
157: !   Input Parameters:
158: !   X    - vector
159: !
160: !   Output Parameters:
161: !   X    - vector
162: !   info - error code
163: !
164:       subroutine FormInitialGuess(X,info)
165:       implicit none

167: ! mx, my are defined in eptorsion2f.h
168: #include "eptorsion2f.h"

170: !  Input/output variables:
171:       Vec              X
172:       integer          info

174: !  Local variables:
175:       PetscInt         i, j, k, xe, ye
176:       PetscScalar      temp, val, hx, hy
177:       PetscInt         xs, ys, xm, ym
178:       PetscInt         gxm, gym, gxs, gys
179:       PetscInt         i1

181:       i1 = 1
182:       hx = 1.0d0/(mx + 1)
183:       hy = 1.0d0/(my + 1)

185: !  Get corner information
186:       call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,               &
187:      &                  PETSC_NULL_INTEGER,info) 
188:       call DAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,              &
189:      &                   gxm,gym,PETSC_NULL_INTEGER,info)



193: !  Compute initial guess over locally owned part of mesh
194:       xe = xs+xm
195:       ye = ys+ym
196:       do j=ys,ye-1
197:          temp = min(j+1,my-j)*hy
198:          do i=xs,xe-1
199:             k   = (j-gys)*gxm + i-gxs
200:             val = min((min(i+1,mx-i))*hx,temp)
201:             call VecSetValuesLocal(X,i1,k,val,ADD_VALUES,info)
202:          end do
203:       end do

205:       return
206:       end


209: ! ---------------------------------------------------------------------
210: !
211: !  FormFunctionGradient - Evaluates gradient G(X). 
212: !
213: !  Input Parameters:
214: !  tao   - the TAO_SOLVER context
215: !  X     - input vector
216: !  dummy - optional user-defined context (not used here)
217: !    
218: !  Output Parameters:
219: !  f     - the function value at X
220: !  G     - vector containing the newly evaluated gradient
221: !  info  - error code
222: !
223: !  Notes:
224: !  This routine serves as a wrapper for the lower-level routine
225: !  "ApplicationGradient", where the actual computations are 
226: !  done using the standard Fortran style of treating the local
227: !  input vector data as an array over the local mesh.
228: !
229:       subroutine FormFunctionGradient(taoapp,X,f,G,dummy,info)
230:       implicit none

232: ! da, mx, my, param, localX declared in eptorsion2f.h
233: #include "eptorsion2f.h"

235: !  Input/output variables:
236:       TAO_APPLICATION  taoapp
237:       Vec              X, G
238:       PetscScalar      f
239:       integer          info
240:       PetscInt         dummy

242: !  Declarations for use with local array:


245: ! PETSc's VecGetArray acts differently in Fortran than it does in C.
246: ! Calling VecGetArray((Vec) X, (PetscScalar) x_array(0:1), (PetscOffset) x_index, info)
247: ! will return an array of doubles referenced by x_array offset by x_index.
248: !  i.e.,  to reference the kth element of X, use x_array(k + x_index).
249: ! Notice that by declaring the arrays with range (0:1), we are using the C 0-indexing practice.
250:       PetscScalar      lx_v(0:1)
251:       PetscOffset      lx_i

253: !  Local variables:
254:       PetscScalar      zero, p5, area, cdiv3
255:       PetscScalar      val, flin, fquad,floc
256:       PetscScalar      v, vb, vl, vr, vt, dvdx
257:       PetscScalar      dvdy, hx, hy
258:       PetscInt         xe, ye, xsm, ysm
259:       PetscInt         xep, yep, i, j, k, ind
260:       PetscInt         xs, ys, xm, ym 
261:       PetscInt         gxs, gys, gxm, gym
262:       PetscInt         i1

264:       i1 = 1
265:       info  = 0
266:       cdiv3 = param/3.0d0
267:       zero = 0.0d0
268:       p5   = 0.5d0
269:       hx = 1.0d0/(mx + 1)
270:       hy = 1.0d0/(my + 1)
271:       fquad = zero
272:       flin = zero

274: !  Initialize gradient to zero
275:       call VecSet(G,zero,info)

277: !  Scatter ghost points to local vector
278:       call DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX,info)
279:       call DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX,info)


282: !  Get corner information
283:       call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,               &
284:      &                  PETSC_NULL_INTEGER,info)
285:       call DAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,              &
286:      &                   gxm,gym,PETSC_NULL_INTEGER,info)

288: !  Get pointer to vector data.
289:       call VecGetArray(localX,lx_v,lx_i,info)


292: !  Set local loop dimensions
293:       xe = xs+xm
294:       ye = ys+ym
295:       if (xs .eq. 0) then
296:          xsm = xs-1
297:       else
298:          xsm = xs
299:       endif
300:       if (ys .eq. 0) then
301:          ysm = ys-1
302:       else
303:          ysm = ys
304:       endif
305:       if (xe .eq. mx) then
306:          xep = xe+1
307:       else
308:          xep = xe
309:       endif
310:       if (ye .eq. my) then
311:          yep = ye+1
312:       else
313:          yep = ye
314:       endif

316: !     Compute local gradient contributions over the lower triangular elements
317:      
318:       do j = ysm, ye-1
319:          do i = xsm, xe-1
320:             k  = (j-gys)*gxm + i-gxs
321:             v  = zero
322:             vr = zero
323:             vt = zero
324:             if (i .ge. 0 .and. j .ge. 0)      v = lx_v(lx_i+k)
325:             if (i .lt. mx-1 .and. j .gt. -1) vr = lx_v(lx_i+k+1)
326:             if (i .gt. -1 .and. j .lt. my-1) vt = lx_v(lx_i+k+gxm)
327:             dvdx = (vr-v)/hx
328:             dvdy = (vt-v)/hy
329:             if (i .ne. -1 .and. j .ne. -1) then
330:                ind = k
331:                val = - dvdx/hx - dvdy/hy - cdiv3
332:                call VecSetValuesLocal(G,i1,k,val,ADD_VALUES,info)
333:             endif
334:             if (i .ne. mx-1 .and. j .ne. -1) then
335:                ind = k+1
336:                val =  dvdx/hx - cdiv3
337:                call VecSetValuesLocal(G,i1,ind,val,ADD_VALUES,info)
338:             endif
339:             if (i .ne. -1 .and. j .ne. my-1) then
340:               ind = k+gxm
341:               val = dvdy/hy - cdiv3
342:               call VecSetValuesLocal(G,i1,ind,val,ADD_VALUES,info)
343:             endif
344:             fquad = fquad + dvdx*dvdx + dvdy*dvdy
345:             flin = flin - cdiv3 * (v+vr+vt)
346:          end do
347:       end do

349: !     Compute local gradient contributions over the upper triangular elements

351:       do j = ys, yep-1
352:          do i = xs, xep-1
353:             k  = (j-gys)*gxm + i-gxs
354:             vb = zero
355:             vl = zero
356:             v  = zero
357:             if (i .lt. mx .and. j .gt. 0) vb = lx_v(lx_i+k-gxm)
358:             if (i .gt. 0 .and. j .lt. my) vl = lx_v(lx_i+k-1)
359:             if (i .lt. mx .and. j .lt. my) v = lx_v(lx_i+k)
360:             dvdx = (v-vl)/hx
361:             dvdy = (v-vb)/hy
362:             if (i .ne. mx .and. j .ne. 0) then
363:                ind = k-gxm
364:                val = - dvdy/hy - cdiv3
365:                call VecSetValuesLocal(G,i1,ind,val,ADD_VALUES,info)
366:             endif
367:             if (i .ne. 0 .and. j .ne. my) then
368:                ind = k-1
369:                val =  - dvdx/hx - cdiv3
370:                call VecSetValuesLocal(G,i1,ind,val,ADD_VALUES,info)
371:             endif
372:             if (i .ne. mx .and. j .ne. my) then
373:                ind = k
374:                val =  dvdx/hx + dvdy/hy - cdiv3
375:                call VecSetValuesLocal(G,i1,ind,val,ADD_VALUES,info)
376:             endif
377:             fquad = fquad + dvdx*dvdx + dvdy*dvdy
378:             flin = flin - cdiv3*(vb + vl + v)
379:          end do
380:       end do

382: !  Restore vector
383:       call VecRestoreArray(localX,lx_v,lx_i,info)
384:       if (info .ne. 0) print *,'VecRestoreArray'
385: !  Assemble gradient vector
386:       call VecAssemblyBegin(G,info)
387:       call VecAssemblyEnd(G,info)

389: ! Scale the gradient      
390:       area = p5*hx*hy
391:       floc = area *(p5*fquad+flin)
392:       call VecScale(G,area,info)


395: !  Sum function contributions from all processes
396:       call MPI_Allreduce(floc,f,1,MPI_DOUBLE_PRECISION,MPI_SUM,           &
397:      &                   MPI_COMM_WORLD,info)
398:       if (info .ne. 0) print *,'MPI_Allreduce'
399:       call PetscLogFlops((ye-ysm)*(xe-xsm)*20+(xep-xs)*(yep-ys)*16,       &
400:      &                    info)



404:       return
405:       end




410:       subroutine ComputeHessian(taoapp, X, H, Hpre, flag, dummy, info)
411:       implicit none
412: #include "eptorsion2f.h"      
413:       TAO_APPLICATION taoapp
414:       Vec             X
415:       Mat             H,Hpre
416:       MatStructure    flag
417:       integer         info
418:       PetscInt        dummy

420:       
421:       PetscInt i,j,k
422:       PetscInt col(0:4),row
423:       PetscInt xs,xm,gxs,gxm
424:       PetscInt ys,ym,gys,gym
425:       PetscScalar v(0:4)
426:       PetscInt i1

428:       i1 = 1

430: !     Get local grid boundaries
431:       call DAGetCorners(da,xs,ys,TAO_NULL_INTEGER,xm,ym,                &
432:      &                TAO_NULL_INTEGER,info)
433:       call DAGetGhostCorners(da,gxs,gys,TAO_NULL_INTEGER,gxm,gym,        &
434:      &                TAO_NULL_INTEGER,info)

436:       do j=ys,ys+ym-1
437:          do i=xs,xs+xm-1
438:             row = (j-gys)*gxm + (i-gxs)

440:             k = 0
441:             if (j .gt. gys) then
442:                v(k) = -1.0d0
443:                col(k) = row-gxm
444:                k = k + 1
445:             endif

447:             if (i .gt. gxs) then
448:                v(k) = -1.0d0
449:                col(k) = row - 1
450:                k = k +1
451:             endif

453:             v(k) = 4.0d0
454:             col(k) = row
455:             k = k + 1

457:             if (i+1 .lt. gxs + gxm) then
458:                v(k) = -1.0d0
459:                col(k) = row + 1
460:                k = k + 1
461:             endif

463:             if (j+1 .lt. gys + gym) then
464:                v(k) = -1.0d0
465:                col(k) = row + gxm
466:                k = k + 1
467:             endif

469:             call MatSetValuesLocal(H,i1,row,k,col,v,INSERT_VALUES,info)
470:          enddo
471:       enddo

473:       
474: !     Assemble matrix
475:       call MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY,info)
476:       call MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY,info)


479: !     Tell the matrix we will never add a new nonzero location to the
480: !     matrix.  If we do it will generate an error.

482:       call MatSetOption(H,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE,info)
483:       call MatSetOption(H,MAT_SYMMETRIC,PETSC_TRUE,info)


486:       call PetscLogFlops(9*xm*ym + 49*xm,info)

488:       info = 0
489:       return 
490:       end
491:       
492:                
493: