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