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