xref: /petsc/src/tao/unconstrained/tutorials/eptorsion2f.F90 (revision 9b88ac225e01f016352a5f4cd90e158abe5f5675)
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! ----------------------------------------------------------------------
25ce78bad3SBarry Smith#include "petsc/finclude/petscdmda.h"
26d8606c27SBarry Smith#include "petsc/finclude/petsctao.h"
27c5e229c2SMartin Diehlmodule eptorsion2fmodule
28d8606c27SBarry Smith  use petscdmda
29d8606c27SBarry Smith  use petsctao
30d8606c27SBarry Smith  implicit none
31d8606c27SBarry Smith
32392a88c0SBlaise Bourdin  type(tVec) localX
33392a88c0SBlaise Bourdin  type(tDM) dm
34d8606c27SBarry Smith  PetscReal param
35d8606c27SBarry Smith  PetscInt mx, my
36d8606c27SBarry Smith
37*e7a95102SMartin Diehlcontains
38d8606c27SBarry Smith! ---------------------------------------------------------------------
39d8606c27SBarry Smith!
40d8606c27SBarry Smith!   FormInitialGuess - Computes an initial approximation to the solution.
41d8606c27SBarry Smith!
42d8606c27SBarry Smith!   Input Parameters:
43d8606c27SBarry Smith!   X    - vector
44d8606c27SBarry Smith!
45d8606c27SBarry Smith!   Output Parameters:
46d8606c27SBarry Smith!   X    - vector
47d8606c27SBarry Smith!   ierr - error code
48d8606c27SBarry Smith!
49d8606c27SBarry Smith  subroutine FormInitialGuess(X, ierr)
50d8606c27SBarry Smith!  Input/output variables:
51d8606c27SBarry Smith    Vec X
52d8606c27SBarry Smith    PetscErrorCode ierr
53d8606c27SBarry Smith
54d8606c27SBarry Smith!  Local variables:
55d8606c27SBarry Smith    PetscInt i, j, k, xe, ye
56d8606c27SBarry Smith    PetscReal temp, val, hx, hy
57d8606c27SBarry Smith    PetscInt xs, ys, xm, ym
58d8606c27SBarry Smith    PetscInt gxm, gym, gxs, gys
59d8606c27SBarry Smith    PetscInt i1
60d8606c27SBarry Smith
61d8606c27SBarry Smith    i1 = 1
62d8606c27SBarry Smith    hx = 1.0/real(mx + 1)
63d8606c27SBarry Smith    hy = 1.0/real(my + 1)
64d8606c27SBarry Smith
65d8606c27SBarry Smith!  Get corner information
66d8606c27SBarry Smith    PetscCall(DMDAGetCorners(dm, xs, ys, PETSC_NULL_INTEGER, xm, ym, PETSC_NULL_INTEGER, ierr))
67d8606c27SBarry Smith    PetscCall(DMDAGetGhostCorners(dm, gxs, gys, PETSC_NULL_INTEGER, gxm, gym, PETSC_NULL_INTEGER, ierr))
68d8606c27SBarry Smith
69d8606c27SBarry Smith!  Compute initial guess over locally owned part of mesh
70d8606c27SBarry Smith    xe = xs + xm
71d8606c27SBarry Smith    ye = ys + ym
72d8606c27SBarry Smith    do j = ys, ye - 1
73d8606c27SBarry Smith      temp = min(j + 1, my - j)*hy
74d8606c27SBarry Smith      do i = xs, xe - 1
75d8606c27SBarry Smith        k = (j - gys)*gxm + i - gxs
76d8606c27SBarry Smith        val = min((min(i + 1, mx - i))*hx, temp)
775d83a8b1SBarry Smith        PetscCall(VecSetValuesLocal(X, i1, [k], [val], ADD_VALUES, ierr))
78d8606c27SBarry Smith      end do
79d8606c27SBarry Smith    end do
80d8606c27SBarry Smith    PetscCall(VecAssemblyBegin(X, ierr))
81d8606c27SBarry Smith    PetscCall(VecAssemblyEnd(X, ierr))
82d8606c27SBarry Smith  end
83d8606c27SBarry Smith
84d8606c27SBarry Smith! ---------------------------------------------------------------------
85d8606c27SBarry Smith!
86d8606c27SBarry Smith!  FormFunctionGradient - Evaluates gradient G(X).
87d8606c27SBarry Smith!
88d8606c27SBarry Smith!  Input Parameters:
89d8606c27SBarry Smith!  tao   - the Tao context
90d8606c27SBarry Smith!  X     - input vector
91d8606c27SBarry Smith!  dummy - optional user-defined context (not used here)
92d8606c27SBarry Smith!
93d8606c27SBarry Smith!  Output Parameters:
94d8606c27SBarry Smith!  f     - the function value at X
95d8606c27SBarry Smith!  G     - vector containing the newly evaluated gradient
96d8606c27SBarry Smith!  ierr  - error code
97d8606c27SBarry Smith!
98d8606c27SBarry Smith!  Notes:
99d8606c27SBarry Smith!  This routine serves as a wrapper for the lower-level routine
100d8606c27SBarry Smith!  "ApplicationGradient", where the actual computations are
101d8606c27SBarry Smith!  done using the standard Fortran style of treating the local
102d8606c27SBarry Smith!  input vector data as an array over the local mesh.
103d8606c27SBarry Smith!
104ce78bad3SBarry Smith  subroutine FormFunctionGradient(ta, X, f, G, dummy, ierr)
105d8606c27SBarry Smith!  Input/output variables:
106ce78bad3SBarry Smith    type(tTao) ta
107392a88c0SBlaise Bourdin    type(tVec) X, G
108d8606c27SBarry Smith    PetscReal f
109d8606c27SBarry Smith    PetscErrorCode ierr
110d8606c27SBarry Smith    PetscInt dummy
111d8606c27SBarry Smith
112d8606c27SBarry Smith!  Declarations for use with local array:
113d8606c27SBarry Smith
11442ce371bSBarry Smith    PetscReal, pointer :: lx_v(:)
115d8606c27SBarry Smith
116d8606c27SBarry Smith!  Local variables:
117d8606c27SBarry Smith    PetscReal zero, p5, area, cdiv3
118d8606c27SBarry Smith    PetscReal val, flin, fquad, floc
119d8606c27SBarry Smith    PetscReal v, vb, vl, vr, vt, dvdx
120d8606c27SBarry Smith    PetscReal dvdy, hx, hy
121d8606c27SBarry Smith    PetscInt xe, ye, xsm, ysm
122d8606c27SBarry Smith    PetscInt xep, yep, i, j, k, ind
123d8606c27SBarry Smith    PetscInt xs, ys, xm, ym
124d8606c27SBarry Smith    PetscInt gxs, gys, gxm, gym
125d8606c27SBarry Smith    PetscInt i1
126d8606c27SBarry Smith
127d8606c27SBarry Smith    i1 = 1
128d8606c27SBarry Smith    ierr = 0
129d8606c27SBarry Smith    cdiv3 = param/3.0
130d8606c27SBarry Smith    zero = 0.0
131d8606c27SBarry Smith    p5 = 0.5
132d8606c27SBarry Smith    hx = 1.0/real(mx + 1)
133d8606c27SBarry Smith    hy = 1.0/real(my + 1)
134d8606c27SBarry Smith    fquad = zero
135d8606c27SBarry Smith    flin = zero
136d8606c27SBarry Smith
137d8606c27SBarry Smith!  Initialize gradient to zero
138d8606c27SBarry Smith    PetscCall(VecSet(G, zero, ierr))
139d8606c27SBarry Smith
140d8606c27SBarry Smith!  Scatter ghost points to local vector
141d8606c27SBarry Smith    PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX, ierr))
142d8606c27SBarry Smith    PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX, ierr))
143d8606c27SBarry Smith
144d8606c27SBarry Smith!  Get corner information
145d8606c27SBarry Smith    PetscCall(DMDAGetCorners(dm, xs, ys, PETSC_NULL_INTEGER, xm, ym, PETSC_NULL_INTEGER, ierr))
146d8606c27SBarry Smith    PetscCall(DMDAGetGhostCorners(dm, gxs, gys, PETSC_NULL_INTEGER, gxm, gym, PETSC_NULL_INTEGER, ierr))
147d8606c27SBarry Smith
148d8606c27SBarry Smith!  Get pointer to vector data.
149ce78bad3SBarry Smith    PetscCall(VecGetArrayRead(localX, lx_v, ierr))
150d8606c27SBarry Smith
151d8606c27SBarry Smith!  Set local loop dimensions
152d8606c27SBarry Smith    xe = xs + xm
153d8606c27SBarry Smith    ye = ys + ym
1544820e4eaSBarry Smith    if (xs == 0) then
155d8606c27SBarry Smith      xsm = xs - 1
156d8606c27SBarry Smith    else
157d8606c27SBarry Smith      xsm = xs
158d8606c27SBarry Smith    end if
1594820e4eaSBarry Smith    if (ys == 0) then
160d8606c27SBarry Smith      ysm = ys - 1
161d8606c27SBarry Smith    else
162d8606c27SBarry Smith      ysm = ys
163d8606c27SBarry Smith    end if
1644820e4eaSBarry Smith    if (xe == mx) then
165d8606c27SBarry Smith      xep = xe + 1
166d8606c27SBarry Smith    else
167d8606c27SBarry Smith      xep = xe
168d8606c27SBarry Smith    end if
1694820e4eaSBarry Smith    if (ye == my) then
170d8606c27SBarry Smith      yep = ye + 1
171d8606c27SBarry Smith    else
172d8606c27SBarry Smith      yep = ye
173d8606c27SBarry Smith    end if
174d8606c27SBarry Smith
175d8606c27SBarry Smith!     Compute local gradient contributions over the lower triangular elements
176d8606c27SBarry Smith
177d8606c27SBarry Smith    do j = ysm, ye - 1
178d8606c27SBarry Smith      do i = xsm, xe - 1
179d8606c27SBarry Smith        k = (j - gys)*gxm + i - gxs
180d8606c27SBarry Smith        v = zero
181d8606c27SBarry Smith        vr = zero
182d8606c27SBarry Smith        vt = zero
1834820e4eaSBarry Smith        if (i >= 0 .and. j >= 0) v = lx_v(k + 1)
1844820e4eaSBarry Smith        if (i < mx - 1 .and. j > -1) vr = lx_v(k + 2)
1854820e4eaSBarry Smith        if (i > -1 .and. j < my - 1) vt = lx_v(k + 1 + gxm)
186d8606c27SBarry Smith        dvdx = (vr - v)/hx
187d8606c27SBarry Smith        dvdy = (vt - v)/hy
1884820e4eaSBarry Smith        if (i /= -1 .and. j /= -1) then
189d8606c27SBarry Smith          ind = k
190d8606c27SBarry Smith          val = -dvdx/hx - dvdy/hy - cdiv3
1915d83a8b1SBarry Smith          PetscCall(VecSetValuesLocal(G, i1, [k], [val], ADD_VALUES, ierr))
192d8606c27SBarry Smith        end if
1934820e4eaSBarry Smith        if (i /= mx - 1 .and. j /= -1) then
194d8606c27SBarry Smith          ind = k + 1
195d8606c27SBarry Smith          val = dvdx/hx - cdiv3
1965d83a8b1SBarry Smith          PetscCall(VecSetValuesLocal(G, i1, [ind], [val], ADD_VALUES, ierr))
197d8606c27SBarry Smith        end if
1984820e4eaSBarry Smith        if (i /= -1 .and. j /= my - 1) then
199d8606c27SBarry Smith          ind = k + gxm
200d8606c27SBarry Smith          val = dvdy/hy - cdiv3
2015d83a8b1SBarry Smith          PetscCall(VecSetValuesLocal(G, i1, [ind], [val], ADD_VALUES, ierr))
202d8606c27SBarry Smith        end if
203d8606c27SBarry Smith        fquad = fquad + dvdx*dvdx + dvdy*dvdy
204d8606c27SBarry Smith        flin = flin - cdiv3*(v + vr + vt)
205d8606c27SBarry Smith      end do
206d8606c27SBarry Smith    end do
207d8606c27SBarry Smith
208d8606c27SBarry Smith!     Compute local gradient contributions over the upper triangular elements
209d8606c27SBarry Smith
210d8606c27SBarry Smith    do j = ys, yep - 1
211d8606c27SBarry Smith      do i = xs, xep - 1
212d8606c27SBarry Smith        k = (j - gys)*gxm + i - gxs
213d8606c27SBarry Smith        vb = zero
214d8606c27SBarry Smith        vl = zero
215d8606c27SBarry Smith        v = zero
2164820e4eaSBarry Smith        if (i < mx .and. j > 0) vb = lx_v(k + 1 - gxm)
2174820e4eaSBarry Smith        if (i > 0 .and. j < my) vl = lx_v(k)
2184820e4eaSBarry Smith        if (i < mx .and. j < my) v = lx_v(1 + k)
219d8606c27SBarry Smith        dvdx = (v - vl)/hx
220d8606c27SBarry Smith        dvdy = (v - vb)/hy
2214820e4eaSBarry Smith        if (i /= mx .and. j /= 0) then
222d8606c27SBarry Smith          ind = k - gxm
223d8606c27SBarry Smith          val = -dvdy/hy - cdiv3
2245d83a8b1SBarry Smith          PetscCall(VecSetValuesLocal(G, i1, [ind], [val], ADD_VALUES, ierr))
225d8606c27SBarry Smith        end if
2264820e4eaSBarry Smith        if (i /= 0 .and. j /= my) then
227d8606c27SBarry Smith          ind = k - 1
228d8606c27SBarry Smith          val = -dvdx/hx - cdiv3
2295d83a8b1SBarry Smith          PetscCall(VecSetValuesLocal(G, i1, [ind], [val], ADD_VALUES, ierr))
230d8606c27SBarry Smith        end if
2314820e4eaSBarry Smith        if (i /= mx .and. j /= my) then
232d8606c27SBarry Smith          ind = k
233d8606c27SBarry Smith          val = dvdx/hx + dvdy/hy - cdiv3
2345d83a8b1SBarry Smith          PetscCall(VecSetValuesLocal(G, i1, [ind], [val], ADD_VALUES, ierr))
235d8606c27SBarry Smith        end if
236d8606c27SBarry Smith        fquad = fquad + dvdx*dvdx + dvdy*dvdy
237d8606c27SBarry Smith        flin = flin - cdiv3*(vb + vl + v)
238d8606c27SBarry Smith      end do
239d8606c27SBarry Smith    end do
240d8606c27SBarry Smith
241d8606c27SBarry Smith!  Restore vector
242ce78bad3SBarry Smith    PetscCall(VecRestoreArrayRead(localX, lx_v, ierr))
243d8606c27SBarry Smith
244d8606c27SBarry Smith!  Assemble gradient vector
245d8606c27SBarry Smith    PetscCall(VecAssemblyBegin(G, ierr))
246d8606c27SBarry Smith    PetscCall(VecAssemblyEnd(G, ierr))
247d8606c27SBarry Smith
248d8606c27SBarry Smith! Scale the gradient
249d8606c27SBarry Smith    area = p5*hx*hy
250d8606c27SBarry Smith    floc = area*(p5*fquad + flin)
251d8606c27SBarry Smith    PetscCall(VecScale(G, area, ierr))
252d8606c27SBarry Smith
253d8606c27SBarry Smith!  Sum function contributions from all processes
254d8606c27SBarry Smith    PetscCallMPI(MPI_Allreduce(floc, f, 1, MPIU_SCALAR, MPIU_SUM, PETSC_COMM_WORLD, ierr))
255d8606c27SBarry Smith    PetscCall(PetscLogFlops(20.0d0*(ye - ysm)*(xe - xsm) + 16.0d0*(xep - xs)*(yep - ys), ierr))
256d8606c27SBarry Smith  end
257d8606c27SBarry Smith
258ce78bad3SBarry Smith  subroutine ComputeHessian(ta, X, H, Hpre, dummy, ierr)
259ce78bad3SBarry Smith    type(tTao) ta
260392a88c0SBlaise Bourdin    type(tVec) X
261392a88c0SBlaise Bourdin    type(tMat) H, Hpre
262d8606c27SBarry Smith    PetscErrorCode ierr
263d8606c27SBarry Smith    PetscInt dummy
264d8606c27SBarry Smith
265d8606c27SBarry Smith    PetscInt i, j, k
266d8606c27SBarry Smith    PetscInt col(0:4), row
267d8606c27SBarry Smith    PetscInt xs, xm, gxs, gxm
268d8606c27SBarry Smith    PetscInt ys, ym, gys, gym
269d8606c27SBarry Smith    PetscReal v(0:4)
270d8606c27SBarry Smith    PetscInt i1
271d8606c27SBarry Smith
272d8606c27SBarry Smith    i1 = 1
273d8606c27SBarry Smith
274d8606c27SBarry Smith!     Get local grid boundaries
275d8606c27SBarry Smith    PetscCall(DMDAGetCorners(dm, xs, ys, PETSC_NULL_INTEGER, xm, ym, PETSC_NULL_INTEGER, ierr))
276d8606c27SBarry Smith    PetscCall(DMDAGetGhostCorners(dm, gxs, gys, PETSC_NULL_INTEGER, gxm, gym, PETSC_NULL_INTEGER, ierr))
277d8606c27SBarry Smith
278d8606c27SBarry Smith    do j = ys, ys + ym - 1
279d8606c27SBarry Smith      do i = xs, xs + xm - 1
280d8606c27SBarry Smith        row = (j - gys)*gxm + (i - gxs)
281d8606c27SBarry Smith
282d8606c27SBarry Smith        k = 0
2834820e4eaSBarry Smith        if (j > gys) then
284d8606c27SBarry Smith          v(k) = -1.0
285d8606c27SBarry Smith          col(k) = row - gxm
286d8606c27SBarry Smith          k = k + 1
287d8606c27SBarry Smith        end if
288d8606c27SBarry Smith
2894820e4eaSBarry Smith        if (i > gxs) then
290d8606c27SBarry Smith          v(k) = -1.0
291d8606c27SBarry Smith          col(k) = row - 1
292d8606c27SBarry Smith          k = k + 1
293d8606c27SBarry Smith        end if
294d8606c27SBarry Smith
295d8606c27SBarry Smith        v(k) = 4.0
296d8606c27SBarry Smith        col(k) = row
297d8606c27SBarry Smith        k = k + 1
298d8606c27SBarry Smith
2994820e4eaSBarry Smith        if (i + 1 < gxs + gxm) then
300d8606c27SBarry Smith          v(k) = -1.0
301d8606c27SBarry Smith          col(k) = row + 1
302d8606c27SBarry Smith          k = k + 1
303d8606c27SBarry Smith        end if
304d8606c27SBarry Smith
3054820e4eaSBarry Smith        if (j + 1 < gys + gym) then
306d8606c27SBarry Smith          v(k) = -1.0
307d8606c27SBarry Smith          col(k) = row + gxm
308d8606c27SBarry Smith          k = k + 1
309d8606c27SBarry Smith        end if
310d8606c27SBarry Smith
3115d83a8b1SBarry Smith        PetscCall(MatSetValuesLocal(H, i1, [row], k, col, v, INSERT_VALUES, ierr))
312d8606c27SBarry Smith      end do
313d8606c27SBarry Smith    end do
314d8606c27SBarry Smith
315d8606c27SBarry Smith!     Assemble matrix
316d8606c27SBarry Smith    PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY, ierr))
317d8606c27SBarry Smith    PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY, ierr))
318d8606c27SBarry Smith
319d8606c27SBarry Smith!     Tell the matrix we will never add a new nonzero location to the
320d8606c27SBarry Smith!     matrix.  If we do it will generate an error.
321d8606c27SBarry Smith
322d8606c27SBarry Smith    PetscCall(MatSetOption(H, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE, ierr))
323d8606c27SBarry Smith    PetscCall(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE, ierr))
324d8606c27SBarry Smith
325d8606c27SBarry Smith    PetscCall(PetscLogFlops(9.0d0*xm*ym + 49.0d0*xm, ierr))
326d8606c27SBarry Smith
327d8606c27SBarry Smith    ierr = 0
328d8606c27SBarry Smith  end
329d8606c27SBarry Smith
330ce78bad3SBarry Smith  subroutine Monitor(ta, dummy, ierr)
331ce78bad3SBarry Smith    type(tTao) ta
332d8606c27SBarry Smith    PetscInt dummy
333d8606c27SBarry Smith    PetscErrorCode ierr
334d8606c27SBarry Smith
335d8606c27SBarry Smith    PetscInt its
336d8606c27SBarry Smith    PetscReal f, gnorm, cnorm, xdiff
337d8606c27SBarry Smith    TaoConvergedReason reason
338d8606c27SBarry Smith
339ce78bad3SBarry Smith    PetscCall(TaoGetSolutionStatus(ta, its, f, gnorm, cnorm, xdiff, reason, ierr))
3404820e4eaSBarry Smith    if (mod(its, 5) /= 0) then
341d8606c27SBarry Smith      PetscCall(PetscPrintf(PETSC_COMM_WORLD, 'iteration multiple of 5\n', ierr))
342d8606c27SBarry Smith    end if
343d8606c27SBarry Smith
344d8606c27SBarry Smith    ierr = 0
345d8606c27SBarry Smith
346d8606c27SBarry Smith  end
347d8606c27SBarry Smith
348ce78bad3SBarry Smith  subroutine ConvergenceTest(ta, dummy, ierr)
349ce78bad3SBarry Smith    type(tTao) ta
350d8606c27SBarry Smith    PetscInt dummy
351d8606c27SBarry Smith    PetscErrorCode ierr
352d8606c27SBarry Smith
353d8606c27SBarry Smith    PetscInt its
354d8606c27SBarry Smith    PetscReal f, gnorm, cnorm, xdiff
355d8606c27SBarry Smith    TaoConvergedReason reason
356d8606c27SBarry Smith
357ce78bad3SBarry Smith    PetscCall(TaoGetSolutionStatus(ta, its, f, gnorm, cnorm, xdiff, reason, ierr))
3584820e4eaSBarry Smith    if (its == 7) then
359ce78bad3SBarry Smith      PetscCall(TaoSetConvergedReason(ta, TAO_DIVERGED_MAXITS, ierr))
360d8606c27SBarry Smith    end if
361d8606c27SBarry Smith
362d8606c27SBarry Smith    ierr = 0
363d8606c27SBarry Smith
364d8606c27SBarry Smith  end
365*e7a95102SMartin Diehlend module
366*e7a95102SMartin Diehl
367*e7a95102SMartin Diehlprogram eptorsion2f
368*e7a95102SMartin Diehl  use eptorsion2fmodule
369*e7a95102SMartin Diehl  implicit none
370*e7a95102SMartin Diehl! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
371*e7a95102SMartin Diehl!                   Variable declarations
372*e7a95102SMartin Diehl! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
373*e7a95102SMartin Diehl!
374*e7a95102SMartin Diehl!  See additional variable declarations in the file eptorsion2f.h
375*e7a95102SMartin Diehl!
376*e7a95102SMartin Diehl  PetscErrorCode ierr           ! used to check for functions returning nonzeros
377*e7a95102SMartin Diehl  type(tVec) x              ! solution vector
378*e7a95102SMartin Diehl  type(tMat) H              ! hessian matrix
379*e7a95102SMartin Diehl  PetscInt Nx, Ny         ! number of processes in x- and y- directions
380*e7a95102SMartin Diehl  type(tTao) ta            ! Tao solver context
381*e7a95102SMartin Diehl  PetscBool flg
382*e7a95102SMartin Diehl  PetscInt i1
383*e7a95102SMartin Diehl  PetscInt dummy
384*e7a95102SMartin Diehl
385*e7a95102SMartin Diehl  i1 = 1
386*e7a95102SMartin Diehl
387*e7a95102SMartin Diehl!     Initialize TAO, PETSc  contexts
388*e7a95102SMartin Diehl  PetscCallA(PetscInitialize(ierr))
389*e7a95102SMartin Diehl
390*e7a95102SMartin Diehl!     Specify default parameters
391*e7a95102SMartin Diehl  param = 5.0
392*e7a95102SMartin Diehl  mx = 10
393*e7a95102SMartin Diehl  my = 10
394*e7a95102SMartin Diehl  Nx = PETSC_DECIDE
395*e7a95102SMartin Diehl  Ny = PETSC_DECIDE
396*e7a95102SMartin Diehl
397*e7a95102SMartin Diehl!     Check for any command line arguments that might override defaults
398*e7a95102SMartin Diehl  PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-mx', mx, flg, ierr))
399*e7a95102SMartin Diehl  PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-my', my, flg, ierr))
400*e7a95102SMartin Diehl  PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-par', param, flg, ierr))
401*e7a95102SMartin Diehl
402*e7a95102SMartin Diehl!     Set up distributed array and vectors
403*e7a95102SMartin Diehl  PetscCallA(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, mx, my, Nx, Ny, i1, i1, PETSC_NULL_INTEGER_ARRAY, PETSC_NULL_INTEGER_ARRAY, dm, ierr))
404*e7a95102SMartin Diehl  PetscCallA(DMSetFromOptions(dm, ierr))
405*e7a95102SMartin Diehl  PetscCallA(DMSetUp(dm, ierr))
406*e7a95102SMartin Diehl
407*e7a95102SMartin Diehl!     Create vectors
408*e7a95102SMartin Diehl  PetscCallA(DMCreateGlobalVector(dm, x, ierr))
409*e7a95102SMartin Diehl  PetscCallA(DMCreateLocalVector(dm, localX, ierr))
410*e7a95102SMartin Diehl
411*e7a95102SMartin Diehl!     Create Hessian
412*e7a95102SMartin Diehl  PetscCallA(DMCreateMatrix(dm, H, ierr))
413*e7a95102SMartin Diehl  PetscCallA(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE, ierr))
414*e7a95102SMartin Diehl
415*e7a95102SMartin Diehl!     The TAO code begins here
416*e7a95102SMartin Diehl
417*e7a95102SMartin Diehl!     Create TAO solver
418*e7a95102SMartin Diehl  PetscCallA(TaoCreate(PETSC_COMM_WORLD, ta, ierr))
419*e7a95102SMartin Diehl  PetscCallA(TaoSetType(ta, TAOCG, ierr))
420*e7a95102SMartin Diehl
421*e7a95102SMartin Diehl!     Set routines for function and gradient evaluation
422*e7a95102SMartin Diehl
423*e7a95102SMartin Diehl  PetscCallA(TaoSetObjectiveAndGradient(ta, PETSC_NULL_VEC, FormFunctionGradient, 0, ierr))
424*e7a95102SMartin Diehl  PetscCallA(TaoSetHessian(ta, H, H, ComputeHessian, 0, ierr))
425*e7a95102SMartin Diehl
426*e7a95102SMartin Diehl!     Set initial guess
427*e7a95102SMartin Diehl  PetscCallA(FormInitialGuess(x, ierr))
428*e7a95102SMartin Diehl  PetscCallA(TaoSetSolution(ta, x, ierr))
429*e7a95102SMartin Diehl
430*e7a95102SMartin Diehl  PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-testmonitor', flg, ierr))
431*e7a95102SMartin Diehl  if (flg) then
432*e7a95102SMartin Diehl    PetscCallA(TaoMonitorSet(ta, Monitor, dummy, PETSC_NULL_FUNCTION, ierr))
433*e7a95102SMartin Diehl  end if
434*e7a95102SMartin Diehl
435*e7a95102SMartin Diehl  PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-testconvergence', flg, ierr))
436*e7a95102SMartin Diehl  if (flg) then
437*e7a95102SMartin Diehl    PetscCallA(TaoSetConvergenceTest(ta, ConvergenceTest, dummy, ierr))
438*e7a95102SMartin Diehl  end if
439*e7a95102SMartin Diehl
440*e7a95102SMartin Diehl!     Check for any TAO command line options
441*e7a95102SMartin Diehl  PetscCallA(TaoSetFromOptions(ta, ierr))
442*e7a95102SMartin Diehl
443*e7a95102SMartin Diehl!     SOLVE THE APPLICATION
444*e7a95102SMartin Diehl  PetscCallA(TaoSolve(ta, ierr))
445*e7a95102SMartin Diehl
446*e7a95102SMartin Diehl!     Free TAO data structures
447*e7a95102SMartin Diehl  PetscCallA(TaoDestroy(ta, ierr))
448*e7a95102SMartin Diehl
449*e7a95102SMartin Diehl!     Free PETSc data structures
450*e7a95102SMartin Diehl  PetscCallA(VecDestroy(x, ierr))
451*e7a95102SMartin Diehl  PetscCallA(VecDestroy(localX, ierr))
452*e7a95102SMartin Diehl  PetscCallA(MatDestroy(H, ierr))
453*e7a95102SMartin Diehl  PetscCallA(DMDestroy(dm, ierr))
454*e7a95102SMartin Diehl
455*e7a95102SMartin Diehl!     Finalize TAO and PETSc
456*e7a95102SMartin Diehl  PetscCallA(PetscFinalize(ierr))
457*e7a95102SMartin Diehlend
458d8606c27SBarry Smith
459d8606c27SBarry Smith!/*TEST
460d8606c27SBarry Smith!
461d8606c27SBarry Smith!   build:
462d8606c27SBarry Smith!      requires: !complex
463d8606c27SBarry Smith!
464d8606c27SBarry Smith!   test:
46510978b7dSBarry Smith!      args: -tao_monitor_short -tao_type nls -tao_gttol 1.e-2
466d8606c27SBarry Smith!
467d8606c27SBarry Smith!   test:
468d8606c27SBarry Smith!      suffix: 2
469d8606c27SBarry Smith!      nsize: 2
47010978b7dSBarry Smith!      args: -tao_monitor_short -tao_type lmvm -tao_gttol 1.e-2
471d8606c27SBarry Smith!
472d8606c27SBarry Smith!   test:
473d8606c27SBarry Smith!      suffix: 3
474d8606c27SBarry Smith!      args: -testmonitor -tao_type lmvm -tao_gttol 1.e-2
475d8606c27SBarry Smith!TEST*/
476