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:       integer          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:       integer          iter           ! iteration information
 58:       PetscScalar      ff,gnorm,cnorm,xdiff
 59:       

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

 64:       external FormInitialGuess,FormFunctionGradient,ComputeHessian

 66: !     Initialize TAO, PETSc  contexts
 67:       call PetscInitialize(PETSC_NULL_CHARACTER,info)
 68:       call TaoInitialize(PETSC_NULL_CHARACTER,info)

 70: !     Specify default parameters 
 71:       param = 5.0d0
 72:       mx = 10
 73:       my = 10
 74:       Nx = PETSC_DECIDE
 75:       Ny = PETSC_DECIDE

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

 83:       
 84: !     Set up distributed array and vectors
 85:       call DACreate2d(MPI_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_BOX,      &
 86:      &     mx,my,Nx,Ny,1,1,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,        &
 87:      &     da,info)

 89: !     Create vectors
 90:       call DACreateGlobalVector(da,x,info)
 91:       call DACreateLocalVector(da,localX,info)

 93: !     Create Hessian
 94:       call DAGetMatrix(da,MATAIJ,H,info)
 95:       call MatSetOption(H,MAT_SYMMETRIC,info)

 97: !     The TAO code begins here

 99: !     Create TAO solver
100:       call TaoCreate(MPI_COMM_WORLD,'tao_cg',tao,info)
101:       call TaoApplicationCreate(MPI_COMM_WORLD,torsionapp,info)

103: !     Set routines for function and gradient evaluation

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

112: !     Set initial guess
113:       call FormInitialGuess(x,info)
114:       call TaoAppSetInitialSolutionVec(torsionapp,x,info)

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


120: !     SOLVE THE APPLICATION
121:       call TaoSolveApplication(torsionapp,tao,info)

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

130:       
131: !     Free TAO data structures
132:       call TaoDestroy(tao,info)
133:       call TaoAppDestroy(torsionapp,info);  

135:   
136: !     Free PETSc data structures 
137:       call VecDestroy(x,info)
138:       call VecDestroy(localX,info)
139:       call MatDestroy(H,info)
140:       call DADestroy(da,info)


143: !     Finalize TAO and PETSc
144:       call PetscFinalize(info)
145:       call TaoFinalize(info)

147:       end


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

164: ! mx, my are defined in eptorsion2f.h
165: #include "eptorsion2f.h"

167: !  Input/output variables:
168:       Vec              X
169:       integer          info

171: !  Local variables:
172:       integer          i, j, k, xe, ye
173:       PetscScalar      temp, val, hx, hy
174:       integer          xs, ys, xm, ym, gxm, gym, gxs, gys

176:       hx = 1.0d0/(mx + 1)
177:       hy = 1.0d0/(my + 1)

179: !  Get corner information
180:       call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,               &
181:      &                  PETSC_NULL_INTEGER,info) 
182:       call DAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,              &
183:      &                   gxm,gym,PETSC_NULL_INTEGER,info)



187: !  Compute initial guess over locally owned part of mesh
188:       xe = xs+xm
189:       ye = ys+ym
190:       do j=ys,ye-1
191:          temp = min(j+1,my-j)*hy
192:          do i=xs,xe-1
193:             k   = (j-gys)*gxm + i-gxs
194:             val = min((min(i+1,mx-i))*hx,temp)
195:             call VecSetValuesLocal(X,1,k,val,ADD_VALUES,info)
196:          end do
197:       end do

199:       return
200:       end


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

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

229: !  Input/output variables:
230:       TAO_APPLICATION  taoapp
231:       Vec              X, G
232:       PetscScalar      f
233:       integer          dummy, info

235: !  Declarations for use with local array:


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

246: !  Local variables:
247:       PetscScalar      zero, p5, area, cdiv3
248:       PetscScalar      val, flin, fquad,floc
249:       PetscScalar      v, vb, vl, vr, vt, dvdx
250:       PetscScalar      dvdy, hx, hy
251:       integer          xe, ye, xsm, ysm, xep, yep, i, j, k, ind
252:       integer          xs, ys, xm, ym, gxs, gys, gxm, gym

254:       info  = 0
255:       cdiv3 = param/3.0d0
256:       zero = 0.0d0
257:       p5   = 0.5d0
258:       hx = 1.0d0/(mx + 1)
259:       hy = 1.0d0/(my + 1)
260:       fquad = zero
261:       flin = zero

263: !  Initialize gradient to zero
264:       call VecSet(G,zero,info)

266: !  Scatter ghost points to local vector
267:       call DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX,info)
268:       call DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX,info)


271: !  Get corner information
272:       call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,               &
273:      &                  PETSC_NULL_INTEGER,info)
274:       call DAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,              &
275:      &                   gxm,gym,PETSC_NULL_INTEGER,info)

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


281: !  Set local loop dimensions
282:       xe = xs+xm
283:       ye = ys+ym
284:       if (xs .eq. 0) then
285:          xsm = xs-1
286:       else
287:          xsm = xs
288:       endif
289:       if (ys .eq. 0) then
290:          ysm = ys-1
291:       else
292:          ysm = ys
293:       endif
294:       if (xe .eq. mx) then
295:          xep = xe+1
296:       else
297:          xep = xe
298:       endif
299:       if (ye .eq. my) then
300:          yep = ye+1
301:       else
302:          yep = ye
303:       endif

305: !     Compute local gradient contributions over the lower triangular elements
306:      
307:       do j = ysm, ye-1
308:          do i = xsm, xe-1
309:             k  = (j-gys)*gxm + i-gxs
310:             v  = zero
311:             vr = zero
312:             vt = zero
313:             if (i .ge. 0 .and. j .ge. 0)      v = lx_v(lx_i+k)
314:             if (i .lt. mx-1 .and. j .gt. -1) vr = lx_v(lx_i+k+1)
315:             if (i .gt. -1 .and. j .lt. my-1) vt = lx_v(lx_i+k+gxm)
316:             dvdx = (vr-v)/hx
317:             dvdy = (vt-v)/hy
318:             if (i .ne. -1 .and. j .ne. -1) then
319:                ind = k
320:                val = - dvdx/hx - dvdy/hy - cdiv3
321:                call VecSetValuesLocal(G,1,k,val,ADD_VALUES,info)
322:             endif
323:             if (i .ne. mx-1 .and. j .ne. -1) then
324:                ind = k+1
325:                val =  dvdx/hx - cdiv3
326:                call VecSetValuesLocal(G,1,ind,val,ADD_VALUES,info)
327:             endif
328:             if (i .ne. -1 .and. j .ne. my-1) then
329:               ind = k+gxm
330:               val = dvdy/hy - cdiv3
331:               call VecSetValuesLocal(G,1,ind,val,ADD_VALUES,info)
332:             endif
333:             fquad = fquad + dvdx*dvdx + dvdy*dvdy
334:             flin = flin - cdiv3 * (v+vr+vt)
335:          end do
336:       end do

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

340:       do j = ys, yep-1
341:          do i = xs, xep-1
342:             k  = (j-gys)*gxm + i-gxs
343:             vb = zero
344:             vl = zero
345:             v  = zero
346:             if (i .lt. mx .and. j .gt. 0) vb = lx_v(lx_i+k-gxm)
347:             if (i .gt. 0 .and. j .lt. my) vl = lx_v(lx_i+k-1)
348:             if (i .lt. mx .and. j .lt. my) v = lx_v(lx_i+k)
349:             dvdx = (v-vl)/hx
350:             dvdy = (v-vb)/hy
351:             if (i .ne. mx .and. j .ne. 0) then
352:                ind = k-gxm
353:                val = - dvdy/hy - cdiv3
354:                call VecSetValuesLocal(G,1,ind,val,ADD_VALUES,info)
355:             endif
356:             if (i .ne. 0 .and. j .ne. my) then
357:                ind = k-1
358:                val =  - dvdx/hx - cdiv3
359:                call VecSetValuesLocal(G,1,ind,val,ADD_VALUES,info)
360:             endif
361:             if (i .ne. mx .and. j .ne. my) then
362:                ind = k
363:                val =  dvdx/hx + dvdy/hy - cdiv3
364:                call VecSetValuesLocal(G,1,ind,val,ADD_VALUES,info)
365:             endif
366:             fquad = fquad + dvdx*dvdx + dvdy*dvdy
367:             flin = flin - cdiv3*(vb + vl + v)
368:          end do
369:       end do

371: !  Restore vector
372:       call VecRestoreArray(localX,lx_v,lx_i,info)

374: !  Assemble gradient vector
375:       call VecAssemblyBegin(G,info)
376:       call VecAssemblyEnd(G,info)

378: ! Scale the gradient      
379:       area = p5*hx*hy
380:       floc = area *(p5*fquad+flin)
381:       call VecScale(G,area,info)


384: !  Sum function contributions from all processes
385:       call MPI_Allreduce(floc,f,1,MPI_DOUBLE_PRECISION,MPI_SUM,           &
386:      &                   MPI_COMM_WORLD,info)

388:       call PetscLogFlops((ye-ysm)*(xe-xsm)*20+(xep-xs)*(yep-ys)*16,       &
389:      &                    info)


392:       return
393:       end




398:       subroutine ComputeHessian(taoapp, X, H, Hpre, flag, dummy, info)
399:       implicit none
400: #include "eptorsion2f.h"      
401:       TAO_APPLICATION taoapp
402:       Vec             X
403:       Mat             H,Hpre
404:       MatStructure    flag
405:       integer         dummy,info

407:       
408:       integer i,j,k
409:       integer col(0:4),row
410:       integer xs,xm,gxs,gxm,ys,ym,gys,gym
411:       PetscScalar v(0:4)

413: !     Get local grid boundaries
414:       call DAGetCorners(da,xs,ys,TAO_NULL_INTEGER,xm,ym,                &
415:      &                TAO_NULL_INTEGER,info)
416:       call DAGetGhostCorners(da,gxs,gys,TAO_NULL_INTEGER,gxm,gym,        &
417:      &                TAO_NULL_INTEGER,info)

419:       do j=ys,ys+ym-1
420:          do i=xs,xs+xm-1
421:             row = (j-gys)*gxm + (i-gxs)

423:             k = 0
424:             if (j .gt. gys) then
425:                v(k) = -1.0d0
426:                col(k) = row-gxm
427:                k = k + 1
428:             endif

430:             if (i .gt. gxs) then
431:                v(k) = -1.0d0
432:                col(k) = row - 1
433:                k = k +1
434:             endif

436:             v(k) = 4.0d0
437:             col(k) = row
438:             k = k + 1

440:             if (i+1 .lt. gxs + gxm) then
441:                v(k) = -1.0d0
442:                col(k) = row + 1
443:                k = k + 1
444:             endif

446:             if (j+1 .lt. gys + gym) then
447:                v(k) = -1.0d0
448:                col(k) = row + gxm
449:                k = k + 1
450:             endif

452:             call MatSetValuesLocal(H,1,row,k,col,v,INSERT_VALUES,info)
453:          enddo
454:       enddo

456:       
457: !     Assemble matrix
458:       call MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY,info)
459:       call MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY,info)


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

465:       call MatSetOption(H,MAT_NEW_NONZERO_LOCATION_ERR,info)
466:       call MatSetOption(H,MAT_SYMMETRIC,info)


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

471:       info = 0
472:       return 
473:       end
474:       
475:                
476: