xref: /petsc/src/tao/unconstrained/tutorials/eptorsion2f.F90 (revision 42ce371b2bd7d45eb85bb2bb31075ac1967f9fc8)
1d8606c27SBarry Smith!  Program usage: mpiexec -n <proc> eptorsion2f [all TAO options]
2d8606c27SBarry Smith!
3d8606c27SBarry Smith!  Description:  This example demonstrates use of the TAO package to solve
4d8606c27SBarry Smith!  unconstrained minimization problems in parallel.  This example is based
5d8606c27SBarry Smith!  on the Elastic-Plastic Torsion (dept) problem from the MINPACK-2 test suite.
6d8606c27SBarry Smith!  The command line options are:
7d8606c27SBarry Smith!    -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction
8d8606c27SBarry Smith!    -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction
9d8606c27SBarry Smith!    -par <param>, where <param> = angle of twist per unit length
10d8606c27SBarry Smith!
11d8606c27SBarry Smith!
12d8606c27SBarry Smith! ----------------------------------------------------------------------
13d8606c27SBarry Smith!
14d8606c27SBarry Smith!  Elastic-plastic torsion problem.
15d8606c27SBarry Smith!
16d8606c27SBarry Smith!  The elastic plastic torsion problem arises from the deconverged
17d8606c27SBarry Smith!  of the stress field on an infinitely long cylindrical bar, which is
18d8606c27SBarry Smith!  equivalent to the solution of the following problem:
19d8606c27SBarry Smith!     min{ .5 * integral(||gradient(v(x))||^2 dx) - C * integral(v(x) dx)}
20d8606c27SBarry Smith!  where C is the torsion angle per unit length.
21d8606c27SBarry Smith!
22d8606c27SBarry Smith!  The C version of this code is eptorsion2.c
23d8606c27SBarry Smith!
24d8606c27SBarry Smith! ----------------------------------------------------------------------
25d8606c27SBarry Smith
2677d968b7SBarry Smith      module eptorsion2fmodule
27d8606c27SBarry Smith#include "petsc/finclude/petsctao.h"
28d8606c27SBarry Smith      use petscdmda
29d8606c27SBarry Smith      use petsctao
30d8606c27SBarry Smith      implicit none
31d8606c27SBarry Smith
32d8606c27SBarry Smith      Vec              localX
33d8606c27SBarry Smith      DM               dm
34d8606c27SBarry Smith      PetscReal      param
35d8606c27SBarry Smith      PetscInt         mx, my
36d8606c27SBarry Smith      end module
37d8606c27SBarry Smith
3877d968b7SBarry Smith      use eptorsion2fmodule
39d8606c27SBarry Smith      implicit none
40d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
41d8606c27SBarry Smith!                   Variable declarations
42d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43d8606c27SBarry Smith!
44d8606c27SBarry Smith!  See additional variable declarations in the file eptorsion2f.h
45d8606c27SBarry Smith!
46d8606c27SBarry Smith      PetscErrorCode   ierr           ! used to check for functions returning nonzeros
47d8606c27SBarry Smith      Vec              x              ! solution vector
48d8606c27SBarry Smith      Mat              H              ! hessian matrix
49d8606c27SBarry Smith      PetscInt         Nx, Ny         ! number of processes in x- and y- directions
50d8606c27SBarry Smith      Tao        tao            ! Tao solver context
51d8606c27SBarry Smith      PetscBool        flg
52d8606c27SBarry Smith      PetscInt         i1
53d8606c27SBarry Smith      PetscInt         dummy
54d8606c27SBarry Smith
55d8606c27SBarry Smith!  Note: Any user-defined Fortran routines (such as FormGradient)
56d8606c27SBarry Smith!  MUST be declared as external.
57d8606c27SBarry Smith
58d8606c27SBarry Smith      external FormInitialGuess,FormFunctionGradient,ComputeHessian
59d8606c27SBarry Smith      external Monitor,ConvergenceTest
60d8606c27SBarry Smith
61d8606c27SBarry Smith      i1 = 1
62d8606c27SBarry Smith
63d8606c27SBarry Smith!     Initialize TAO, PETSc  contexts
64d8606c27SBarry Smith      PetscCallA(PetscInitialize(ierr))
65d8606c27SBarry Smith
66d8606c27SBarry Smith!     Specify default parameters
67d8606c27SBarry Smith      param = 5.0
68d8606c27SBarry Smith      mx = 10
69d8606c27SBarry Smith      my = 10
70d8606c27SBarry Smith      Nx = PETSC_DECIDE
71d8606c27SBarry Smith      Ny = PETSC_DECIDE
72d8606c27SBarry Smith
73d8606c27SBarry Smith!     Check for any command line arguments that might override defaults
74d8606c27SBarry Smith      PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-mx',mx,flg,ierr))
75d8606c27SBarry Smith      PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-my',my,flg,ierr))
76d8606c27SBarry Smith      PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-par',param,flg,ierr))
77d8606c27SBarry Smith
78d8606c27SBarry Smith!     Set up distributed array and vectors
79d8606c27SBarry Smith      PetscCallA(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,mx,my,Nx,Ny,i1,i1,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,dm,ierr))
80d8606c27SBarry Smith      PetscCallA(DMSetFromOptions(dm,ierr))
81d8606c27SBarry Smith      PetscCallA(DMSetUp(dm,ierr))
82d8606c27SBarry Smith
83d8606c27SBarry Smith!     Create vectors
84d8606c27SBarry Smith      PetscCallA(DMCreateGlobalVector(dm,x,ierr))
85d8606c27SBarry Smith      PetscCallA(DMCreateLocalVector(dm,localX,ierr))
86d8606c27SBarry Smith
87d8606c27SBarry Smith!     Create Hessian
88d8606c27SBarry Smith      PetscCallA(DMCreateMatrix(dm,H,ierr))
89d8606c27SBarry Smith      PetscCallA(MatSetOption(H,MAT_SYMMETRIC,PETSC_TRUE,ierr))
90d8606c27SBarry Smith
91d8606c27SBarry Smith!     The TAO code begins here
92d8606c27SBarry Smith
93d8606c27SBarry Smith!     Create TAO solver
94d8606c27SBarry Smith      PetscCallA(TaoCreate(PETSC_COMM_WORLD,tao,ierr))
95d8606c27SBarry Smith      PetscCallA(TaoSetType(tao,TAOCG,ierr))
96d8606c27SBarry Smith
97d8606c27SBarry Smith!     Set routines for function and gradient evaluation
98d8606c27SBarry Smith
99d8606c27SBarry Smith      PetscCallA(TaoSetObjectiveAndGradient(tao,PETSC_NULL_VEC,FormFunctionGradient,0,ierr))
100d8606c27SBarry Smith      PetscCallA(TaoSetHessian(tao,H,H,ComputeHessian,0,ierr))
101d8606c27SBarry Smith
102d8606c27SBarry Smith!     Set initial guess
103d8606c27SBarry Smith      PetscCallA(FormInitialGuess(x,ierr))
104d8606c27SBarry Smith      PetscCallA(TaoSetSolution(tao,x,ierr))
105d8606c27SBarry Smith
106d8606c27SBarry Smith      PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-testmonitor',flg,ierr))
107d8606c27SBarry Smith      if (flg) then
108d8606c27SBarry Smith         PetscCallA(TaoSetMonitor(tao,Monitor,dummy,PETSC_NULL_FUNCTION,ierr))
109d8606c27SBarry Smith      endif
110d8606c27SBarry Smith
111d8606c27SBarry Smith      PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-testconvergence',flg, ierr))
112d8606c27SBarry Smith      if (flg) then
113d8606c27SBarry Smith         PetscCallA(TaoSetConvergenceTest(tao,ConvergenceTest,dummy,ierr))
114d8606c27SBarry Smith      endif
115d8606c27SBarry Smith
116d8606c27SBarry Smith!     Check for any TAO command line options
117d8606c27SBarry Smith      PetscCallA(TaoSetFromOptions(tao,ierr))
118d8606c27SBarry Smith
119d8606c27SBarry Smith!     SOLVE THE APPLICATION
120d8606c27SBarry Smith      PetscCallA(TaoSolve(tao,ierr))
121d8606c27SBarry Smith
122d8606c27SBarry Smith!     Free TAO data structures
123d8606c27SBarry Smith      PetscCallA(TaoDestroy(tao,ierr))
124d8606c27SBarry Smith
125d8606c27SBarry Smith!     Free PETSc data structures
126d8606c27SBarry Smith      PetscCallA(VecDestroy(x,ierr))
127d8606c27SBarry Smith      PetscCallA(VecDestroy(localX,ierr))
128d8606c27SBarry Smith      PetscCallA(MatDestroy(H,ierr))
129d8606c27SBarry Smith      PetscCallA(DMDestroy(dm,ierr))
130d8606c27SBarry Smith
131d8606c27SBarry Smith!     Finalize TAO and PETSc
132d8606c27SBarry Smith      PetscCallA(PetscFinalize(ierr))
133d8606c27SBarry Smith      end
134d8606c27SBarry Smith
135d8606c27SBarry Smith! ---------------------------------------------------------------------
136d8606c27SBarry Smith!
137d8606c27SBarry Smith!   FormInitialGuess - Computes an initial approximation to the solution.
138d8606c27SBarry Smith!
139d8606c27SBarry Smith!   Input Parameters:
140d8606c27SBarry Smith!   X    - vector
141d8606c27SBarry Smith!
142d8606c27SBarry Smith!   Output Parameters:
143d8606c27SBarry Smith!   X    - vector
144d8606c27SBarry Smith!   ierr - error code
145d8606c27SBarry Smith!
146d8606c27SBarry Smith      subroutine FormInitialGuess(X,ierr)
14777d968b7SBarry Smith      use eptorsion2fmodule
148d8606c27SBarry Smith      implicit none
149d8606c27SBarry Smith
150d8606c27SBarry Smith!  Input/output variables:
151d8606c27SBarry Smith      Vec              X
152d8606c27SBarry Smith      PetscErrorCode   ierr
153d8606c27SBarry Smith
154d8606c27SBarry Smith!  Local variables:
155d8606c27SBarry Smith      PetscInt         i, j, k, xe, ye
156d8606c27SBarry Smith      PetscReal      temp, val, hx, hy
157d8606c27SBarry Smith      PetscInt         xs, ys, xm, ym
158d8606c27SBarry Smith      PetscInt         gxm, gym, gxs, gys
159d8606c27SBarry Smith      PetscInt         i1
160d8606c27SBarry Smith
161d8606c27SBarry Smith      i1 = 1
162d8606c27SBarry Smith      hx = 1.0/real(mx + 1)
163d8606c27SBarry Smith      hy = 1.0/real(my + 1)
164d8606c27SBarry Smith
165d8606c27SBarry Smith!  Get corner information
166d8606c27SBarry Smith      PetscCall(DMDAGetCorners(dm,xs,ys,PETSC_NULL_INTEGER,xm,ym,PETSC_NULL_INTEGER,ierr))
167d8606c27SBarry Smith      PetscCall(DMDAGetGhostCorners(dm,gxs,gys,PETSC_NULL_INTEGER,gxm,gym,PETSC_NULL_INTEGER,ierr))
168d8606c27SBarry Smith
169d8606c27SBarry Smith!  Compute initial guess over locally owned part of mesh
170d8606c27SBarry Smith      xe = xs+xm
171d8606c27SBarry Smith      ye = ys+ym
172d8606c27SBarry Smith      do j=ys,ye-1
173d8606c27SBarry Smith         temp = min(j+1,my-j)*hy
174d8606c27SBarry Smith         do i=xs,xe-1
175d8606c27SBarry Smith            k   = (j-gys)*gxm + i-gxs
176d8606c27SBarry Smith            val = min((min(i+1,mx-i))*hx,temp)
177d8606c27SBarry Smith            PetscCall(VecSetValuesLocal(X,i1,k,val,ADD_VALUES,ierr))
178d8606c27SBarry Smith         end do
179d8606c27SBarry Smith      end do
180d8606c27SBarry Smith      PetscCall(VecAssemblyBegin(X,ierr))
181d8606c27SBarry Smith      PetscCall(VecAssemblyEnd(X,ierr))
182d8606c27SBarry Smith      return
183d8606c27SBarry Smith      end
184d8606c27SBarry Smith
185d8606c27SBarry Smith! ---------------------------------------------------------------------
186d8606c27SBarry Smith!
187d8606c27SBarry Smith!  FormFunctionGradient - Evaluates gradient G(X).
188d8606c27SBarry Smith!
189d8606c27SBarry Smith!  Input Parameters:
190d8606c27SBarry Smith!  tao   - the Tao context
191d8606c27SBarry Smith!  X     - input vector
192d8606c27SBarry Smith!  dummy - optional user-defined context (not used here)
193d8606c27SBarry Smith!
194d8606c27SBarry Smith!  Output Parameters:
195d8606c27SBarry Smith!  f     - the function value at X
196d8606c27SBarry Smith!  G     - vector containing the newly evaluated gradient
197d8606c27SBarry Smith!  ierr  - error code
198d8606c27SBarry Smith!
199d8606c27SBarry Smith!  Notes:
200d8606c27SBarry Smith!  This routine serves as a wrapper for the lower-level routine
201d8606c27SBarry Smith!  "ApplicationGradient", where the actual computations are
202d8606c27SBarry Smith!  done using the standard Fortran style of treating the local
203d8606c27SBarry Smith!  input vector data as an array over the local mesh.
204d8606c27SBarry Smith!
205d8606c27SBarry Smith      subroutine FormFunctionGradient(tao,X,f,G,dummy,ierr)
20677d968b7SBarry Smith      use eptorsion2fmodule
207d8606c27SBarry Smith      implicit none
208d8606c27SBarry Smith
209d8606c27SBarry Smith!  Input/output variables:
210d8606c27SBarry Smith      Tao        tao
211d8606c27SBarry Smith      Vec              X, G
212d8606c27SBarry Smith      PetscReal      f
213d8606c27SBarry Smith      PetscErrorCode   ierr
214d8606c27SBarry Smith      PetscInt         dummy
215d8606c27SBarry Smith
216d8606c27SBarry Smith!  Declarations for use with local array:
217d8606c27SBarry Smith
218*42ce371bSBarry Smith      PetscReal, pointer :: lx_v(:)
219d8606c27SBarry Smith
220d8606c27SBarry Smith!  Local variables:
221d8606c27SBarry Smith      PetscReal      zero, p5, area, cdiv3
222d8606c27SBarry Smith      PetscReal      val, flin, fquad,floc
223d8606c27SBarry Smith      PetscReal      v, vb, vl, vr, vt, dvdx
224d8606c27SBarry Smith      PetscReal      dvdy, hx, hy
225d8606c27SBarry Smith      PetscInt         xe, ye, xsm, ysm
226d8606c27SBarry Smith      PetscInt         xep, yep, i, j, k, ind
227d8606c27SBarry Smith      PetscInt         xs, ys, xm, ym
228d8606c27SBarry Smith      PetscInt         gxs, gys, gxm, gym
229d8606c27SBarry Smith      PetscInt         i1
230d8606c27SBarry Smith
231d8606c27SBarry Smith      i1 = 1
232d8606c27SBarry Smith      ierr  = 0
233d8606c27SBarry Smith      cdiv3 = param/3.0
234d8606c27SBarry Smith      zero = 0.0
235d8606c27SBarry Smith      p5   = 0.5
236d8606c27SBarry Smith      hx = 1.0/real(mx + 1)
237d8606c27SBarry Smith      hy = 1.0/real(my + 1)
238d8606c27SBarry Smith      fquad = zero
239d8606c27SBarry Smith      flin = zero
240d8606c27SBarry Smith
241d8606c27SBarry Smith!  Initialize gradient to zero
242d8606c27SBarry Smith      PetscCall(VecSet(G,zero,ierr))
243d8606c27SBarry Smith
244d8606c27SBarry Smith!  Scatter ghost points to local vector
245d8606c27SBarry Smith      PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,localX,ierr))
246d8606c27SBarry Smith      PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,localX,ierr))
247d8606c27SBarry Smith
248d8606c27SBarry Smith!  Get corner information
249d8606c27SBarry Smith      PetscCall(DMDAGetCorners(dm,xs,ys,PETSC_NULL_INTEGER,xm,ym,PETSC_NULL_INTEGER,ierr))
250d8606c27SBarry Smith      PetscCall(DMDAGetGhostCorners(dm,gxs,gys,PETSC_NULL_INTEGER,gxm,gym,PETSC_NULL_INTEGER,ierr))
251d8606c27SBarry Smith
252d8606c27SBarry Smith!  Get pointer to vector data.
253*42ce371bSBarry Smith      PetscCall(VecGetArrayReadF90(localX,lx_v,ierr))
254d8606c27SBarry Smith
255d8606c27SBarry Smith!  Set local loop dimensions
256d8606c27SBarry Smith      xe = xs+xm
257d8606c27SBarry Smith      ye = ys+ym
258d8606c27SBarry Smith      if (xs .eq. 0) then
259d8606c27SBarry Smith         xsm = xs-1
260d8606c27SBarry Smith      else
261d8606c27SBarry Smith         xsm = xs
262d8606c27SBarry Smith      endif
263d8606c27SBarry Smith      if (ys .eq. 0) then
264d8606c27SBarry Smith         ysm = ys-1
265d8606c27SBarry Smith      else
266d8606c27SBarry Smith         ysm = ys
267d8606c27SBarry Smith      endif
268d8606c27SBarry Smith      if (xe .eq. mx) then
269d8606c27SBarry Smith         xep = xe+1
270d8606c27SBarry Smith      else
271d8606c27SBarry Smith         xep = xe
272d8606c27SBarry Smith      endif
273d8606c27SBarry Smith      if (ye .eq. my) then
274d8606c27SBarry Smith         yep = ye+1
275d8606c27SBarry Smith      else
276d8606c27SBarry Smith         yep = ye
277d8606c27SBarry Smith      endif
278d8606c27SBarry Smith
279d8606c27SBarry Smith!     Compute local gradient contributions over the lower triangular elements
280d8606c27SBarry Smith
281d8606c27SBarry Smith      do j = ysm, ye-1
282d8606c27SBarry Smith         do i = xsm, xe-1
283d8606c27SBarry Smith            k  = (j-gys)*gxm + i-gxs
284d8606c27SBarry Smith            v  = zero
285d8606c27SBarry Smith            vr = zero
286d8606c27SBarry Smith            vt = zero
287*42ce371bSBarry Smith            if (i .ge. 0 .and. j .ge. 0)      v = lx_v(k+1)
288*42ce371bSBarry Smith            if (i .lt. mx-1 .and. j .gt. -1) vr = lx_v(k+2)
289*42ce371bSBarry Smith            if (i .gt. -1 .and. j .lt. my-1) vt = lx_v(k+1+gxm)
290d8606c27SBarry Smith            dvdx = (vr-v)/hx
291d8606c27SBarry Smith            dvdy = (vt-v)/hy
292d8606c27SBarry Smith            if (i .ne. -1 .and. j .ne. -1) then
293d8606c27SBarry Smith               ind = k
294d8606c27SBarry Smith               val = - dvdx/hx - dvdy/hy - cdiv3
295d8606c27SBarry Smith               PetscCall(VecSetValuesLocal(G,i1,k,val,ADD_VALUES,ierr))
296d8606c27SBarry Smith            endif
297d8606c27SBarry Smith            if (i .ne. mx-1 .and. j .ne. -1) then
298d8606c27SBarry Smith               ind = k+1
299d8606c27SBarry Smith               val =  dvdx/hx - cdiv3
300d8606c27SBarry Smith               PetscCall(VecSetValuesLocal(G,i1,ind,val,ADD_VALUES,ierr))
301d8606c27SBarry Smith            endif
302d8606c27SBarry Smith            if (i .ne. -1 .and. j .ne. my-1) then
303d8606c27SBarry Smith              ind = k+gxm
304d8606c27SBarry Smith              val = dvdy/hy - cdiv3
305d8606c27SBarry Smith              PetscCall(VecSetValuesLocal(G,i1,ind,val,ADD_VALUES,ierr))
306d8606c27SBarry Smith            endif
307d8606c27SBarry Smith            fquad = fquad + dvdx*dvdx + dvdy*dvdy
308d8606c27SBarry Smith            flin = flin - cdiv3 * (v+vr+vt)
309d8606c27SBarry Smith         end do
310d8606c27SBarry Smith      end do
311d8606c27SBarry Smith
312d8606c27SBarry Smith!     Compute local gradient contributions over the upper triangular elements
313d8606c27SBarry Smith
314d8606c27SBarry Smith      do j = ys, yep-1
315d8606c27SBarry Smith         do i = xs, xep-1
316d8606c27SBarry Smith            k  = (j-gys)*gxm + i-gxs
317d8606c27SBarry Smith            vb = zero
318d8606c27SBarry Smith            vl = zero
319d8606c27SBarry Smith            v  = zero
320*42ce371bSBarry Smith            if (i .lt. mx .and. j .gt. 0) vb = lx_v(k+1-gxm)
321*42ce371bSBarry Smith            if (i .gt. 0 .and. j .lt. my) vl = lx_v(k)
322*42ce371bSBarry Smith            if (i .lt. mx .and. j .lt. my) v = lx_v(1+k)
323d8606c27SBarry Smith            dvdx = (v-vl)/hx
324d8606c27SBarry Smith            dvdy = (v-vb)/hy
325d8606c27SBarry Smith            if (i .ne. mx .and. j .ne. 0) then
326d8606c27SBarry Smith               ind = k-gxm
327d8606c27SBarry Smith               val = - dvdy/hy - cdiv3
328d8606c27SBarry Smith               PetscCall(VecSetValuesLocal(G,i1,ind,val,ADD_VALUES,ierr))
329d8606c27SBarry Smith            endif
330d8606c27SBarry Smith            if (i .ne. 0 .and. j .ne. my) then
331d8606c27SBarry Smith               ind = k-1
332d8606c27SBarry Smith               val =  - dvdx/hx - cdiv3
333d8606c27SBarry Smith               PetscCall(VecSetValuesLocal(G,i1,ind,val,ADD_VALUES,ierr))
334d8606c27SBarry Smith            endif
335d8606c27SBarry Smith            if (i .ne. mx .and. j .ne. my) then
336d8606c27SBarry Smith               ind = k
337d8606c27SBarry Smith               val =  dvdx/hx + dvdy/hy - cdiv3
338d8606c27SBarry Smith               PetscCall(VecSetValuesLocal(G,i1,ind,val,ADD_VALUES,ierr))
339d8606c27SBarry Smith            endif
340d8606c27SBarry Smith            fquad = fquad + dvdx*dvdx + dvdy*dvdy
341d8606c27SBarry Smith            flin = flin - cdiv3*(vb + vl + v)
342d8606c27SBarry Smith         end do
343d8606c27SBarry Smith      end do
344d8606c27SBarry Smith
345d8606c27SBarry Smith!  Restore vector
346*42ce371bSBarry Smith      PetscCall(VecRestoreArrayReadF90(localX,lx_v,ierr))
347d8606c27SBarry Smith
348d8606c27SBarry Smith!  Assemble gradient vector
349d8606c27SBarry Smith      PetscCall(VecAssemblyBegin(G,ierr))
350d8606c27SBarry Smith      PetscCall(VecAssemblyEnd(G,ierr))
351d8606c27SBarry Smith
352d8606c27SBarry Smith! Scale the gradient
353d8606c27SBarry Smith      area = p5*hx*hy
354d8606c27SBarry Smith      floc = area *(p5*fquad+flin)
355d8606c27SBarry Smith      PetscCall(VecScale(G,area,ierr))
356d8606c27SBarry Smith
357d8606c27SBarry Smith!  Sum function contributions from all processes
358d8606c27SBarry Smith      PetscCallMPI(MPI_Allreduce(floc,f,1,MPIU_SCALAR,MPIU_SUM,PETSC_COMM_WORLD,ierr))
359d8606c27SBarry Smith      PetscCall(PetscLogFlops(20.0d0*(ye-ysm)*(xe-xsm)+16.0d0*(xep-xs)*(yep-ys),ierr))
360d8606c27SBarry Smith      return
361d8606c27SBarry Smith      end
362d8606c27SBarry Smith
363d8606c27SBarry Smith      subroutine ComputeHessian(tao, X, H, Hpre, dummy, ierr)
36477d968b7SBarry Smith      use eptorsion2fmodule
365d8606c27SBarry Smith      implicit none
366d8606c27SBarry Smith
367d8606c27SBarry Smith      Tao       tao
368d8606c27SBarry Smith      Vec             X
369d8606c27SBarry Smith      Mat             H,Hpre
370d8606c27SBarry Smith      PetscErrorCode  ierr
371d8606c27SBarry Smith      PetscInt        dummy
372d8606c27SBarry Smith
373d8606c27SBarry Smith      PetscInt i,j,k
374d8606c27SBarry Smith      PetscInt col(0:4),row
375d8606c27SBarry Smith      PetscInt xs,xm,gxs,gxm
376d8606c27SBarry Smith      PetscInt ys,ym,gys,gym
377d8606c27SBarry Smith      PetscReal v(0:4)
378d8606c27SBarry Smith      PetscInt i1
379d8606c27SBarry Smith
380d8606c27SBarry Smith      i1 = 1
381d8606c27SBarry Smith
382d8606c27SBarry Smith!     Get local grid boundaries
383d8606c27SBarry Smith      PetscCall(DMDAGetCorners(dm,xs,ys,PETSC_NULL_INTEGER,xm,ym,PETSC_NULL_INTEGER,ierr))
384d8606c27SBarry Smith      PetscCall(DMDAGetGhostCorners(dm,gxs,gys,PETSC_NULL_INTEGER,gxm,gym,PETSC_NULL_INTEGER,ierr))
385d8606c27SBarry Smith
386d8606c27SBarry Smith      do j=ys,ys+ym-1
387d8606c27SBarry Smith         do i=xs,xs+xm-1
388d8606c27SBarry Smith            row = (j-gys)*gxm + (i-gxs)
389d8606c27SBarry Smith
390d8606c27SBarry Smith            k = 0
391d8606c27SBarry Smith            if (j .gt. gys) then
392d8606c27SBarry Smith               v(k) = -1.0
393d8606c27SBarry Smith               col(k) = row-gxm
394d8606c27SBarry Smith               k = k + 1
395d8606c27SBarry Smith            endif
396d8606c27SBarry Smith
397d8606c27SBarry Smith            if (i .gt. gxs) then
398d8606c27SBarry Smith               v(k) = -1.0
399d8606c27SBarry Smith               col(k) = row - 1
400d8606c27SBarry Smith               k = k +1
401d8606c27SBarry Smith            endif
402d8606c27SBarry Smith
403d8606c27SBarry Smith            v(k) = 4.0
404d8606c27SBarry Smith            col(k) = row
405d8606c27SBarry Smith            k = k + 1
406d8606c27SBarry Smith
407d8606c27SBarry Smith            if (i+1 .lt. gxs + gxm) then
408d8606c27SBarry Smith               v(k) = -1.0
409d8606c27SBarry Smith               col(k) = row + 1
410d8606c27SBarry Smith               k = k + 1
411d8606c27SBarry Smith            endif
412d8606c27SBarry Smith
413d8606c27SBarry Smith            if (j+1 .lt. gys + gym) then
414d8606c27SBarry Smith               v(k) = -1.0
415d8606c27SBarry Smith               col(k) = row + gxm
416d8606c27SBarry Smith               k = k + 1
417d8606c27SBarry Smith            endif
418d8606c27SBarry Smith
419d8606c27SBarry Smith            PetscCall(MatSetValuesLocal(H,i1,row,k,col,v,INSERT_VALUES,ierr))
420d8606c27SBarry Smith         enddo
421d8606c27SBarry Smith      enddo
422d8606c27SBarry Smith
423d8606c27SBarry Smith!     Assemble matrix
424d8606c27SBarry Smith      PetscCall(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY,ierr))
425d8606c27SBarry Smith      PetscCall(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY,ierr))
426d8606c27SBarry Smith
427d8606c27SBarry Smith!     Tell the matrix we will never add a new nonzero location to the
428d8606c27SBarry Smith!     matrix.  If we do it will generate an error.
429d8606c27SBarry Smith
430d8606c27SBarry Smith      PetscCall(MatSetOption(H,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE,ierr))
431d8606c27SBarry Smith      PetscCall(MatSetOption(H,MAT_SYMMETRIC,PETSC_TRUE,ierr))
432d8606c27SBarry Smith
433d8606c27SBarry Smith      PetscCall(PetscLogFlops(9.0d0*xm*ym + 49.0d0*xm,ierr))
434d8606c27SBarry Smith
435d8606c27SBarry Smith      ierr = 0
436d8606c27SBarry Smith      return
437d8606c27SBarry Smith      end
438d8606c27SBarry Smith
439d8606c27SBarry Smith      subroutine Monitor(tao, dummy, ierr)
44077d968b7SBarry Smith      use eptorsion2fmodule
441d8606c27SBarry Smith      implicit none
442d8606c27SBarry Smith
443d8606c27SBarry Smith      Tao tao
444d8606c27SBarry Smith      PetscInt dummy
445d8606c27SBarry Smith      PetscErrorCode ierr
446d8606c27SBarry Smith
447d8606c27SBarry Smith      PetscInt its
448d8606c27SBarry Smith      PetscReal f,gnorm,cnorm,xdiff
449d8606c27SBarry Smith      TaoConvergedReason reason
450d8606c27SBarry Smith
451d8606c27SBarry Smith      PetscCall(TaoGetSolutionStatus(tao,its,f,gnorm,cnorm,xdiff,reason,ierr))
452d8606c27SBarry Smith      if (mod(its,5) .ne. 0) then
453d8606c27SBarry Smith         PetscCall(PetscPrintf(PETSC_COMM_WORLD,'iteration multiple of 5\n',ierr))
454d8606c27SBarry Smith      endif
455d8606c27SBarry Smith
456d8606c27SBarry Smith      ierr = 0
457d8606c27SBarry Smith
458d8606c27SBarry Smith      return
459d8606c27SBarry Smith      end
460d8606c27SBarry Smith
461d8606c27SBarry Smith      subroutine ConvergenceTest(tao, dummy, ierr)
46277d968b7SBarry Smith      use eptorsion2fmodule
463d8606c27SBarry Smith      implicit none
464d8606c27SBarry Smith
465d8606c27SBarry Smith      Tao tao
466d8606c27SBarry Smith      PetscInt dummy
467d8606c27SBarry Smith      PetscErrorCode ierr
468d8606c27SBarry Smith
469d8606c27SBarry Smith      PetscInt its
470d8606c27SBarry Smith      PetscReal f,gnorm,cnorm,xdiff
471d8606c27SBarry Smith      TaoConvergedReason reason
472d8606c27SBarry Smith
473d8606c27SBarry Smith      PetscCall(TaoGetSolutionStatus(tao,its,f,gnorm,cnorm,xdiff,reason,ierr))
474d8606c27SBarry Smith      if (its .eq. 7) then
475d8606c27SBarry Smith       PetscCall(TaoSetConvergedReason(tao,TAO_DIVERGED_MAXITS,ierr))
476d8606c27SBarry Smith      endif
477d8606c27SBarry Smith
478d8606c27SBarry Smith      ierr = 0
479d8606c27SBarry Smith
480d8606c27SBarry Smith      return
481d8606c27SBarry Smith      end
482d8606c27SBarry Smith
483d8606c27SBarry Smith!/*TEST
484d8606c27SBarry Smith!
485d8606c27SBarry Smith!   build:
486d8606c27SBarry Smith!      requires: !complex
487d8606c27SBarry Smith!
488d8606c27SBarry Smith!   test:
489d8606c27SBarry Smith!      args: -tao_smonitor -tao_type nls -tao_gttol 1.e-2
490d8606c27SBarry Smith!
491d8606c27SBarry Smith!   test:
492d8606c27SBarry Smith!      suffix: 2
493d8606c27SBarry Smith!      nsize: 2
494d8606c27SBarry Smith!      args: -tao_smonitor -tao_type lmvm -tao_gttol 1.e-2
495d8606c27SBarry Smith!
496d8606c27SBarry Smith!   test:
497d8606c27SBarry Smith!      suffix: 3
498d8606c27SBarry Smith!      args: -testmonitor -tao_type lmvm -tao_gttol 1.e-2
499d8606c27SBarry Smith!TEST*/
500