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