1c4762a1bSJed Brown! Program usage: mpiexec -n 1 rosenbrock1f [-help] [all TAO options] 2c4762a1bSJed Brown! 3c4762a1bSJed Brown! Description: This example demonstrates use of the TAO package to solve an 4c4762a1bSJed Brown! unconstrained minimization problem on a single processor. We minimize the 5c4762a1bSJed Brown! extended Rosenbrock function: 6c4762a1bSJed Brown! sum_{i=0}^{n/2-1} (alpha*(x_{2i+1}-x_{2i}^2)^2 + (1-x_{2i})^2) 7c4762a1bSJed Brown! 8c4762a1bSJed Brown! The C version of this code is rosenbrock1.c 9c4762a1bSJed Brown! 10c4762a1bSJed Brown 11c4762a1bSJed Brown! 12c4762a1bSJed Brown 13c4762a1bSJed Brown! ---------------------------------------------------------------------- 14c4762a1bSJed Brown! 158ba8ee0cSBarry Smith#include "petsc/finclude/petsctao.h" 16*e7a95102SMartin Diehlmodule rosenbrock1f_mod 178ba8ee0cSBarry Smith use petsctao 188ba8ee0cSBarry Smith implicit none 19*e7a95102SMartin Diehlcontains 20*e7a95102SMartin Diehl! -------------------------------------------------------------------- 21*e7a95102SMartin Diehl! FormFunctionGradient - Evaluates the function f(X) and gradient G(X) 22*e7a95102SMartin Diehl! 23*e7a95102SMartin Diehl! Input Parameters: 24*e7a95102SMartin Diehl! tao - the Tao context 25*e7a95102SMartin Diehl! X - input vector 26*e7a95102SMartin Diehl! dummy - not used 27*e7a95102SMartin Diehl! 28*e7a95102SMartin Diehl! Output Parameters: 29*e7a95102SMartin Diehl! G - vector containing the newly evaluated gradient 30*e7a95102SMartin Diehl! f - function value 31*e7a95102SMartin Diehl 32*e7a95102SMartin Diehl subroutine FormFunctionGradient(ta, X, f, G, dummy, ierr) 33*e7a95102SMartin Diehl type(tTao) ta 34*e7a95102SMartin Diehl type(tVec) X, G 35*e7a95102SMartin Diehl PetscReal f 36*e7a95102SMartin Diehl PetscErrorCode ierr 37*e7a95102SMartin Diehl PetscInt dummy 38*e7a95102SMartin Diehl 39*e7a95102SMartin Diehl PetscReal ff, t1, t2 40*e7a95102SMartin Diehl PetscInt i, nn 41*e7a95102SMartin Diehl PetscReal, pointer :: g_v(:), x_v(:) 42*e7a95102SMartin Diehl PetscReal alpha 43*e7a95102SMartin Diehl PetscInt n 44*e7a95102SMartin Diehl common/params/alpha, n 45*e7a95102SMartin Diehl 46*e7a95102SMartin Diehl ierr = 0 47*e7a95102SMartin Diehl nn = n/2 48*e7a95102SMartin Diehl ff = 0 49*e7a95102SMartin Diehl 50*e7a95102SMartin Diehl! Get pointers to vector data 51*e7a95102SMartin Diehl PetscCall(VecGetArrayRead(X, x_v, ierr)) 52*e7a95102SMartin Diehl PetscCall(VecGetArray(G, g_v, ierr)) 53*e7a95102SMartin Diehl 54*e7a95102SMartin Diehl! Compute G(X) 55*e7a95102SMartin Diehl do i = 0, nn - 1 56*e7a95102SMartin Diehl t1 = x_v(1 + 2*i + 1) - x_v(1 + 2*i)*x_v(1 + 2*i) 57*e7a95102SMartin Diehl t2 = 1.0 - x_v(1 + 2*i) 58*e7a95102SMartin Diehl ff = ff + alpha*t1*t1 + t2*t2 59*e7a95102SMartin Diehl g_v(1 + 2*i) = -4*alpha*t1*x_v(1 + 2*i) - 2.0*t2 60*e7a95102SMartin Diehl g_v(1 + 2*i + 1) = 2.0*alpha*t1 61*e7a95102SMartin Diehl end do 62*e7a95102SMartin Diehl 63*e7a95102SMartin Diehl! Restore vectors 64*e7a95102SMartin Diehl PetscCall(VecRestoreArrayRead(X, x_v, ierr)) 65*e7a95102SMartin Diehl PetscCall(VecRestoreArray(G, g_v, ierr)) 66*e7a95102SMartin Diehl 67*e7a95102SMartin Diehl f = ff 68*e7a95102SMartin Diehl PetscCall(PetscLogFlops(15.0d0*nn, ierr)) 69*e7a95102SMartin Diehl 70*e7a95102SMartin Diehl end 71*e7a95102SMartin Diehl 72*e7a95102SMartin Diehl! 73*e7a95102SMartin Diehl! --------------------------------------------------------------------- 74*e7a95102SMartin Diehl! 75*e7a95102SMartin Diehl! FormHessian - Evaluates Hessian matrix. 76*e7a95102SMartin Diehl! 77*e7a95102SMartin Diehl! Input Parameters: 78*e7a95102SMartin Diehl! tao - the Tao context 79*e7a95102SMartin Diehl! X - input vector 80*e7a95102SMartin Diehl! dummy - optional user-defined context, as set by SNESSetHessian() 81*e7a95102SMartin Diehl! (not used here) 82*e7a95102SMartin Diehl! 83*e7a95102SMartin Diehl! Output Parameters: 84*e7a95102SMartin Diehl! H - Hessian matrix 85*e7a95102SMartin Diehl! PrecH - optionally different matrix used to compute the preconditioner (not used here) 86*e7a95102SMartin Diehl! ierr - error code 87*e7a95102SMartin Diehl! 88*e7a95102SMartin Diehl! Note: Providing the Hessian may not be necessary. Only some solvers 89*e7a95102SMartin Diehl! require this matrix. 90*e7a95102SMartin Diehl 91*e7a95102SMartin Diehl subroutine FormHessian(ta, X, H, PrecH, dummy, ierr) 92*e7a95102SMartin Diehl 93*e7a95102SMartin Diehl! Input/output variables: 94*e7a95102SMartin Diehl type(tTao) ta 95*e7a95102SMartin Diehl type(tVec) X 96*e7a95102SMartin Diehl type(tMat) H, PrecH 97*e7a95102SMartin Diehl PetscErrorCode ierr 98*e7a95102SMartin Diehl PetscInt dummy 99*e7a95102SMartin Diehl 100*e7a95102SMartin Diehl PetscReal v(0:1, 0:1) 101*e7a95102SMartin Diehl PetscBool assembled 102*e7a95102SMartin Diehl 103*e7a95102SMartin Diehl! PETSc's VecGetArray acts differently in Fortran than it does in C. 104*e7a95102SMartin Diehl! Calling VecGetArray((Vec) X, (PetscReal) x_array(0:1), (PetscOffset) x_index, ierr) 105*e7a95102SMartin Diehl! will return an array of doubles referenced by x_array offset by x_index. 106*e7a95102SMartin Diehl! i.e., to reference the kth element of X, use x_array(k + x_index). 107*e7a95102SMartin Diehl! Notice that by declaring the arrays with range (0:1), we are using the C 0-indexing practice. 108*e7a95102SMartin Diehl PetscReal, pointer :: x_v(:) 109*e7a95102SMartin Diehl PetscInt i, nn, ind(0:1), i2 110*e7a95102SMartin Diehl PetscReal alpha 111*e7a95102SMartin Diehl PetscInt n 112*e7a95102SMartin Diehl common/params/alpha, n 113*e7a95102SMartin Diehl 114*e7a95102SMartin Diehl ierr = 0 115*e7a95102SMartin Diehl nn = n/2 116*e7a95102SMartin Diehl i2 = 2 117*e7a95102SMartin Diehl 118*e7a95102SMartin Diehl! Zero existing matrix entries 119*e7a95102SMartin Diehl PetscCall(MatAssembled(H, assembled, ierr)) 120*e7a95102SMartin Diehl if (assembled .eqv. PETSC_TRUE) PetscCall(MatZeroEntries(H, ierr)) 121*e7a95102SMartin Diehl 122*e7a95102SMartin Diehl! Get a pointer to vector data 123*e7a95102SMartin Diehl 124*e7a95102SMartin Diehl PetscCall(VecGetArrayRead(X, x_v, ierr)) 125*e7a95102SMartin Diehl 126*e7a95102SMartin Diehl! Compute Hessian entries 127*e7a95102SMartin Diehl 128*e7a95102SMartin Diehl do i = 0, nn - 1 129*e7a95102SMartin Diehl v(1, 1) = 2.0*alpha 130*e7a95102SMartin Diehl v(0, 0) = -4.0*alpha*(x_v(1 + 2*i + 1) - 3*x_v(1 + 2*i)*x_v(1 + 2*i)) + 2 131*e7a95102SMartin Diehl v(1, 0) = -4.0*alpha*x_v(1 + 2*i) 132*e7a95102SMartin Diehl v(0, 1) = v(1, 0) 133*e7a95102SMartin Diehl ind(0) = 2*i 134*e7a95102SMartin Diehl ind(1) = 2*i + 1 135*e7a95102SMartin Diehl PetscCall(MatSetValues(H, i2, ind, i2, ind, reshape(v, [i2*i2]), INSERT_VALUES, ierr)) 136*e7a95102SMartin Diehl end do 137*e7a95102SMartin Diehl 138*e7a95102SMartin Diehl! Restore vector 139*e7a95102SMartin Diehl 140*e7a95102SMartin Diehl PetscCall(VecRestoreArrayRead(X, x_v, ierr)) 141*e7a95102SMartin Diehl 142*e7a95102SMartin Diehl! Assemble matrix 143*e7a95102SMartin Diehl 144*e7a95102SMartin Diehl PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY, ierr)) 145*e7a95102SMartin Diehl PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY, ierr)) 146*e7a95102SMartin Diehl 147*e7a95102SMartin Diehl PetscCall(PetscLogFlops(9.0d0*nn, ierr)) 148*e7a95102SMartin Diehl 149*e7a95102SMartin Diehl end 150*e7a95102SMartin Diehlend module rosenbrock1f_mod 151*e7a95102SMartin Diehl 152*e7a95102SMartin Diehlprogram rosenbrock1f 153*e7a95102SMartin Diehl use petsctao 154*e7a95102SMartin Diehl use rosenbrock1f_mod 155*e7a95102SMartin Diehl implicit none 156c4762a1bSJed Brown 157c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 158c4762a1bSJed Brown! Variable declarations 159c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 160c4762a1bSJed Brown! 161c4762a1bSJed Brown! See additional variable declarations in the file rosenbrock1f.h 162c4762a1bSJed Brown 163c4762a1bSJed Brown PetscErrorCode ierr ! used to check for functions returning nonzeros 164392a88c0SBlaise Bourdin type(tVec) x ! solution vector 165392a88c0SBlaise Bourdin type(tMat) H ! hessian matrix 166ce78bad3SBarry Smith type(tTao) ta ! TAO_SOVER context 167c4762a1bSJed Brown PetscBool flg 168c4762a1bSJed Brown PetscInt i2, i1 169c4762a1bSJed Brown PetscMPIInt size 170c4762a1bSJed Brown PetscReal zero 1718ba8ee0cSBarry Smith PetscReal alpha 1728ba8ee0cSBarry Smith PetscInt n 1738ba8ee0cSBarry Smith common/params/alpha, n 174c4762a1bSJed Brown 175c4762a1bSJed Brown zero = 0.0d0 176c4762a1bSJed Brown i2 = 2 177c4762a1bSJed Brown i1 = 1 178c4762a1bSJed Brown 179c4762a1bSJed Brown! Initialize TAO and PETSc 180d8606c27SBarry Smith PetscCallA(PetscInitialize(ierr)) 181c4762a1bSJed Brown 182d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, size, ierr)) 1834820e4eaSBarry Smith PetscCheckA(size == 1, PETSC_COMM_SELF, PETSC_ERR_WRONG_MPI_SIZE, 'This is a uniprocessor example only') 184c4762a1bSJed Brown 185c4762a1bSJed Brown! Initialize problem parameters 186c4762a1bSJed Brown n = 2 187c4762a1bSJed Brown alpha = 99.0d0 188c4762a1bSJed Brown 189c4762a1bSJed Brown! Check for command line arguments to override defaults 190d8606c27SBarry Smith PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-n', n, flg, ierr)) 191d8606c27SBarry Smith PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-alpha', alpha, flg, ierr)) 192c4762a1bSJed Brown 193c4762a1bSJed Brown! Allocate vectors for the solution and gradient 194d8606c27SBarry Smith PetscCallA(VecCreateSeq(PETSC_COMM_SELF, n, x, ierr)) 195c4762a1bSJed Brown 196ccfd86f1SBarry Smith! Allocate storage space for Hessian 1975d83a8b1SBarry Smith PetscCallA(MatCreateSeqBAIJ(PETSC_COMM_SELF, i2, n, n, i1, PETSC_NULL_INTEGER_ARRAY, H, ierr)) 198c4762a1bSJed Brown 199d8606c27SBarry Smith PetscCallA(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE, ierr)) 200c4762a1bSJed Brown 201c4762a1bSJed Brown! The TAO code begins here 202c4762a1bSJed Brown 203c4762a1bSJed Brown! Create TAO solver 204ce78bad3SBarry Smith PetscCallA(TaoCreate(PETSC_COMM_SELF, ta, ierr)) 205ce78bad3SBarry Smith PetscCallA(TaoSetType(ta, TAOLMVM, ierr)) 206c4762a1bSJed Brown 207c4762a1bSJed Brown! Set routines for function, gradient, and hessian evaluation 208ce78bad3SBarry Smith PetscCallA(TaoSetObjectiveAndGradient(ta, PETSC_NULL_VEC, FormFunctionGradient, 0, ierr)) 209ce78bad3SBarry Smith PetscCallA(TaoSetHessian(ta, H, H, FormHessian, 0, ierr)) 210c4762a1bSJed Brown 211c4762a1bSJed Brown! Optional: Set initial guess 212d8606c27SBarry Smith PetscCallA(VecSet(x, zero, ierr)) 213ce78bad3SBarry Smith PetscCallA(TaoSetSolution(ta, x, ierr)) 214c4762a1bSJed Brown 215c4762a1bSJed Brown! Check for TAO command line options 216ce78bad3SBarry Smith PetscCallA(TaoSetFromOptions(ta, ierr)) 217c4762a1bSJed Brown 218c4762a1bSJed Brown! SOLVE THE APPLICATION 219ce78bad3SBarry Smith PetscCallA(TaoSolve(ta, ierr)) 220c4762a1bSJed Brown 221c4762a1bSJed Brown! TaoView() prints ierr about the TAO solver; the option 222c4762a1bSJed Brown! -tao_view 223c4762a1bSJed Brown! can alternatively be used to activate this at runtime. 224ce78bad3SBarry Smith! PetscCallA(TaoView(ta,PETSC_VIEWER_STDOUT_SELF,ierr)) 225c4762a1bSJed Brown 226c4762a1bSJed Brown! Free TAO data structures 227ce78bad3SBarry Smith PetscCallA(TaoDestroy(ta, ierr)) 228c4762a1bSJed Brown 229c4762a1bSJed Brown! Free PETSc data structures 230d8606c27SBarry Smith PetscCallA(VecDestroy(x, ierr)) 231d8606c27SBarry Smith PetscCallA(MatDestroy(H, ierr)) 232c4762a1bSJed Brown 233d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr)) 234*e7a95102SMartin Diehlend program rosenbrock1f 235c4762a1bSJed Brown! 236c4762a1bSJed Brown!/*TEST 237c4762a1bSJed Brown! 238c4762a1bSJed Brown! build: 239c4762a1bSJed Brown! requires: !complex 240c4762a1bSJed Brown! 241c4762a1bSJed Brown! test: 24210978b7dSBarry Smith! args: -tao_monitor_short -tao_type ntr -tao_gatol 1.e-5 243c4762a1bSJed Brown! requires: !single 244c4762a1bSJed Brown! 245c4762a1bSJed Brown!TEST*/ 246