xref: /petsc/src/tao/bound/tutorials/plate2f.F90 (revision 4820e4ea99a084ae862a8c395f732bc7c0e1a6d0)
1c4762a1bSJed Brown!  Program usage: mpiexec -n <proc> plate2f [all TAO options]
2c4762a1bSJed Brown!
3c4762a1bSJed Brown!  This example demonstrates use of the TAO package to solve a bound constrained
4c4762a1bSJed Brown!  minimization problem.  This example is based on a problem from the
5c4762a1bSJed Brown!  MINPACK-2 test suite.  Given a rectangular 2-D domain and boundary values
6c4762a1bSJed Brown!  along the edges of the domain, the objective is to find the surface
7c4762a1bSJed Brown!  with the minimal area that satisfies the boundary conditions.
8c4762a1bSJed Brown!  The command line options are:
9c4762a1bSJed Brown!    -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction
10c4762a1bSJed Brown!    -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction
11c4762a1bSJed Brown!    -bmx <bxg>, where <bxg> = number of grid points under plate in 1st direction
12c4762a1bSJed Brown!    -bmy <byg>, where <byg> = number of grid points under plate in 2nd direction
13c4762a1bSJed Brown!    -bheight <ht>, where <ht> = height of the plate
14c4762a1bSJed Brown!
15c4762a1bSJed Brown
1677d968b7SBarry Smithmodule plate2fmodule
17c4762a1bSJed Brown#include "petsc/finclude/petscdmda.h"
18c4762a1bSJed Brown#include "petsc/finclude/petsctao.h"
19c4762a1bSJed Brown  use petscdmda
20c4762a1bSJed Brown  use petsctao
21c4762a1bSJed Brown
22c4762a1bSJed Brown  Vec localX, localV
23c4762a1bSJed Brown  Vec Top, Left
24c4762a1bSJed Brown  Vec Right, Bottom
25c4762a1bSJed Brown  DM dm
26c4762a1bSJed Brown  PetscReal bheight
27c4762a1bSJed Brown  PetscInt bmx, bmy
28c4762a1bSJed Brown  PetscInt mx, my, Nx, Ny, N
29c4762a1bSJed Brownend module
30c4762a1bSJed Brown
31c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
32c4762a1bSJed Brown!                   Variable declarations
33c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34c4762a1bSJed Brown!
35c4762a1bSJed Brown!  Variables:
36c4762a1bSJed Brown!    (common from plate2f.h):
37c4762a1bSJed Brown!    Nx, Ny           number of processors in x- and y- directions
38c4762a1bSJed Brown!    mx, my           number of grid points in x,y directions
39c4762a1bSJed Brown!    N    global dimension of vector
4077d968b7SBarry Smithuse plate2fmodule
41c4762a1bSJed Brownimplicit none
42c4762a1bSJed Brown
43c4762a1bSJed BrownPetscErrorCode ierr          ! used to check for functions returning nonzeros
44c4762a1bSJed BrownVec x             ! solution vector
45c4762a1bSJed BrownPetscInt m             ! number of local elements in vector
46ce78bad3SBarry SmithTao ta           ! Tao solver context
47c4762a1bSJed BrownMat H             ! Hessian matrix
48c4762a1bSJed BrownISLocalToGlobalMapping isltog  ! local to global mapping object
49c4762a1bSJed BrownPetscBool flg
50c4762a1bSJed BrownPetscInt i1, i3, i7
51c4762a1bSJed Brown
52c4762a1bSJed Brownexternal FormFunctionGradient
53c4762a1bSJed Brownexternal FormHessian
54c4762a1bSJed Brownexternal MSA_BoundaryConditions
55c4762a1bSJed Brownexternal MSA_Plate
56c4762a1bSJed Brownexternal MSA_InitialPoint
57c4762a1bSJed Brown! Initialize Tao
58c4762a1bSJed Brown
59c4762a1bSJed Browni1 = 1
60c4762a1bSJed Browni3 = 3
61c4762a1bSJed Browni7 = 7
62c4762a1bSJed Brown
63d8606c27SBarry SmithPetscCallA(PetscInitialize(ierr))
64c4762a1bSJed Brown
65c4762a1bSJed Brown! Specify default dimensions of the problem
66c4762a1bSJed Brownmx = 10
67c4762a1bSJed Brownmy = 10
68c4762a1bSJed Brownbheight = 0.1
69c4762a1bSJed Brown
70c4762a1bSJed Brown! Check for any command line arguments that override defaults
71c4762a1bSJed Brown
72d8606c27SBarry SmithPetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-mx', mx, flg, ierr))
73d8606c27SBarry SmithPetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-my', my, flg, ierr))
74c4762a1bSJed Brown
75c4762a1bSJed Brownbmx = mx/2
76c4762a1bSJed Brownbmy = my/2
77c4762a1bSJed Brown
78d8606c27SBarry SmithPetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-bmx', bmx, flg, ierr))
79d8606c27SBarry SmithPetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-bmy', bmy, flg, ierr))
80d8606c27SBarry SmithPetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-bheight', bheight, flg, ierr))
81c4762a1bSJed Brown
82c4762a1bSJed Brown! Calculate any derived values from parameters
83c4762a1bSJed BrownN = mx*my
84c4762a1bSJed Brown
85f0b74427SPierre Jolivet! Let PETSc determine the dimensions of the local vectors
86c4762a1bSJed BrownNx = PETSC_DECIDE
87c4762a1bSJed BrownNY = PETSC_DECIDE
88c4762a1bSJed Brown
89c4762a1bSJed Brown! A two dimensional distributed array will help define this problem, which
90c4762a1bSJed Brown! derives from an elliptic PDE on a two-dimensional domain.  From the
91c4762a1bSJed Brown! distributed array, create the vectors
92c4762a1bSJed Brown
935d83a8b1SBarry SmithPetscCallA(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))
94d8606c27SBarry SmithPetscCallA(DMSetFromOptions(dm, ierr))
95d8606c27SBarry SmithPetscCallA(DMSetUp(dm, ierr))
96c4762a1bSJed Brown
97c4762a1bSJed Brown! Extract global and local vectors from DM; The local vectors are
98c4762a1bSJed Brown! used solely as work space for the evaluation of the function,
99c4762a1bSJed Brown! gradient, and Hessian.  Duplicate for remaining vectors that are
100c4762a1bSJed Brown! the same types.
101c4762a1bSJed Brown
102d8606c27SBarry SmithPetscCallA(DMCreateGlobalVector(dm, x, ierr))
103d8606c27SBarry SmithPetscCallA(DMCreateLocalVector(dm, localX, ierr))
104d8606c27SBarry SmithPetscCallA(VecDuplicate(localX, localV, ierr))
105c4762a1bSJed Brown
106c4762a1bSJed Brown! Create a matrix data structure to store the Hessian.
107c4762a1bSJed Brown! Here we (optionally) also associate the local numbering scheme
108c4762a1bSJed Brown! with the matrix so that later we can use local indices for matrix
109c4762a1bSJed Brown! assembly
110c4762a1bSJed Brown
111d8606c27SBarry SmithPetscCallA(VecGetLocalSize(x, m, ierr))
1125d83a8b1SBarry SmithPetscCallA(MatCreateAIJ(PETSC_COMM_WORLD, m, m, N, N, i7, PETSC_NULL_INTEGER_ARRAY, i3, PETSC_NULL_INTEGER_ARRAY, H, ierr))
113c4762a1bSJed Brown
114d8606c27SBarry SmithPetscCallA(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE, ierr))
115d8606c27SBarry SmithPetscCallA(DMGetLocalToGlobalMapping(dm, isltog, ierr))
116d8606c27SBarry SmithPetscCallA(MatSetLocalToGlobalMapping(H, isltog, isltog, ierr))
117c4762a1bSJed Brown
118c4762a1bSJed Brown! The Tao code begins here
119c4762a1bSJed Brown! Create TAO solver and set desired solution method.
120c4762a1bSJed Brown! This problems uses bounded variables, so the
121c4762a1bSJed Brown! method must either be 'tao_tron' or 'tao_blmvm'
122c4762a1bSJed Brown
123ce78bad3SBarry SmithPetscCallA(TaoCreate(PETSC_COMM_WORLD, ta, ierr))
124ce78bad3SBarry SmithPetscCallA(TaoSetType(ta, TAOBLMVM, ierr))
125c4762a1bSJed Brown
126c4762a1bSJed Brown!     Set minimization function and gradient, hessian evaluation functions
127c4762a1bSJed Brown
128ce78bad3SBarry SmithPetscCallA(TaoSetObjectiveAndGradient(ta, PETSC_NULL_VEC, FormFunctionGradient, 0, ierr))
129c4762a1bSJed Brown
130ce78bad3SBarry SmithPetscCallA(TaoSetHessian(ta, H, H, FormHessian, 0, ierr))
131c4762a1bSJed Brown
132c4762a1bSJed Brown! Set Variable bounds
133d8606c27SBarry SmithPetscCallA(MSA_BoundaryConditions(ierr))
134ce78bad3SBarry SmithPetscCallA(TaoSetVariableBoundsRoutine(ta, MSA_Plate, 0, ierr))
135c4762a1bSJed Brown
136c4762a1bSJed Brown! Set the initial solution guess
137d8606c27SBarry SmithPetscCallA(MSA_InitialPoint(x, ierr))
138ce78bad3SBarry SmithPetscCallA(TaoSetSolution(ta, x, ierr))
139c4762a1bSJed Brown
140c4762a1bSJed Brown! Check for any tao command line options
141ce78bad3SBarry SmithPetscCallA(TaoSetFromOptions(ta, ierr))
142c4762a1bSJed Brown
143c4762a1bSJed Brown! Solve the application
144ce78bad3SBarry SmithPetscCallA(TaoSolve(ta, ierr))
145c4762a1bSJed Brown
146c4762a1bSJed Brown! Free TAO data structures
147ce78bad3SBarry SmithPetscCallA(TaoDestroy(ta, ierr))
148c4762a1bSJed Brown
149c4762a1bSJed Brown! Free PETSc data structures
150d8606c27SBarry SmithPetscCallA(VecDestroy(x, ierr))
151d8606c27SBarry SmithPetscCallA(VecDestroy(Top, ierr))
152d8606c27SBarry SmithPetscCallA(VecDestroy(Bottom, ierr))
153d8606c27SBarry SmithPetscCallA(VecDestroy(Left, ierr))
154d8606c27SBarry SmithPetscCallA(VecDestroy(Right, ierr))
155d8606c27SBarry SmithPetscCallA(MatDestroy(H, ierr))
156d8606c27SBarry SmithPetscCallA(VecDestroy(localX, ierr))
157d8606c27SBarry SmithPetscCallA(VecDestroy(localV, ierr))
158d8606c27SBarry SmithPetscCallA(DMDestroy(dm, ierr))
159c4762a1bSJed Brown
160c4762a1bSJed Brown! Finalize TAO
161c4762a1bSJed Brown
162d8606c27SBarry SmithPetscCallA(PetscFinalize(ierr))
163c4762a1bSJed Brown
164c4762a1bSJed Brownend
165c4762a1bSJed Brown
166c4762a1bSJed Brown! ---------------------------------------------------------------------
167c4762a1bSJed Brown!
168c4762a1bSJed Brown!  FormFunctionGradient - Evaluates function f(X).
169c4762a1bSJed Brown!
170c4762a1bSJed Brown!  Input Parameters:
171c4762a1bSJed Brown!  tao   - the Tao context
172c4762a1bSJed Brown!  X     - the input vector
173c4762a1bSJed Brown!  dummy - optional user-defined context, as set by TaoSetFunction()
174c4762a1bSJed Brown!          (not used here)
175c4762a1bSJed Brown!
176c4762a1bSJed Brown!  Output Parameters:
177c4762a1bSJed Brown!  fcn     - the newly evaluated function
178c4762a1bSJed Brown!  G       - the gradient vector
179c4762a1bSJed Brown!  info  - error code
180c4762a1bSJed Brown!
181c4762a1bSJed Brown
182ce78bad3SBarry Smithsubroutine FormFunctionGradient(ta, X, fcn, G, dummy, ierr)
18377d968b7SBarry Smith  use plate2fmodule
184c4762a1bSJed Brown  implicit none
185c4762a1bSJed Brown
186c4762a1bSJed Brown! Input/output variables
187c4762a1bSJed Brown
188ce78bad3SBarry Smith  Tao ta
189c4762a1bSJed Brown  PetscReal fcn
190c4762a1bSJed Brown  Vec X, G
191c4762a1bSJed Brown  PetscErrorCode ierr
192c4762a1bSJed Brown  PetscInt dummy
193c4762a1bSJed Brown
194c4762a1bSJed Brown  PetscInt i, j, row
195c4762a1bSJed Brown  PetscInt xs, xm
196c4762a1bSJed Brown  PetscInt gxs, gxm
197c4762a1bSJed Brown  PetscInt ys, ym
198c4762a1bSJed Brown  PetscInt gys, gym
199c4762a1bSJed Brown  PetscReal ft, zero, hx, hy, hydhx, hxdhy
200c4762a1bSJed Brown  PetscReal area, rhx, rhy
201c4762a1bSJed Brown  PetscReal f1, f2, f3, f4, f5, f6, d1, d2, d3
202c4762a1bSJed Brown  PetscReal d4, d5, d6, d7, d8
203c4762a1bSJed Brown  PetscReal df1dxc, df2dxc, df3dxc, df4dxc
204c4762a1bSJed Brown  PetscReal df5dxc, df6dxc
205c4762a1bSJed Brown  PetscReal xc, xl, xr, xt, xb, xlt, xrb
206c4762a1bSJed Brown
20742ce371bSBarry Smith  PetscReal, pointer :: g_v(:), x_v(:)
20842ce371bSBarry Smith  PetscReal, pointer :: top_v(:), left_v(:)
20942ce371bSBarry Smith  PetscReal, pointer :: right_v(:), bottom_v(:)
210c4762a1bSJed Brown
211c4762a1bSJed Brown  ft = 0.0
212c4762a1bSJed Brown  zero = 0.0
213c4762a1bSJed Brown  hx = 1.0/real(mx + 1)
214c4762a1bSJed Brown  hy = 1.0/real(my + 1)
215c4762a1bSJed Brown  hydhx = hy/hx
216c4762a1bSJed Brown  hxdhy = hx/hy
217c4762a1bSJed Brown  area = 0.5*hx*hy
218c4762a1bSJed Brown  rhx = real(mx) + 1.0
219c4762a1bSJed Brown  rhy = real(my) + 1.0
220c4762a1bSJed Brown
221c4762a1bSJed Brown! Get local mesh boundaries
222d8606c27SBarry Smith  PetscCall(DMDAGetCorners(dm, xs, ys, PETSC_NULL_INTEGER, xm, ym, PETSC_NULL_INTEGER, ierr))
223d8606c27SBarry Smith  PetscCall(DMDAGetGhostCorners(dm, gxs, gys, PETSC_NULL_INTEGER, gxm, gym, PETSC_NULL_INTEGER, ierr))
224c4762a1bSJed Brown
225c4762a1bSJed Brown! Scatter ghost points to local vector
226d8606c27SBarry Smith  PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX, ierr))
227d8606c27SBarry Smith  PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX, ierr))
228c4762a1bSJed Brown
229c4762a1bSJed Brown! Initialize the vector to zero
230d8606c27SBarry Smith  PetscCall(VecSet(localV, zero, ierr))
231c4762a1bSJed Brown
232ce78bad3SBarry Smith! Get arrays to vector data (See note above about using VecGetArray in Fortran)
233ce78bad3SBarry Smith  PetscCall(VecGetArray(localX, x_v, ierr))
234ce78bad3SBarry Smith  PetscCall(VecGetArray(localV, g_v, ierr))
235ce78bad3SBarry Smith  PetscCall(VecGetArray(Top, top_v, ierr))
236ce78bad3SBarry Smith  PetscCall(VecGetArray(Bottom, bottom_v, ierr))
237ce78bad3SBarry Smith  PetscCall(VecGetArray(Left, left_v, ierr))
238ce78bad3SBarry Smith  PetscCall(VecGetArray(Right, right_v, ierr))
239c4762a1bSJed Brown
240c4762a1bSJed Brown! Compute function over the locally owned part of the mesh
241c4762a1bSJed Brown  do j = ys, ys + ym - 1
242c4762a1bSJed Brown    do i = xs, xs + xm - 1
243c4762a1bSJed Brown      row = (j - gys)*gxm + (i - gxs)
24442ce371bSBarry Smith      xc = x_v(1 + row)
245c4762a1bSJed Brown      xt = xc
246c4762a1bSJed Brown      xb = xc
247c4762a1bSJed Brown      xr = xc
248c4762a1bSJed Brown      xl = xc
249c4762a1bSJed Brown      xrb = xc
250c4762a1bSJed Brown      xlt = xc
251c4762a1bSJed Brown
252*4820e4eaSBarry Smith      if (i == 0) then !left side
25342ce371bSBarry Smith        xl = left_v(1 + j - ys + 1)
25442ce371bSBarry Smith        xlt = left_v(1 + j - ys + 2)
255c4762a1bSJed Brown      else
25642ce371bSBarry Smith        xl = x_v(1 + row - 1)
257c4762a1bSJed Brown      end if
258c4762a1bSJed Brown
259*4820e4eaSBarry Smith      if (j == 0) then !bottom side
26042ce371bSBarry Smith        xb = bottom_v(1 + i - xs + 1)
26142ce371bSBarry Smith        xrb = bottom_v(1 + i - xs + 2)
262c4762a1bSJed Brown      else
26342ce371bSBarry Smith        xb = x_v(1 + row - gxm)
264c4762a1bSJed Brown      end if
265c4762a1bSJed Brown
266*4820e4eaSBarry Smith      if (i + 1 == gxs + gxm) then !right side
26742ce371bSBarry Smith        xr = right_v(1 + j - ys + 1)
26842ce371bSBarry Smith        xrb = right_v(1 + j - ys)
269c4762a1bSJed Brown      else
27042ce371bSBarry Smith        xr = x_v(1 + row + 1)
271c4762a1bSJed Brown      end if
272c4762a1bSJed Brown
273*4820e4eaSBarry Smith      if (j + 1 == gys + gym) then !top side
27442ce371bSBarry Smith        xt = top_v(1 + i - xs + 1)
27542ce371bSBarry Smith        xlt = top_v(1 + i - xs)
276c4762a1bSJed Brown      else
27742ce371bSBarry Smith        xt = x_v(1 + row + gxm)
278c4762a1bSJed Brown      end if
279c4762a1bSJed Brown
280*4820e4eaSBarry Smith      if ((i > gxs) .and. (j + 1 < gys + gym)) then
28142ce371bSBarry Smith        xlt = x_v(1 + row - 1 + gxm)
282c4762a1bSJed Brown      end if
283c4762a1bSJed Brown
284*4820e4eaSBarry Smith      if ((j > gys) .and. (i + 1 < gxs + gxm)) then
28542ce371bSBarry Smith        xrb = x_v(1 + row + 1 - gxm)
286c4762a1bSJed Brown      end if
287c4762a1bSJed Brown
288c4762a1bSJed Brown      d1 = xc - xl
289c4762a1bSJed Brown      d2 = xc - xr
290c4762a1bSJed Brown      d3 = xc - xt
291c4762a1bSJed Brown      d4 = xc - xb
292c4762a1bSJed Brown      d5 = xr - xrb
293c4762a1bSJed Brown      d6 = xrb - xb
294c4762a1bSJed Brown      d7 = xlt - xl
295c4762a1bSJed Brown      d8 = xt - xlt
296c4762a1bSJed Brown
297c4762a1bSJed Brown      df1dxc = d1*hydhx
298c4762a1bSJed Brown      df2dxc = d1*hydhx + d4*hxdhy
299c4762a1bSJed Brown      df3dxc = d3*hxdhy
300c4762a1bSJed Brown      df4dxc = d2*hydhx + d3*hxdhy
301c4762a1bSJed Brown      df5dxc = d2*hydhx
302c4762a1bSJed Brown      df6dxc = d4*hxdhy
303c4762a1bSJed Brown
304c4762a1bSJed Brown      d1 = d1*rhx
305c4762a1bSJed Brown      d2 = d2*rhx
306c4762a1bSJed Brown      d3 = d3*rhy
307c4762a1bSJed Brown      d4 = d4*rhy
308c4762a1bSJed Brown      d5 = d5*rhy
309c4762a1bSJed Brown      d6 = d6*rhx
310c4762a1bSJed Brown      d7 = d7*rhy
311c4762a1bSJed Brown      d8 = d8*rhx
312c4762a1bSJed Brown
313c4762a1bSJed Brown      f1 = sqrt(1.0 + d1*d1 + d7*d7)
314c4762a1bSJed Brown      f2 = sqrt(1.0 + d1*d1 + d4*d4)
315c4762a1bSJed Brown      f3 = sqrt(1.0 + d3*d3 + d8*d8)
316c4762a1bSJed Brown      f4 = sqrt(1.0 + d3*d3 + d2*d2)
317c4762a1bSJed Brown      f5 = sqrt(1.0 + d2*d2 + d5*d5)
318c4762a1bSJed Brown      f6 = sqrt(1.0 + d4*d4 + d6*d6)
319c4762a1bSJed Brown
320c4762a1bSJed Brown      ft = ft + f2 + f4
321c4762a1bSJed Brown
322c4762a1bSJed Brown      df1dxc = df1dxc/f1
323c4762a1bSJed Brown      df2dxc = df2dxc/f2
324c4762a1bSJed Brown      df3dxc = df3dxc/f3
325c4762a1bSJed Brown      df4dxc = df4dxc/f4
326c4762a1bSJed Brown      df5dxc = df5dxc/f5
327c4762a1bSJed Brown      df6dxc = df6dxc/f6
328c4762a1bSJed Brown
32942ce371bSBarry Smith      g_v(1 + row) = 0.5*(df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc)
330c4762a1bSJed Brown    end do
331c4762a1bSJed Brown  end do
332c4762a1bSJed Brown
333c4762a1bSJed Brown! Compute triangular areas along the border of the domain.
334*4820e4eaSBarry Smith  if (xs == 0) then  ! left side
335c4762a1bSJed Brown    do j = ys, ys + ym - 1
33642ce371bSBarry Smith      d3 = (left_v(1 + j - ys + 1) - left_v(1 + j - ys + 2))*rhy
33742ce371bSBarry Smith      d2 = (left_v(1 + j - ys + 1) - x_v(1 + (j - gys)*gxm))*rhx
338c4762a1bSJed Brown      ft = ft + sqrt(1.0 + d3*d3 + d2*d2)
339c4762a1bSJed Brown    end do
340c4762a1bSJed Brown  end if
341c4762a1bSJed Brown
342*4820e4eaSBarry Smith  if (ys == 0) then !bottom side
343c4762a1bSJed Brown    do i = xs, xs + xm - 1
34442ce371bSBarry Smith      d2 = (bottom_v(1 + i + 1 - xs) - bottom_v(1 + i - xs + 2))*rhx
34542ce371bSBarry Smith      d3 = (bottom_v(1 + i - xs + 1) - x_v(1 + i - gxs))*rhy
346c4762a1bSJed Brown      ft = ft + sqrt(1.0 + d3*d3 + d2*d2)
347c4762a1bSJed Brown    end do
348c4762a1bSJed Brown  end if
349c4762a1bSJed Brown
350*4820e4eaSBarry Smith  if (xs + xm == mx) then ! right side
351c4762a1bSJed Brown    do j = ys, ys + ym - 1
35242ce371bSBarry Smith      d1 = (x_v(1 + (j + 1 - gys)*gxm - 1) - right_v(1 + j - ys + 1))*rhx
35342ce371bSBarry Smith      d4 = (right_v(1 + j - ys) - right_v(1 + j - ys + 1))*rhy
354c4762a1bSJed Brown      ft = ft + sqrt(1.0 + d1*d1 + d4*d4)
355c4762a1bSJed Brown    end do
356c4762a1bSJed Brown  end if
357c4762a1bSJed Brown
358*4820e4eaSBarry Smith  if (ys + ym == my) then
359c4762a1bSJed Brown    do i = xs, xs + xm - 1
36042ce371bSBarry Smith      d1 = (x_v(1 + (gym - 1)*gxm + i - gxs) - top_v(1 + i - xs + 1))*rhy
36142ce371bSBarry Smith      d4 = (top_v(1 + i - xs + 1) - top_v(1 + i - xs))*rhx
362c4762a1bSJed Brown      ft = ft + sqrt(1.0 + d1*d1 + d4*d4)
363c4762a1bSJed Brown    end do
364c4762a1bSJed Brown  end if
365c4762a1bSJed Brown
366*4820e4eaSBarry Smith  if ((ys == 0) .and. (xs == 0)) then
36742ce371bSBarry Smith    d1 = (left_v(1 + 0) - left_v(1 + 1))*rhy
36842ce371bSBarry Smith    d2 = (bottom_v(1 + 0) - bottom_v(1 + 1))*rhx
369c4762a1bSJed Brown    ft = ft + sqrt(1.0 + d1*d1 + d2*d2)
370c4762a1bSJed Brown  end if
371c4762a1bSJed Brown
372*4820e4eaSBarry Smith  if ((ys + ym == my) .and. (xs + xm == mx)) then
37342ce371bSBarry Smith    d1 = (right_v(1 + ym + 1) - right_v(1 + ym))*rhy
37442ce371bSBarry Smith    d2 = (top_v(1 + xm + 1) - top_v(1 + xm))*rhx
375c4762a1bSJed Brown    ft = ft + sqrt(1.0 + d1*d1 + d2*d2)
376c4762a1bSJed Brown  end if
377c4762a1bSJed Brown
378c4762a1bSJed Brown  ft = ft*area
379d8606c27SBarry Smith  PetscCallMPI(MPI_Allreduce(ft, fcn, 1, MPIU_SCALAR, MPIU_SUM, PETSC_COMM_WORLD, ierr))
380c4762a1bSJed Brown
381c4762a1bSJed Brown! Restore vectors
382ce78bad3SBarry Smith  PetscCall(VecRestoreArray(localX, x_v, ierr))
383ce78bad3SBarry Smith  PetscCall(VecRestoreArray(localV, g_v, ierr))
384ce78bad3SBarry Smith  PetscCall(VecRestoreArray(Left, left_v, ierr))
385ce78bad3SBarry Smith  PetscCall(VecRestoreArray(Top, top_v, ierr))
386ce78bad3SBarry Smith  PetscCall(VecRestoreArray(Bottom, bottom_v, ierr))
387ce78bad3SBarry Smith  PetscCall(VecRestoreArray(Right, right_v, ierr))
388c4762a1bSJed Brown
389c4762a1bSJed Brown! Scatter values to global vector
390d8606c27SBarry Smith  PetscCall(DMLocalToGlobalBegin(dm, localV, INSERT_VALUES, G, ierr))
391d8606c27SBarry Smith  PetscCall(DMLocalToGlobalEnd(dm, localV, INSERT_VALUES, G, ierr))
392c4762a1bSJed Brown
393d8606c27SBarry Smith  PetscCall(PetscLogFlops(70.0d0*xm*ym, ierr))
394c4762a1bSJed Brown
395c4762a1bSJed Brownend  !FormFunctionGradient
396c4762a1bSJed Brown
397c4762a1bSJed Brown! ----------------------------------------------------------------------------
398c4762a1bSJed Brown!
399c4762a1bSJed Brown!
400c4762a1bSJed Brown!   FormHessian - Evaluates Hessian matrix.
401c4762a1bSJed Brown!
402c4762a1bSJed Brown!   Input Parameters:
403c4762a1bSJed Brown!.  tao  - the Tao context
404c4762a1bSJed Brown!.  X    - input vector
405c4762a1bSJed Brown!.  dummy  - not used
406c4762a1bSJed Brown!
407c4762a1bSJed Brown!   Output Parameters:
408c4762a1bSJed Brown!.  Hessian    - Hessian matrix
4097addb90fSBarry Smith!.  Hpc    - optionally different matrix used to construct the preconditioner
410c4762a1bSJed Brown!
411c4762a1bSJed Brown!   Notes:
412c4762a1bSJed Brown!   Due to mesh point reordering with DMs, we must always work
413c4762a1bSJed Brown!   with the local mesh points, and then transform them to the new
414c4762a1bSJed Brown!   global numbering with the local-to-global mapping.  We cannot work
415c4762a1bSJed Brown!   directly with the global numbers for the original uniprocessor mesh!
416c4762a1bSJed Brown!
417c4762a1bSJed Brown!      MatSetValuesLocal(), using the local ordering (including
418c4762a1bSJed Brown!         ghost points!)
419c4762a1bSJed Brown!         - Then set matrix entries using the local ordering
420c4762a1bSJed Brown!           by calling MatSetValuesLocal()
421c4762a1bSJed Brown
422ce78bad3SBarry Smithsubroutine FormHessian(ta, X, Hessian, Hpc, dummy, ierr)
42377d968b7SBarry Smith  use plate2fmodule
424c4762a1bSJed Brown  implicit none
425c4762a1bSJed Brown
426ce78bad3SBarry Smith  Tao ta
427c4762a1bSJed Brown  Vec X
428c4762a1bSJed Brown  Mat Hessian, Hpc
429c4762a1bSJed Brown  PetscInt dummy
430c4762a1bSJed Brown  PetscErrorCode ierr
431c4762a1bSJed Brown
432c4762a1bSJed Brown  PetscInt i, j, k, row
433c4762a1bSJed Brown  PetscInt xs, xm, gxs, gxm
434c4762a1bSJed Brown  PetscInt ys, ym, gys, gym
435c4762a1bSJed Brown  PetscInt col(0:6)
436c4762a1bSJed Brown  PetscReal hx, hy, hydhx, hxdhy, rhx, rhy
437c4762a1bSJed Brown  PetscReal f1, f2, f3, f4, f5, f6, d1, d2, d3
438c4762a1bSJed Brown  PetscReal d4, d5, d6, d7, d8
439c4762a1bSJed Brown  PetscReal xc, xl, xr, xt, xb, xlt, xrb
440c4762a1bSJed Brown  PetscReal hl, hr, ht, hb, hc, htl, hbr
441c4762a1bSJed Brown
44242ce371bSBarry Smith  PetscReal, pointer ::  right_v(:), left_v(:)
44342ce371bSBarry Smith  PetscReal, pointer ::  bottom_v(:), top_v(:)
44442ce371bSBarry Smith  PetscReal, pointer ::  x_v(:)
445c4762a1bSJed Brown  PetscReal v(0:6)
446c4762a1bSJed Brown  PetscBool assembled
447c4762a1bSJed Brown  PetscInt i1
448c4762a1bSJed Brown
449c4762a1bSJed Brown  i1 = 1
450c4762a1bSJed Brown
451c4762a1bSJed Brown! Set various matrix options
452d8606c27SBarry Smith  PetscCall(MatSetOption(Hessian, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE, ierr))
453c4762a1bSJed Brown
454c4762a1bSJed Brown! Get local mesh boundaries
455d8606c27SBarry Smith  PetscCall(DMDAGetCorners(dm, xs, ys, PETSC_NULL_INTEGER, xm, ym, PETSC_NULL_INTEGER, ierr))
456d8606c27SBarry Smith  PetscCall(DMDAGetGhostCorners(dm, gxs, gys, PETSC_NULL_INTEGER, gxm, gym, PETSC_NULL_INTEGER, ierr))
457c4762a1bSJed Brown
458c4762a1bSJed Brown! Scatter ghost points to local vectors
459d8606c27SBarry Smith  PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX, ierr))
460d8606c27SBarry Smith  PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX, ierr))
461c4762a1bSJed Brown
462c4762a1bSJed Brown! Get pointers to vector data (see note on Fortran arrays above)
463ce78bad3SBarry Smith  PetscCall(VecGetArray(localX, x_v, ierr))
464ce78bad3SBarry Smith  PetscCall(VecGetArray(Top, top_v, ierr))
465ce78bad3SBarry Smith  PetscCall(VecGetArray(Bottom, bottom_v, ierr))
466ce78bad3SBarry Smith  PetscCall(VecGetArray(Left, left_v, ierr))
467ce78bad3SBarry Smith  PetscCall(VecGetArray(Right, right_v, ierr))
468c4762a1bSJed Brown
469c4762a1bSJed Brown! Initialize matrix entries to zero
470d8606c27SBarry Smith  PetscCall(MatAssembled(Hessian, assembled, ierr))
471d8606c27SBarry Smith  if (assembled .eqv. PETSC_TRUE) PetscCall(MatZeroEntries(Hessian, ierr))
472c4762a1bSJed Brown
473c4762a1bSJed Brown  rhx = real(mx) + 1.0
474c4762a1bSJed Brown  rhy = real(my) + 1.0
475c4762a1bSJed Brown  hx = 1.0/rhx
476c4762a1bSJed Brown  hy = 1.0/rhy
477c4762a1bSJed Brown  hydhx = hy/hx
478c4762a1bSJed Brown  hxdhy = hx/hy
479c4762a1bSJed Brown! compute Hessian over the locally owned part of the mesh
480c4762a1bSJed Brown
481c4762a1bSJed Brown  do i = xs, xs + xm - 1
482c4762a1bSJed Brown    do j = ys, ys + ym - 1
483c4762a1bSJed Brown      row = (j - gys)*gxm + (i - gxs)
484c4762a1bSJed Brown
48542ce371bSBarry Smith      xc = x_v(1 + row)
486c4762a1bSJed Brown      xt = xc
487c4762a1bSJed Brown      xb = xc
488c4762a1bSJed Brown      xr = xc
489c4762a1bSJed Brown      xl = xc
490c4762a1bSJed Brown      xrb = xc
491c4762a1bSJed Brown      xlt = xc
492c4762a1bSJed Brown
493*4820e4eaSBarry Smith      if (i == gxs) then   ! Left side
49442ce371bSBarry Smith        xl = left_v(1 + j - ys + 1)
49542ce371bSBarry Smith        xlt = left_v(1 + j - ys + 2)
496c4762a1bSJed Brown      else
49742ce371bSBarry Smith        xl = x_v(1 + row - 1)
498c4762a1bSJed Brown      end if
499c4762a1bSJed Brown
500*4820e4eaSBarry Smith      if (j == gys) then ! bottom side
50142ce371bSBarry Smith        xb = bottom_v(1 + i - xs + 1)
50242ce371bSBarry Smith        xrb = bottom_v(1 + i - xs + 2)
503c4762a1bSJed Brown      else
50442ce371bSBarry Smith        xb = x_v(1 + row - gxm)
505c4762a1bSJed Brown      end if
506c4762a1bSJed Brown
507*4820e4eaSBarry Smith      if (i + 1 == gxs + gxm) then !right side
50842ce371bSBarry Smith        xr = right_v(1 + j - ys + 1)
50942ce371bSBarry Smith        xrb = right_v(1 + j - ys)
510c4762a1bSJed Brown      else
51142ce371bSBarry Smith        xr = x_v(1 + row + 1)
512c4762a1bSJed Brown      end if
513c4762a1bSJed Brown
514*4820e4eaSBarry Smith      if (j + 1 == gym + gys) then !top side
51542ce371bSBarry Smith        xt = top_v(1 + i - xs + 1)
51642ce371bSBarry Smith        xlt = top_v(1 + i - xs)
517c4762a1bSJed Brown      else
51842ce371bSBarry Smith        xt = x_v(1 + row + gxm)
519c4762a1bSJed Brown      end if
520c4762a1bSJed Brown
521*4820e4eaSBarry Smith      if ((i > gxs) .and. (j + 1 < gys + gym)) then
52242ce371bSBarry Smith        xlt = x_v(1 + row - 1 + gxm)
523c4762a1bSJed Brown      end if
524c4762a1bSJed Brown
525*4820e4eaSBarry Smith      if ((i + 1 < gxs + gxm) .and. (j > gys)) then
52642ce371bSBarry Smith        xrb = x_v(1 + row + 1 - gxm)
527c4762a1bSJed Brown      end if
528c4762a1bSJed Brown
529c4762a1bSJed Brown      d1 = (xc - xl)*rhx
530c4762a1bSJed Brown      d2 = (xc - xr)*rhx
531c4762a1bSJed Brown      d3 = (xc - xt)*rhy
532c4762a1bSJed Brown      d4 = (xc - xb)*rhy
533c4762a1bSJed Brown      d5 = (xrb - xr)*rhy
534c4762a1bSJed Brown      d6 = (xrb - xb)*rhx
535c4762a1bSJed Brown      d7 = (xlt - xl)*rhy
536c4762a1bSJed Brown      d8 = (xlt - xt)*rhx
537c4762a1bSJed Brown
538c4762a1bSJed Brown      f1 = sqrt(1.0 + d1*d1 + d7*d7)
539c4762a1bSJed Brown      f2 = sqrt(1.0 + d1*d1 + d4*d4)
540c4762a1bSJed Brown      f3 = sqrt(1.0 + d3*d3 + d8*d8)
541c4762a1bSJed Brown      f4 = sqrt(1.0 + d3*d3 + d2*d2)
542c4762a1bSJed Brown      f5 = sqrt(1.0 + d2*d2 + d5*d5)
543c4762a1bSJed Brown      f6 = sqrt(1.0 + d4*d4 + d6*d6)
544c4762a1bSJed Brown
545d8606c27SBarry Smith      hl = (-hydhx*(1.0 + d7*d7) + d1*d7)/(f1*f1*f1) + (-hydhx*(1.0 + d4*d4) + d1*d4)/(f2*f2*f2)
546c4762a1bSJed Brown
547d8606c27SBarry Smith      hr = (-hydhx*(1.0 + d5*d5) + d2*d5)/(f5*f5*f5) + (-hydhx*(1.0 + d3*d3) + d2*d3)/(f4*f4*f4)
548c4762a1bSJed Brown
549d8606c27SBarry Smith      ht = (-hxdhy*(1.0 + d8*d8) + d3*d8)/(f3*f3*f3) + (-hxdhy*(1.0 + d2*d2) + d2*d3)/(f4*f4*f4)
550c4762a1bSJed Brown
551d8606c27SBarry Smith      hb = (-hxdhy*(1.0 + d6*d6) + d4*d6)/(f6*f6*f6) + (-hxdhy*(1.0 + d1*d1) + d1*d4)/(f2*f2*f2)
552c4762a1bSJed Brown
553c4762a1bSJed Brown      hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6)
554c4762a1bSJed Brown      htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3)
555c4762a1bSJed Brown
556d8606c27SBarry Smith      hc = hydhx*(1.0 + d7*d7)/(f1*f1*f1) + hxdhy*(1.0 + d8*d8)/(f3*f3*f3) + hydhx*(1.0 + d5*d5)/(f5*f5*f5) +                      &
557d8606c27SBarry Smith&           hxdhy*(1.0 + d6*d6)/(f6*f6*f6) + (hxdhy*(1.0 + d1*d1) + hydhx*(1.0 + d4*d4) - 2*d1*d4)/(f2*f2*f2) + (hxdhy*(1.0 + d2*d2) +   &
558c4762a1bSJed Brown&           hydhx*(1.0 + d3*d3) - 2*d2*d3)/(f4*f4*f4)
559c4762a1bSJed Brown
560c4762a1bSJed Brown      hl = hl*0.5
561c4762a1bSJed Brown      hr = hr*0.5
562c4762a1bSJed Brown      ht = ht*0.5
563c4762a1bSJed Brown      hb = hb*0.5
564c4762a1bSJed Brown      hbr = hbr*0.5
565c4762a1bSJed Brown      htl = htl*0.5
566c4762a1bSJed Brown      hc = hc*0.5
567c4762a1bSJed Brown
568c4762a1bSJed Brown      k = 0
569c4762a1bSJed Brown
570*4820e4eaSBarry Smith      if (j > 0) then
571c4762a1bSJed Brown        v(k) = hb
572c4762a1bSJed Brown        col(k) = row - gxm
573c4762a1bSJed Brown        k = k + 1
574c4762a1bSJed Brown      end if
575c4762a1bSJed Brown
576*4820e4eaSBarry Smith      if ((j > 0) .and. (i < mx - 1)) then
577c4762a1bSJed Brown        v(k) = hbr
578c4762a1bSJed Brown        col(k) = row - gxm + 1
579c4762a1bSJed Brown        k = k + 1
580c4762a1bSJed Brown      end if
581c4762a1bSJed Brown
582*4820e4eaSBarry Smith      if (i > 0) then
583c4762a1bSJed Brown        v(k) = hl
584c4762a1bSJed Brown        col(k) = row - 1
585c4762a1bSJed Brown        k = k + 1
586c4762a1bSJed Brown      end if
587c4762a1bSJed Brown
588c4762a1bSJed Brown      v(k) = hc
589c4762a1bSJed Brown      col(k) = row
590c4762a1bSJed Brown      k = k + 1
591c4762a1bSJed Brown
592*4820e4eaSBarry Smith      if (i < mx - 1) then
593c4762a1bSJed Brown        v(k) = hr
594c4762a1bSJed Brown        col(k) = row + 1
595c4762a1bSJed Brown        k = k + 1
596c4762a1bSJed Brown      end if
597c4762a1bSJed Brown
598*4820e4eaSBarry Smith      if ((i > 0) .and. (j < my - 1)) then
599c4762a1bSJed Brown        v(k) = htl
600c4762a1bSJed Brown        col(k) = row + gxm - 1
601c4762a1bSJed Brown        k = k + 1
602c4762a1bSJed Brown      end if
603c4762a1bSJed Brown
604*4820e4eaSBarry Smith      if (j < my - 1) then
605c4762a1bSJed Brown        v(k) = ht
606c4762a1bSJed Brown        col(k) = row + gxm
607c4762a1bSJed Brown        k = k + 1
608c4762a1bSJed Brown      end if
609c4762a1bSJed Brown
610c4762a1bSJed Brown! Set matrix values using local numbering, defined earlier in main routine
6115d83a8b1SBarry Smith      PetscCall(MatSetValuesLocal(Hessian, i1, [row], k, col, v, INSERT_VALUES, ierr))
612c4762a1bSJed Brown
613c4762a1bSJed Brown    end do
614c4762a1bSJed Brown  end do
615c4762a1bSJed Brown
616c4762a1bSJed Brown! restore vectors
617ce78bad3SBarry Smith  PetscCall(VecRestoreArray(localX, x_v, ierr))
618ce78bad3SBarry Smith  PetscCall(VecRestoreArray(Left, left_v, ierr))
619ce78bad3SBarry Smith  PetscCall(VecRestoreArray(Right, right_v, ierr))
620ce78bad3SBarry Smith  PetscCall(VecRestoreArray(Top, top_v, ierr))
621ce78bad3SBarry Smith  PetscCall(VecRestoreArray(Bottom, bottom_v, ierr))
622c4762a1bSJed Brown
623c4762a1bSJed Brown! Assemble the matrix
624d8606c27SBarry Smith  PetscCall(MatAssemblyBegin(Hessian, MAT_FINAL_ASSEMBLY, ierr))
625d8606c27SBarry Smith  PetscCall(MatAssemblyEnd(Hessian, MAT_FINAL_ASSEMBLY, ierr))
626c4762a1bSJed Brown
627d8606c27SBarry Smith  PetscCall(PetscLogFlops(199.0d0*xm*ym, ierr))
628c4762a1bSJed Brown
629c4762a1bSJed Brownend
630c4762a1bSJed Brown
631c4762a1bSJed Brown! Top,Left,Right,Bottom,bheight,mx,my,bmx,bmy,H, defined in plate2f.h
632c4762a1bSJed Brown
633c4762a1bSJed Brown! ----------------------------------------------------------------------------
634c4762a1bSJed Brown!
635c4762a1bSJed Brown!/*
636c4762a1bSJed Brown!     MSA_BoundaryConditions - calculates the boundary conditions for the region
637c4762a1bSJed Brown!
638c4762a1bSJed Brown!
639c4762a1bSJed Brown!*/
640c4762a1bSJed Brown
641c4762a1bSJed Brownsubroutine MSA_BoundaryConditions(ierr)
64277d968b7SBarry Smith  use plate2fmodule
643c4762a1bSJed Brown  implicit none
644c4762a1bSJed Brown
645c4762a1bSJed Brown  PetscErrorCode ierr
646c4762a1bSJed Brown  PetscInt i, j, k, limit, maxits
647c4762a1bSJed Brown  PetscInt xs, xm, gxs, gxm
648c4762a1bSJed Brown  PetscInt ys, ym, gys, gym
649c4762a1bSJed Brown  PetscInt bsize, lsize
650c4762a1bSJed Brown  PetscInt tsize, rsize
651c4762a1bSJed Brown  PetscReal one, two, three, tol
652c4762a1bSJed Brown  PetscReal scl, fnorm, det, xt
653c4762a1bSJed Brown  PetscReal yt, hx, hy, u1, u2, nf1, nf2
654c4762a1bSJed Brown  PetscReal njac11, njac12, njac21, njac22
655c4762a1bSJed Brown  PetscReal b, t, l, r
65642ce371bSBarry Smith  PetscReal, pointer :: boundary_v(:)
65742ce371bSBarry Smith
658c4762a1bSJed Brown  logical exitloop
659c4762a1bSJed Brown  PetscBool flg
660c4762a1bSJed Brown
661c4762a1bSJed Brown  limit = 0
662c4762a1bSJed Brown  maxits = 5
663c4762a1bSJed Brown  tol = 1e-10
664c4762a1bSJed Brown  b = -0.5
665c4762a1bSJed Brown  t = 0.5
666c4762a1bSJed Brown  l = -0.5
667c4762a1bSJed Brown  r = 0.5
668c4762a1bSJed Brown  xt = 0
669c4762a1bSJed Brown  yt = 0
670c4762a1bSJed Brown  one = 1.0
671c4762a1bSJed Brown  two = 2.0
672c4762a1bSJed Brown  three = 3.0
673c4762a1bSJed Brown
674d8606c27SBarry Smith  PetscCall(DMDAGetCorners(dm, xs, ys, PETSC_NULL_INTEGER, xm, ym, PETSC_NULL_INTEGER, ierr))
675d8606c27SBarry Smith  PetscCall(DMDAGetGhostCorners(dm, gxs, gys, PETSC_NULL_INTEGER, gxm, gym, PETSC_NULL_INTEGER, ierr))
676c4762a1bSJed Brown
677c4762a1bSJed Brown  bsize = xm + 2
678c4762a1bSJed Brown  lsize = ym + 2
679c4762a1bSJed Brown  rsize = ym + 2
680c4762a1bSJed Brown  tsize = xm + 2
681c4762a1bSJed Brown
682d8606c27SBarry Smith  PetscCall(VecCreateMPI(PETSC_COMM_WORLD, bsize, PETSC_DECIDE, Bottom, ierr))
683d8606c27SBarry Smith  PetscCall(VecCreateMPI(PETSC_COMM_WORLD, tsize, PETSC_DECIDE, Top, ierr))
684d8606c27SBarry Smith  PetscCall(VecCreateMPI(PETSC_COMM_WORLD, lsize, PETSC_DECIDE, Left, ierr))
685d8606c27SBarry Smith  PetscCall(VecCreateMPI(PETSC_COMM_WORLD, rsize, PETSC_DECIDE, Right, ierr))
686c4762a1bSJed Brown
687c4762a1bSJed Brown  hx = (r - l)/(mx + 1)
688c4762a1bSJed Brown  hy = (t - b)/(my + 1)
689c4762a1bSJed Brown
690c4762a1bSJed Brown  do j = 0, 3
691c4762a1bSJed Brown
692*4820e4eaSBarry Smith    if (j == 0) then
693c4762a1bSJed Brown      yt = b
694c4762a1bSJed Brown      xt = l + hx*xs
695c4762a1bSJed Brown      limit = bsize
696ce78bad3SBarry Smith      PetscCall(VecGetArray(Bottom, boundary_v, ierr))
697c4762a1bSJed Brown
698*4820e4eaSBarry Smith    elseif (j == 1) then
699c4762a1bSJed Brown      yt = t
700c4762a1bSJed Brown      xt = l + hx*xs
701c4762a1bSJed Brown      limit = tsize
702ce78bad3SBarry Smith      PetscCall(VecGetArray(Top, boundary_v, ierr))
703c4762a1bSJed Brown
704*4820e4eaSBarry Smith    elseif (j == 2) then
705c4762a1bSJed Brown      yt = b + hy*ys
706c4762a1bSJed Brown      xt = l
707c4762a1bSJed Brown      limit = lsize
708ce78bad3SBarry Smith      PetscCall(VecGetArray(Left, boundary_v, ierr))
709c4762a1bSJed Brown
710*4820e4eaSBarry Smith    elseif (j == 3) then
711c4762a1bSJed Brown      yt = b + hy*ys
712c4762a1bSJed Brown      xt = r
713c4762a1bSJed Brown      limit = rsize
714ce78bad3SBarry Smith      PetscCall(VecGetArray(Right, boundary_v, ierr))
715c4762a1bSJed Brown    end if
716c4762a1bSJed Brown
717c4762a1bSJed Brown    do i = 0, limit - 1
718c4762a1bSJed Brown
719c4762a1bSJed Brown      u1 = xt
720c4762a1bSJed Brown      u2 = -yt
721c4762a1bSJed Brown      k = 0
722c4762a1bSJed Brown      exitloop = .false.
723*4820e4eaSBarry Smith      do while (k < maxits .and. (.not. exitloop))
724c4762a1bSJed Brown
725c4762a1bSJed Brown        nf1 = u1 + u1*u2*u2 - u1*u1*u1/three - xt
726c4762a1bSJed Brown        nf2 = -u2 - u1*u1*u2 + u2*u2*u2/three - yt
727c4762a1bSJed Brown        fnorm = sqrt(nf1*nf1 + nf2*nf2)
728*4820e4eaSBarry Smith        if (fnorm > tol) then
729c4762a1bSJed Brown          njac11 = one + u2*u2 - u1*u1
730c4762a1bSJed Brown          njac12 = two*u1*u2
731c4762a1bSJed Brown          njac21 = -two*u1*u2
732c4762a1bSJed Brown          njac22 = -one - u1*u1 + u2*u2
733c4762a1bSJed Brown          det = njac11*njac22 - njac21*njac12
734c4762a1bSJed Brown          u1 = u1 - (njac22*nf1 - njac12*nf2)/det
735c4762a1bSJed Brown          u2 = u2 - (njac11*nf2 - njac21*nf1)/det
736c4762a1bSJed Brown        else
737c4762a1bSJed Brown          exitloop = .true.
738c4762a1bSJed Brown        end if
739c4762a1bSJed Brown        k = k + 1
740c4762a1bSJed Brown      end do
741c4762a1bSJed Brown
74242ce371bSBarry Smith      boundary_v(1 + i) = u1*u1 - u2*u2
743*4820e4eaSBarry Smith      if ((j == 0) .or. (j == 1)) then
744c4762a1bSJed Brown        xt = xt + hx
745c4762a1bSJed Brown      else
746c4762a1bSJed Brown        yt = yt + hy
747c4762a1bSJed Brown      end if
748c4762a1bSJed Brown
749c4762a1bSJed Brown    end do
750c4762a1bSJed Brown
751*4820e4eaSBarry Smith    if (j == 0) then
752ce78bad3SBarry Smith      PetscCall(VecRestoreArray(Bottom, boundary_v, ierr))
753*4820e4eaSBarry Smith    elseif (j == 1) then
754ce78bad3SBarry Smith      PetscCall(VecRestoreArray(Top, boundary_v, ierr))
755*4820e4eaSBarry Smith    elseif (j == 2) then
756ce78bad3SBarry Smith      PetscCall(VecRestoreArray(Left, boundary_v, ierr))
757*4820e4eaSBarry Smith    elseif (j == 3) then
758ce78bad3SBarry Smith      PetscCall(VecRestoreArray(Right, boundary_v, ierr))
759c4762a1bSJed Brown    end if
760c4762a1bSJed Brown
761c4762a1bSJed Brown  end do
762c4762a1bSJed Brown
763c4762a1bSJed Brown! Scale the boundary if desired
764d8606c27SBarry Smith  PetscCall(PetscOptionsGetReal(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-bottom', scl, flg, ierr))
765c4762a1bSJed Brown  if (flg .eqv. PETSC_TRUE) then
766d8606c27SBarry Smith    PetscCall(VecScale(Bottom, scl, ierr))
767c4762a1bSJed Brown  end if
768c4762a1bSJed Brown
769d8606c27SBarry Smith  PetscCall(PetscOptionsGetReal(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-top', scl, flg, ierr))
770c4762a1bSJed Brown  if (flg .eqv. PETSC_TRUE) then
771d8606c27SBarry Smith    PetscCall(VecScale(Top, scl, ierr))
772c4762a1bSJed Brown  end if
773c4762a1bSJed Brown
774d8606c27SBarry Smith  PetscCall(PetscOptionsGetReal(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-right', scl, flg, ierr))
775c4762a1bSJed Brown  if (flg .eqv. PETSC_TRUE) then
776d8606c27SBarry Smith    PetscCall(VecScale(Right, scl, ierr))
777c4762a1bSJed Brown  end if
778c4762a1bSJed Brown
779d8606c27SBarry Smith  PetscCall(PetscOptionsGetReal(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-left', scl, flg, ierr))
780c4762a1bSJed Brown  if (flg .eqv. PETSC_TRUE) then
781d8606c27SBarry Smith    PetscCall(VecScale(Left, scl, ierr))
782c4762a1bSJed Brown  end if
783c4762a1bSJed Brown
784c4762a1bSJed Brownend
785c4762a1bSJed Brown
786c4762a1bSJed Brown! ----------------------------------------------------------------------------
787c4762a1bSJed Brown!
788c4762a1bSJed Brown!/*
789c4762a1bSJed Brown!     MSA_Plate - Calculates an obstacle for surface to stretch over
790c4762a1bSJed Brown!
791c4762a1bSJed Brown!     Output Parameter:
792c4762a1bSJed Brown!.    xl - lower bound vector
793c4762a1bSJed Brown!.    xu - upper bound vector
794c4762a1bSJed Brown!
795c4762a1bSJed Brown!*/
796c4762a1bSJed Brown
797ce78bad3SBarry Smithsubroutine MSA_Plate(ta, xl, xu, dummy, ierr)
79877d968b7SBarry Smith  use plate2fmodule
799c4762a1bSJed Brown  implicit none
800c4762a1bSJed Brown
801ce78bad3SBarry Smith  Tao ta
802c4762a1bSJed Brown  Vec xl, xu
803c4762a1bSJed Brown  PetscErrorCode ierr
804c4762a1bSJed Brown  PetscInt i, j, row
805c4762a1bSJed Brown  PetscInt xs, xm, ys, ym
806c4762a1bSJed Brown  PetscReal lb, ub
807c4762a1bSJed Brown  PetscInt dummy
80842ce371bSBarry Smith  PetscReal, pointer :: xl_v(:)
809c4762a1bSJed Brown
810c4762a1bSJed Brown  lb = PETSC_NINFINITY
811c4762a1bSJed Brown  ub = PETSC_INFINITY
812c4762a1bSJed Brown
813*4820e4eaSBarry Smith  if (bmy < 0) bmy = 0
814*4820e4eaSBarry Smith  if (bmy > my) bmy = my
815*4820e4eaSBarry Smith  if (bmx < 0) bmx = 0
816*4820e4eaSBarry Smith  if (bmx > mx) bmx = mx
817c4762a1bSJed Brown
818d8606c27SBarry Smith  PetscCall(DMDAGetCorners(dm, xs, ys, PETSC_NULL_INTEGER, xm, ym, PETSC_NULL_INTEGER, ierr))
819c4762a1bSJed Brown
820d8606c27SBarry Smith  PetscCall(VecSet(xl, lb, ierr))
821d8606c27SBarry Smith  PetscCall(VecSet(xu, ub, ierr))
822c4762a1bSJed Brown
823ce78bad3SBarry Smith  PetscCall(VecGetArray(xl, xl_v, ierr))
824c4762a1bSJed Brown
825c4762a1bSJed Brown  do i = xs, xs + xm - 1
826c4762a1bSJed Brown
827c4762a1bSJed Brown    do j = ys, ys + ym - 1
828c4762a1bSJed Brown
829c4762a1bSJed Brown      row = (j - ys)*xm + (i - xs)
830c4762a1bSJed Brown
831*4820e4eaSBarry Smith      if (i >= ((mx - bmx)/2) .and. i < (mx - (mx - bmx)/2) .and.           &
832*4820e4eaSBarry Smith&          j >= ((my - bmy)/2) .and. j < (my - (my - bmy)/2)) then
83342ce371bSBarry Smith        xl_v(1 + row) = bheight
834c4762a1bSJed Brown
835c4762a1bSJed Brown      end if
836c4762a1bSJed Brown
837c4762a1bSJed Brown    end do
838c4762a1bSJed Brown  end do
839c4762a1bSJed Brown
840ce78bad3SBarry Smith  PetscCall(VecRestoreArray(xl, xl_v, ierr))
841c4762a1bSJed Brown
842c4762a1bSJed Brownend
843c4762a1bSJed Brown
844c4762a1bSJed Brown! ----------------------------------------------------------------------------
845c4762a1bSJed Brown!
846c4762a1bSJed Brown!/*
847c4762a1bSJed Brown!     MSA_InitialPoint - Calculates an obstacle for surface to stretch over
848c4762a1bSJed Brown!
849c4762a1bSJed Brown!     Output Parameter:
850c4762a1bSJed Brown!.    X - vector for initial guess
851c4762a1bSJed Brown!
852c4762a1bSJed Brown!*/
853c4762a1bSJed Brown
854c4762a1bSJed Brownsubroutine MSA_InitialPoint(X, ierr)
85577d968b7SBarry Smith  use plate2fmodule
856c4762a1bSJed Brown  implicit none
857c4762a1bSJed Brown
858c4762a1bSJed Brown  Vec X
859c4762a1bSJed Brown  PetscErrorCode ierr
860c4762a1bSJed Brown  PetscInt start, i, j
861c4762a1bSJed Brown  PetscInt row
862c4762a1bSJed Brown  PetscInt xs, xm, gxs, gxm
863c4762a1bSJed Brown  PetscInt ys, ym, gys, gym
864c4762a1bSJed Brown  PetscReal zero, np5
865c4762a1bSJed Brown
86642ce371bSBarry Smith  PetscReal, pointer :: left_v(:), right_v(:)
86742ce371bSBarry Smith  PetscReal, pointer :: bottom_v(:), top_v(:)
86842ce371bSBarry Smith  PetscReal, pointer :: x_v(:)
869c4762a1bSJed Brown  PetscBool flg
870c4762a1bSJed Brown  PetscRandom rctx
871c4762a1bSJed Brown
872c4762a1bSJed Brown  zero = 0.0
873c4762a1bSJed Brown  np5 = -0.5
874c4762a1bSJed Brown
875d8606c27SBarry Smith  PetscCall(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-start', start, flg, ierr))
876c4762a1bSJed Brown
877*4820e4eaSBarry Smith  if ((flg .eqv. PETSC_TRUE) .and. (start == 0)) then  ! the zero vector is reasonable
878d8606c27SBarry Smith    PetscCall(VecSet(X, zero, ierr))
879c4762a1bSJed Brown
880*4820e4eaSBarry Smith  elseif ((flg .eqv. PETSC_TRUE) .and. (start > 0)) then  ! random start -0.5 < xi < 0.5
881d8606c27SBarry Smith    PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, rctx, ierr))
882c4762a1bSJed Brown    do i = 0, start - 1
883d8606c27SBarry Smith      PetscCall(VecSetRandom(X, rctx, ierr))
884c4762a1bSJed Brown    end do
885c4762a1bSJed Brown
886d8606c27SBarry Smith    PetscCall(PetscRandomDestroy(rctx, ierr))
887d8606c27SBarry Smith    PetscCall(VecShift(X, np5, ierr))
888c4762a1bSJed Brown
889c4762a1bSJed Brown  else   ! average of boundary conditions
890c4762a1bSJed Brown
891c4762a1bSJed Brown!        Get Local mesh boundaries
892d8606c27SBarry Smith    PetscCall(DMDAGetCorners(dm, xs, ys, PETSC_NULL_INTEGER, xm, ym, PETSC_NULL_INTEGER, ierr))
893d8606c27SBarry Smith    PetscCall(DMDAGetGhostCorners(dm, gxs, gys, PETSC_NULL_INTEGER, gxm, gym, PETSC_NULL_INTEGER, ierr))
894c4762a1bSJed Brown
895c4762a1bSJed Brown!        Get pointers to vector data
896ce78bad3SBarry Smith    PetscCall(VecGetArray(Top, top_v, ierr))
897ce78bad3SBarry Smith    PetscCall(VecGetArray(Bottom, bottom_v, ierr))
898ce78bad3SBarry Smith    PetscCall(VecGetArray(Left, left_v, ierr))
899ce78bad3SBarry Smith    PetscCall(VecGetArray(Right, right_v, ierr))
900c4762a1bSJed Brown
901ce78bad3SBarry Smith    PetscCall(VecGetArray(localX, x_v, ierr))
902c4762a1bSJed Brown
903c4762a1bSJed Brown!        Perform local computations
904c4762a1bSJed Brown    do j = ys, ys + ym - 1
905c4762a1bSJed Brown      do i = xs, xs + xm - 1
906c4762a1bSJed Brown        row = (j - gys)*gxm + (i - gxs)
90742ce371bSBarry Smith        x_v(1 + row) = ((j + 1)*bottom_v(1 + i - xs + 1)/my + (my - j + 1)*top_v(1 + i - xs + 1)/(my + 2) +                  &
90842ce371bSBarry Smith&                          (i + 1)*left_v(1 + j - ys + 1)/mx + (mx - i + 1)*right_v(1 + j - ys + 1)/(mx + 2))*0.5
909c4762a1bSJed Brown      end do
910c4762a1bSJed Brown    end do
911c4762a1bSJed Brown
912c4762a1bSJed Brown!        Restore vectors
913ce78bad3SBarry Smith    PetscCall(VecRestoreArray(localX, x_v, ierr))
914c4762a1bSJed Brown
915ce78bad3SBarry Smith    PetscCall(VecRestoreArray(Left, left_v, ierr))
916ce78bad3SBarry Smith    PetscCall(VecRestoreArray(Top, top_v, ierr))
917ce78bad3SBarry Smith    PetscCall(VecRestoreArray(Bottom, bottom_v, ierr))
918ce78bad3SBarry Smith    PetscCall(VecRestoreArray(Right, right_v, ierr))
919c4762a1bSJed Brown
920d8606c27SBarry Smith    PetscCall(DMLocalToGlobalBegin(dm, localX, INSERT_VALUES, X, ierr))
921d8606c27SBarry Smith    PetscCall(DMLocalToGlobalEnd(dm, localX, INSERT_VALUES, X, ierr))
922c4762a1bSJed Brown
923c4762a1bSJed Brown  end if
924c4762a1bSJed Brown
925c4762a1bSJed Brownend
926c4762a1bSJed Brown
927c4762a1bSJed Brown!
928c4762a1bSJed Brown!/*TEST
929c4762a1bSJed Brown!
930c4762a1bSJed Brown!   build:
931c4762a1bSJed Brown!      requires: !complex
932c4762a1bSJed Brown!
933c4762a1bSJed Brown!   test:
93410978b7dSBarry Smith!      args: -tao_monitor_short -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bqnls -tao_gatol 1.e-4
935c4762a1bSJed Brown!      filter: sort -b
936c4762a1bSJed Brown!      filter_output: sort -b
937c4762a1bSJed Brown!      requires: !single
938c4762a1bSJed Brown!
939c4762a1bSJed Brown!   test:
940c4762a1bSJed Brown!      suffix: 2
941c4762a1bSJed Brown!      nsize: 2
94210978b7dSBarry Smith!      args: -tao_monitor_short -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bqnls -tao_gatol 1.e-4
943c4762a1bSJed Brown!      filter: sort -b
944c4762a1bSJed Brown!      filter_output: sort -b
945c4762a1bSJed Brown!      requires: !single
946c4762a1bSJed Brown!
947c4762a1bSJed Brown!TEST*/
948