xref: /petsc/src/tao/unconstrained/tutorials/rosenbrock1f.F90 (revision ccfd86f17c20321558100f6af55b03dc7cd752d2)
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"
168ba8ee0cSBarry Smith      use petsctao
178ba8ee0cSBarry 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
26392a88c0SBlaise Bourdin      type(tVec)      x       ! solution vector
27392a88c0SBlaise Bourdin      type(tMat)      H       ! hessian matrix
28ce78bad3SBarry Smith      type(tTao)      ta     ! TAO_SOVER context
29c4762a1bSJed Brown      PetscBool       flg
30c4762a1bSJed Brown      PetscInt        i2,i1
31c4762a1bSJed Brown      PetscMPIInt     size
32c4762a1bSJed Brown      PetscReal       zero
338ba8ee0cSBarry Smith      PetscReal       alpha
348ba8ee0cSBarry Smith      PetscInt        n
358ba8ee0cSBarry 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))
50dcb3e689SBarry Smith      PetscCheckA(size .eq. 1,PETSC_COMM_SELF,PETSC_ERR_WRONG_MPI_SIZE,'This is a uniprocessor example only')
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
63*ccfd86f1SBarry Smith!  Allocate storage space for Hessian
645d83a8b1SBarry Smith      PetscCallA(MatCreateSeqBAIJ(PETSC_COMM_SELF,i2,n,n,i1,PETSC_NULL_INTEGER_ARRAY, 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
71ce78bad3SBarry Smith      PetscCallA(TaoCreate(PETSC_COMM_SELF,ta,ierr))
72ce78bad3SBarry Smith      PetscCallA(TaoSetType(ta,TAOLMVM,ierr))
73c4762a1bSJed Brown
74c4762a1bSJed Brown!  Set routines for function, gradient, and hessian evaluation
75ce78bad3SBarry Smith      PetscCallA(TaoSetObjectiveAndGradient(ta,PETSC_NULL_VEC,FormFunctionGradient,0,ierr))
76ce78bad3SBarry Smith      PetscCallA(TaoSetHessian(ta,H,H,FormHessian,0,ierr))
77c4762a1bSJed Brown
78c4762a1bSJed Brown!  Optional: Set initial guess
79d8606c27SBarry Smith      PetscCallA(VecSet(x, zero, ierr))
80ce78bad3SBarry Smith      PetscCallA(TaoSetSolution(ta, x, ierr))
81c4762a1bSJed Brown
82c4762a1bSJed Brown!  Check for TAO command line options
83ce78bad3SBarry Smith      PetscCallA(TaoSetFromOptions(ta,ierr))
84c4762a1bSJed Brown
85c4762a1bSJed Brown!  SOLVE THE APPLICATION
86ce78bad3SBarry Smith      PetscCallA(TaoSolve(ta,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.
91ce78bad3SBarry Smith!      PetscCallA(TaoView(ta,PETSC_VIEWER_STDOUT_SELF,ierr))
92c4762a1bSJed Brown
93c4762a1bSJed Brown!  Free TAO data structures
94ce78bad3SBarry Smith      PetscCallA(TaoDestroy(ta,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
115ce78bad3SBarry Smith      subroutine FormFunctionGradient(ta, X, f, G, dummy, ierr)
1168ba8ee0cSBarry Smith      use petsctao
1178ba8ee0cSBarry Smith      implicit none
118c4762a1bSJed Brown
119ce78bad3SBarry Smith      type(tTao)       ta
120392a88c0SBlaise Bourdin      type(tVec)       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
12742ce371bSBarry Smith      PetscReal, pointer :: g_v(:),x_v(:)
1288ba8ee0cSBarry Smith      PetscReal        alpha
1298ba8ee0cSBarry Smith      PetscInt         n
1308ba8ee0cSBarry Smith      common /params/ alpha, n
131c4762a1bSJed Brown
132c4762a1bSJed Brown      ierr = 0
133c4762a1bSJed Brown      nn = n/2
134c4762a1bSJed Brown      ff = 0
135c4762a1bSJed Brown
136c4762a1bSJed Brown!     Get pointers to vector data
137ce78bad3SBarry Smith      PetscCall(VecGetArrayRead(X,x_v,ierr))
138ce78bad3SBarry Smith      PetscCall(VecGetArray(G,g_v,ierr))
139c4762a1bSJed Brown
140c4762a1bSJed Brown!     Compute G(X)
141c4762a1bSJed Brown      do i=0,nn-1
14242ce371bSBarry Smith         t1 = x_v(1+2*i+1) - x_v(1+2*i)*x_v(1+2*i)
14342ce371bSBarry Smith         t2 = 1.0 - x_v(1 + 2*i)
144c4762a1bSJed Brown         ff = ff + alpha*t1*t1 + t2*t2
14542ce371bSBarry Smith         g_v(1 + 2*i) = -4*alpha*t1*x_v(1 + 2*i) - 2.0*t2
14642ce371bSBarry Smith         g_v(1 + 2*i + 1) = 2.0*alpha*t1
147c4762a1bSJed Brown      enddo
148c4762a1bSJed Brown
149c4762a1bSJed Brown!     Restore vectors
150ce78bad3SBarry Smith      PetscCall(VecRestoreArrayRead(X,x_v,ierr))
151ce78bad3SBarry Smith      PetscCall(VecRestoreArray(G,g_v,ierr))
152c4762a1bSJed Brown
153c4762a1bSJed Brown      f = ff
154d8606c27SBarry Smith      PetscCall(PetscLogFlops(15.0d0*nn,ierr))
155c4762a1bSJed Brown
156c4762a1bSJed Brown      end
157c4762a1bSJed Brown
158c4762a1bSJed Brown!
159c4762a1bSJed Brown! ---------------------------------------------------------------------
160c4762a1bSJed Brown!
161c4762a1bSJed Brown!  FormHessian - Evaluates Hessian matrix.
162c4762a1bSJed Brown!
163c4762a1bSJed Brown!  Input Parameters:
164c4762a1bSJed Brown!  tao     - the Tao context
165c4762a1bSJed Brown!  X       - input vector
166c4762a1bSJed Brown!  dummy   - optional user-defined context, as set by SNESSetHessian()
167c4762a1bSJed Brown!            (not used here)
168c4762a1bSJed Brown!
169c4762a1bSJed Brown!  Output Parameters:
170c4762a1bSJed Brown!  H      - Hessian matrix
1717addb90fSBarry Smith!  PrecH  - optionally different matrix used to compute the preconditioner (not used here)
172c4762a1bSJed Brown!  ierr   - error code
173c4762a1bSJed Brown!
174c4762a1bSJed Brown!  Note: Providing the Hessian may not be necessary.  Only some solvers
175c4762a1bSJed Brown!  require this matrix.
176c4762a1bSJed Brown
177ce78bad3SBarry Smith      subroutine FormHessian(ta,X,H,PrecH,dummy,ierr)
1788ba8ee0cSBarry Smith      use petsctao
1798ba8ee0cSBarry Smith      implicit none
180c4762a1bSJed Brown
181c4762a1bSJed Brown!  Input/output variables:
182ce78bad3SBarry Smith      type(tTao)       ta
183392a88c0SBlaise Bourdin      type(tVec)       X
184392a88c0SBlaise Bourdin      type(tMat)       H, PrecH
185c4762a1bSJed Brown      PetscErrorCode   ierr
186c4762a1bSJed Brown      PetscInt         dummy
187c4762a1bSJed Brown
188c4762a1bSJed Brown      PetscReal        v(0:1,0:1)
189c4762a1bSJed Brown      PetscBool        assembled
190c4762a1bSJed Brown
191c4762a1bSJed Brown! PETSc's VecGetArray acts differently in Fortran than it does in C.
192c4762a1bSJed Brown! Calling VecGetArray((Vec) X, (PetscReal) x_array(0:1), (PetscOffset) x_index, ierr)
193c4762a1bSJed Brown! will return an array of doubles referenced by x_array offset by x_index.
194c4762a1bSJed Brown!  i.e.,  to reference the kth element of X, use x_array(k + x_index).
195c4762a1bSJed Brown! Notice that by declaring the arrays with range (0:1), we are using the C 0-indexing practice.
19642ce371bSBarry Smith      PetscReal, pointer :: x_v(:)
197c4762a1bSJed Brown      PetscInt         i,nn,ind(0:1),i2
1988ba8ee0cSBarry Smith      PetscReal        alpha
1998ba8ee0cSBarry Smith      PetscInt         n
2008ba8ee0cSBarry Smith      common /params/ alpha, n
201c4762a1bSJed Brown
202c4762a1bSJed Brown      ierr = 0
203c4762a1bSJed Brown      nn= n/2
204c4762a1bSJed Brown      i2 = 2
205c4762a1bSJed Brown
206c4762a1bSJed Brown!  Zero existing matrix entries
207d8606c27SBarry Smith      PetscCall(MatAssembled(H,assembled,ierr))
208d8606c27SBarry Smith      if (assembled .eqv. PETSC_TRUE) PetscCall(MatZeroEntries(H,ierr))
209c4762a1bSJed Brown
210c4762a1bSJed Brown!  Get a pointer to vector data
211c4762a1bSJed Brown
212ce78bad3SBarry Smith      PetscCall(VecGetArrayRead(X,x_v,ierr))
213c4762a1bSJed Brown
214c4762a1bSJed Brown!  Compute Hessian entries
215c4762a1bSJed Brown
216c4762a1bSJed Brown      do i=0,nn-1
217c4762a1bSJed Brown         v(1,1) = 2.0*alpha
21842ce371bSBarry Smith         v(0,0) = -4.0*alpha*(x_v(1+2*i+1) - 3*x_v(1+2*i)*x_v(1+2*i))+2
21942ce371bSBarry Smith         v(1,0) = -4.0*alpha*x_v(1+2*i)
220c4762a1bSJed Brown         v(0,1) = v(1,0)
221c4762a1bSJed Brown         ind(0) = 2*i
222c4762a1bSJed Brown         ind(1) = 2*i + 1
223ce78bad3SBarry Smith         PetscCall(MatSetValues(H,i2,ind,i2,ind,reshape(v,[i2*i2]),INSERT_VALUES,ierr))
224c4762a1bSJed Brown      enddo
225c4762a1bSJed Brown
226c4762a1bSJed Brown!  Restore vector
227c4762a1bSJed Brown
228ce78bad3SBarry Smith      PetscCall(VecRestoreArrayRead(X,x_v,ierr))
229c4762a1bSJed Brown
230c4762a1bSJed Brown!  Assemble matrix
231c4762a1bSJed Brown
232d8606c27SBarry Smith      PetscCall(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY,ierr))
233d8606c27SBarry Smith      PetscCall(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY,ierr))
234c4762a1bSJed Brown
235d8606c27SBarry Smith      PetscCall(PetscLogFlops(9.0d0*nn,ierr))
236c4762a1bSJed Brown
237c4762a1bSJed Brown      end
238c4762a1bSJed Brown
239c4762a1bSJed Brown!
240c4762a1bSJed Brown!/*TEST
241c4762a1bSJed Brown!
242c4762a1bSJed Brown!   build:
243c4762a1bSJed Brown!      requires: !complex
244c4762a1bSJed Brown!
245c4762a1bSJed Brown!   test:
24610978b7dSBarry Smith!      args: -tao_monitor_short -tao_type ntr -tao_gatol 1.e-5
247c4762a1bSJed Brown!      requires: !single
248c4762a1bSJed Brown!
249c4762a1bSJed Brown!TEST*/
250