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