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