xref: /petsc/src/tao/leastsquares/tutorials/chwirut1f.F90 (revision dfbbaf821b4c49d07b4ce746493b0d955783fdf9)
1c4762a1bSJed Brown!  Program usage: mpiexec -n 1 chwirut1f [-help] [all TAO options]
2c4762a1bSJed Brown!
3c4762a1bSJed Brown!  Description:  This example demonstrates use of the TAO package to solve a
4c4762a1bSJed Brown!  nonlinear least-squares problem on a single processor.  We minimize the
5c4762a1bSJed Brown!  Chwirut 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 test_chwirut1.c
9c4762a1bSJed Brown!
10c4762a1bSJed Brown
11c4762a1bSJed Brown!
12c4762a1bSJed Brown! ----------------------------------------------------------------------
13c4762a1bSJed Brown!
14*dfbbaf82SBarry Smith      module chwirut1fmodule
15*dfbbaf82SBarry Smith      use petsctao
16*dfbbaf82SBarry Smith#include <petsc/finclude/petsctao.h>
17*dfbbaf82SBarry Smith      PetscReal t(0:213)
18*dfbbaf82SBarry Smith      PetscReal y(0:213)
19*dfbbaf82SBarry Smith      PetscInt  m,n
20*dfbbaf82SBarry Smith      end module chwirut1fmodule
21c4762a1bSJed Brown
22*dfbbaf82SBarry Smith      program main
23*dfbbaf82SBarry Smith      use chwirut1fmodule
24c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
25c4762a1bSJed Brown!                   Variable declarations
26c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
27c4762a1bSJed Brown!
28c4762a1bSJed Brown!  See additional variable declarations in the file chwirut1f.h
29c4762a1bSJed Brown
30c4762a1bSJed Brown      PetscErrorCode   ierr    ! used to check for functions returning nonzeros
31c4762a1bSJed Brown      Vec              x       ! solution vector
32c4762a1bSJed Brown      Vec              f       ! vector of functions
33c4762a1bSJed Brown      Tao        tao     ! Tao context
34c4762a1bSJed Brown      PetscInt         nhist
35c4762a1bSJed Brown      PetscMPIInt  size,rank    ! number of processes running
36c4762a1bSJed Brown      PetscReal      hist(100) ! objective value history
37c4762a1bSJed Brown      PetscReal      resid(100)! residual history
38c4762a1bSJed Brown      PetscReal      cnorm(100)! cnorm history
39c4762a1bSJed Brown      PetscInt      lits(100)   ! #ksp history
40c4762a1bSJed Brown      PetscInt oh
41c4762a1bSJed Brown      TaoConvergedReason reason
42c4762a1bSJed Brown
43c4762a1bSJed Brown!  Note: Any user-defined Fortran routines (such as FormGradient)
44c4762a1bSJed Brown!  MUST be declared as external.
45c4762a1bSJed Brown
46c4762a1bSJed Brown      external FormFunction
47c4762a1bSJed Brown
48c4762a1bSJed Brown!  Initialize TAO and PETSc
49d8606c27SBarry Smith      PetscCallA(PetscInitialize(ierr))
50c4762a1bSJed Brown
51d8606c27SBarry Smith      PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD,size,ierr))
52d8606c27SBarry Smith      PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
53c4762a1bSJed Brown      if (size .ne. 1) then; SETERRA(PETSC_COMM_SELF,PETSC_ERR_WRONG_MPI_SIZE,'This is a uniprocessor example only '); endif
54c4762a1bSJed Brown
55c4762a1bSJed Brown!  Initialize problem parameters
56c4762a1bSJed Brown      m = 214
57c4762a1bSJed Brown      n = 3
58c4762a1bSJed Brown
59c4762a1bSJed Brown!  Allocate vectors for the solution and gradient
60d8606c27SBarry Smith      PetscCallA(VecCreateSeq(PETSC_COMM_SELF,n,x,ierr))
61d8606c27SBarry Smith      PetscCallA(VecCreateSeq(PETSC_COMM_SELF,m,f,ierr))
62c4762a1bSJed Brown
63c4762a1bSJed Brown!  The TAO code begins here
64c4762a1bSJed Brown
65c4762a1bSJed Brown!  Create TAO solver
66d8606c27SBarry Smith      PetscCallA(TaoCreate(PETSC_COMM_SELF,tao,ierr))
67d8606c27SBarry Smith      PetscCallA(TaoSetType(tao,TAOPOUNDERS,ierr))
68c4762a1bSJed Brown!  Set routines for function, gradient, and hessian evaluation
69c4762a1bSJed Brown
70d8606c27SBarry Smith      PetscCallA(TaoSetResidualRoutine(tao,f,FormFunction,0,ierr))
71c4762a1bSJed Brown
72c4762a1bSJed Brown!  Optional: Set initial guess
73c4762a1bSJed Brown      call InitializeData()
74c4762a1bSJed Brown      call FormStartingPoint(x)
75d8606c27SBarry Smith      PetscCallA(TaoSetSolution(tao, x, ierr))
76c4762a1bSJed Brown
77c4762a1bSJed Brown!  Check for TAO command line options
78d8606c27SBarry Smith      PetscCallA(TaoSetFromOptions(tao,ierr))
79c4762a1bSJed Brown      oh = 100
80d8606c27SBarry Smith      PetscCallA(TaoSetConvergenceHistory(tao,hist,resid,cnorm,lits,oh,PETSC_TRUE,ierr))
81c4762a1bSJed Brown!  SOLVE THE APPLICATION
82d8606c27SBarry Smith      PetscCallA(TaoSolve(tao,ierr))
83d8606c27SBarry Smith      PetscCallA(TaoGetConvergenceHistory(tao,nhist,ierr))
84d8606c27SBarry Smith      PetscCallA(TaoGetConvergedReason(tao, reason, ierr))
85c4762a1bSJed Brown      if (reason .le. 0) then
86c4762a1bSJed Brown         print *,'Tao failed.'
87c4762a1bSJed Brown         print *,'Try a different TAO method, adjust some parameters,'
88c4762a1bSJed Brown         print *,'or check the function evaluation routines.'
89c4762a1bSJed Brown      endif
90c4762a1bSJed Brown
91c4762a1bSJed Brown!  Free TAO data structures
92d8606c27SBarry Smith      PetscCallA(TaoDestroy(tao,ierr))
93c4762a1bSJed Brown
94c4762a1bSJed Brown!  Free PETSc data structures
95d8606c27SBarry Smith      PetscCallA(VecDestroy(x,ierr))
96d8606c27SBarry Smith      PetscCallA(VecDestroy(f,ierr))
97c4762a1bSJed Brown
98d8606c27SBarry Smith      PetscCallA(PetscFinalize(ierr))
99c4762a1bSJed Brown
100c4762a1bSJed Brown      end
101c4762a1bSJed Brown
102c4762a1bSJed Brown! --------------------------------------------------------------------
103c4762a1bSJed Brown!  FormFunction - Evaluates the function f(X) and gradient G(X)
104c4762a1bSJed Brown!
105c4762a1bSJed Brown!  Input Parameters:
106c4762a1bSJed Brown!  tao - the Tao context
107c4762a1bSJed Brown!  X   - input vector
108c4762a1bSJed Brown!  dummy - not used
109c4762a1bSJed Brown!
110c4762a1bSJed Brown!  Output Parameters:
111c4762a1bSJed Brown!  f - function vector
112c4762a1bSJed Brown
113c4762a1bSJed Brown      subroutine FormFunction(tao, x, f, dummy, ierr)
114*dfbbaf82SBarry Smith      use chwirut1fmodule
115c4762a1bSJed Brown
116c4762a1bSJed Brown      Tao        tao
117c4762a1bSJed Brown      Vec              x,f
118c4762a1bSJed Brown      PetscErrorCode   ierr
119c4762a1bSJed Brown      PetscInt         dummy
120c4762a1bSJed Brown
121c4762a1bSJed Brown      PetscInt         i
122a00960a2SBarry Smith      PetscScalar, pointer, dimension(:)  :: x_v,f_v
123c4762a1bSJed Brown
124c4762a1bSJed Brown      ierr = 0
125c4762a1bSJed Brown
126c4762a1bSJed Brown!     Get pointers to vector data
127d8606c27SBarry Smith      PetscCall(VecGetArrayF90(x,x_v,ierr))
128d8606c27SBarry Smith      PetscCall(VecGetArrayF90(f,f_v,ierr))
129c4762a1bSJed Brown
130c4762a1bSJed Brown!     Compute F(X)
131c4762a1bSJed Brown      do i=0,m-1
132a00960a2SBarry Smith         f_v(i+1) = y(i) - exp(-x_v(1)*t(i))/(x_v(2) + x_v(3)*t(i))
133c4762a1bSJed Brown      enddo
134c4762a1bSJed Brown
135c4762a1bSJed Brown!     Restore vectors
136d8606c27SBarry Smith      PetscCall(VecRestoreArrayF90(X,x_v,ierr))
137d8606c27SBarry Smith      PetscCall(VecRestoreArrayF90(F,f_v,ierr))
138c4762a1bSJed Brown
139c4762a1bSJed Brown      return
140c4762a1bSJed Brown      end
141c4762a1bSJed Brown
142c4762a1bSJed Brown      subroutine FormStartingPoint(x)
143*dfbbaf82SBarry Smith      use chwirut1fmodule
144c4762a1bSJed Brown
145c4762a1bSJed Brown      Vec             x
146a00960a2SBarry Smith      PetscScalar, pointer, dimension(:)  :: x_v
147c4762a1bSJed Brown      PetscErrorCode  ierr
148c4762a1bSJed Brown
149d8606c27SBarry Smith      PetscCall(VecGetArrayF90(x,x_v,ierr))
150a00960a2SBarry Smith      x_v(1) = 0.15
151a00960a2SBarry Smith      x_v(2) = 0.008
152a00960a2SBarry Smith      x_v(3) = 0.01
153d8606c27SBarry Smith      PetscCall(VecRestoreArrayF90(x,x_v,ierr))
154c4762a1bSJed Brown      return
155c4762a1bSJed Brown      end
156c4762a1bSJed Brown
157c4762a1bSJed Brown      subroutine InitializeData()
158*dfbbaf82SBarry Smith      use chwirut1fmodule
159c4762a1bSJed Brown
160c4762a1bSJed Brown      integer i
161c4762a1bSJed Brown      i=0
162c4762a1bSJed Brown      y(i) =    92.9000;  t(i) =  0.5000; i=i+1
163c4762a1bSJed Brown      y(i) =    78.7000;  t(i) =   0.6250; i=i+1
164c4762a1bSJed Brown      y(i) =    64.2000;  t(i) =   0.7500; i=i+1
165c4762a1bSJed Brown      y(i) =    64.9000;  t(i) =   0.8750; i=i+1
166c4762a1bSJed Brown      y(i) =    57.1000;  t(i) =   1.0000; i=i+1
167c4762a1bSJed Brown      y(i) =    43.3000;  t(i) =   1.2500; i=i+1
168c4762a1bSJed Brown      y(i) =    31.1000;  t(i) =  1.7500; i=i+1
169c4762a1bSJed Brown      y(i) =    23.6000;  t(i) =  2.2500; i=i+1
170c4762a1bSJed Brown      y(i) =    31.0500;  t(i) =  1.7500; i=i+1
171c4762a1bSJed Brown      y(i) =    23.7750;  t(i) =  2.2500; i=i+1
172c4762a1bSJed Brown      y(i) =    17.7375;  t(i) =  2.7500; i=i+1
173c4762a1bSJed Brown      y(i) =    13.8000;  t(i) =  3.2500; i=i+1
174c4762a1bSJed Brown      y(i) =    11.5875;  t(i) =  3.7500; i=i+1
175c4762a1bSJed Brown      y(i) =     9.4125;  t(i) =  4.2500; i=i+1
176c4762a1bSJed Brown      y(i) =     7.7250;  t(i) =  4.7500; i=i+1
177c4762a1bSJed Brown      y(i) =     7.3500;  t(i) =  5.2500; i=i+1
178c4762a1bSJed Brown      y(i) =     8.0250;  t(i) =  5.7500; i=i+1
179c4762a1bSJed Brown      y(i) =    90.6000;  t(i) =  0.5000; i=i+1
180c4762a1bSJed Brown      y(i) =    76.9000;  t(i) =  0.6250; i=i+1
181c4762a1bSJed Brown      y(i) =    71.6000;  t(i) = 0.7500; i=i+1
182c4762a1bSJed Brown      y(i) =    63.6000;  t(i) =  0.8750; i=i+1
183c4762a1bSJed Brown      y(i) =    54.0000;  t(i) =  1.0000; i=i+1
184c4762a1bSJed Brown      y(i) =    39.2000;  t(i) =  1.2500; i=i+1
185c4762a1bSJed Brown      y(i) =    29.3000;  t(i) = 1.7500; i=i+1
186c4762a1bSJed Brown      y(i) =    21.4000;  t(i) =  2.2500; i=i+1
187c4762a1bSJed Brown      y(i) =    29.1750;  t(i) =  1.7500; i=i+1
188c4762a1bSJed Brown      y(i) =    22.1250;  t(i) =  2.2500; i=i+1
189c4762a1bSJed Brown      y(i) =    17.5125;  t(i) =  2.7500; i=i+1
190c4762a1bSJed Brown      y(i) =    14.2500;  t(i) =  3.2500; i=i+1
191c4762a1bSJed Brown      y(i) =     9.4500;  t(i) =  3.7500; i=i+1
192c4762a1bSJed Brown      y(i) =     9.1500;  t(i) =  4.2500; i=i+1
193c4762a1bSJed Brown      y(i) =     7.9125;  t(i) =  4.7500; i=i+1
194c4762a1bSJed Brown      y(i) =     8.4750;  t(i) =  5.2500; i=i+1
195c4762a1bSJed Brown      y(i) =     6.1125;  t(i) =  5.7500; i=i+1
196c4762a1bSJed Brown      y(i) =    80.0000;  t(i) =  0.5000; i=i+1
197c4762a1bSJed Brown      y(i) =    79.0000;  t(i) =  0.6250; i=i+1
198c4762a1bSJed Brown      y(i) =    63.8000;  t(i) =  0.7500; i=i+1
199c4762a1bSJed Brown      y(i) =    57.2000;  t(i) =  0.8750; i=i+1
200c4762a1bSJed Brown      y(i) =    53.2000;  t(i) =  1.0000; i=i+1
201c4762a1bSJed Brown      y(i) =    42.5000;  t(i) =  1.2500; i=i+1
202c4762a1bSJed Brown      y(i) =    26.8000;  t(i) =  1.7500; i=i+1
203c4762a1bSJed Brown      y(i) =    20.4000;  t(i) =  2.2500; i=i+1
204c4762a1bSJed Brown      y(i) =    26.8500;  t(i) =   1.7500; i=i+1
205c4762a1bSJed Brown      y(i) =    21.0000;  t(i) =   2.2500; i=i+1
206c4762a1bSJed Brown      y(i) =    16.4625;  t(i) =   2.7500; i=i+1
207c4762a1bSJed Brown      y(i) =    12.5250;  t(i) =   3.2500; i=i+1
208c4762a1bSJed Brown      y(i) =    10.5375;  t(i) =   3.7500; i=i+1
209c4762a1bSJed Brown      y(i) =     8.5875;  t(i) =   4.2500; i=i+1
210c4762a1bSJed Brown      y(i) =     7.1250;  t(i) =   4.7500; i=i+1
211c4762a1bSJed Brown      y(i) =     6.1125;  t(i) =   5.2500; i=i+1
212c4762a1bSJed Brown      y(i) =     5.9625;  t(i) =   5.7500; i=i+1
213c4762a1bSJed Brown      y(i) =    74.1000;  t(i) =   0.5000; i=i+1
214c4762a1bSJed Brown      y(i) =    67.3000;  t(i) =   0.6250; i=i+1
215c4762a1bSJed Brown      y(i) =    60.8000;  t(i) =   0.7500; i=i+1
216c4762a1bSJed Brown      y(i) =    55.5000;  t(i) =   0.8750; i=i+1
217c4762a1bSJed Brown      y(i) =    50.3000;  t(i) =   1.0000; i=i+1
218c4762a1bSJed Brown      y(i) =    41.0000;  t(i) =   1.2500; i=i+1
219c4762a1bSJed Brown      y(i) =    29.4000;  t(i) =   1.7500; i=i+1
220c4762a1bSJed Brown      y(i) =    20.4000;  t(i) =   2.2500; i=i+1
221c4762a1bSJed Brown      y(i) =    29.3625;  t(i) =   1.7500; i=i+1
222c4762a1bSJed Brown      y(i) =    21.1500;  t(i) =   2.2500; i=i+1
223c4762a1bSJed Brown      y(i) =    16.7625;  t(i) =   2.7500; i=i+1
224c4762a1bSJed Brown      y(i) =    13.2000;  t(i) =   3.2500; i=i+1
225c4762a1bSJed Brown      y(i) =    10.8750;  t(i) =   3.7500; i=i+1
226c4762a1bSJed Brown      y(i) =     8.1750;  t(i) =   4.2500; i=i+1
227c4762a1bSJed Brown      y(i) =     7.3500;  t(i) =   4.7500; i=i+1
228c4762a1bSJed Brown      y(i) =     5.9625;  t(i) =  5.2500; i=i+1
229c4762a1bSJed Brown      y(i) =     5.6250;  t(i) =   5.7500; i=i+1
230c4762a1bSJed Brown      y(i) =    81.5000;  t(i) =    .5000; i=i+1
231c4762a1bSJed Brown      y(i) =    62.4000;  t(i) =    .7500; i=i+1
232c4762a1bSJed Brown      y(i) =    32.5000;  t(i) =   1.5000; i=i+1
233c4762a1bSJed Brown      y(i) =    12.4100;  t(i) =   3.0000; i=i+1
234c4762a1bSJed Brown      y(i) =    13.1200;  t(i) =   3.0000; i=i+1
235c4762a1bSJed Brown      y(i) =    15.5600;  t(i) =   3.0000; i=i+1
236c4762a1bSJed Brown      y(i) =     5.6300;  t(i) =   6.0000; i=i+1
237c4762a1bSJed Brown      y(i) =    78.0000;  t(i) =   .5000; i=i+1
238c4762a1bSJed Brown      y(i) =    59.9000;  t(i) =    .7500; i=i+1
239c4762a1bSJed Brown      y(i) =    33.2000;  t(i) =   1.5000; i=i+1
240c4762a1bSJed Brown      y(i) =    13.8400;  t(i) =   3.0000; i=i+1
241c4762a1bSJed Brown      y(i) =    12.7500;  t(i) =   3.0000; i=i+1
242c4762a1bSJed Brown      y(i) =    14.6200;  t(i) =   3.0000; i=i+1
243c4762a1bSJed Brown      y(i) =     3.9400;  t(i) =   6.0000; i=i+1
244c4762a1bSJed Brown      y(i) =    76.8000;  t(i) =    .5000; i=i+1
245c4762a1bSJed Brown      y(i) =    61.0000;  t(i) =    .7500; i=i+1
246c4762a1bSJed Brown      y(i) =    32.9000;  t(i) =   1.5000; i=i+1
247c4762a1bSJed Brown      y(i) =    13.8700;  t(i) = 3.0000; i=i+1
248c4762a1bSJed Brown      y(i) =    11.8100;  t(i) =   3.0000; i=i+1
249c4762a1bSJed Brown      y(i) =    13.3100;  t(i) =   3.0000; i=i+1
250c4762a1bSJed Brown      y(i) =     5.4400;  t(i) =   6.0000; i=i+1
251c4762a1bSJed Brown      y(i) =    78.0000;  t(i) =    .5000; i=i+1
252c4762a1bSJed Brown      y(i) =    63.5000;  t(i) =    .7500; i=i+1
253c4762a1bSJed Brown      y(i) =    33.8000;  t(i) =   1.5000; i=i+1
254c4762a1bSJed Brown      y(i) =    12.5600;  t(i) =   3.0000; i=i+1
255c4762a1bSJed Brown      y(i) =     5.6300;  t(i) =   6.0000; i=i+1
256c4762a1bSJed Brown      y(i) =    12.7500;  t(i) =   3.0000; i=i+1
257c4762a1bSJed Brown      y(i) =    13.1200;  t(i) =   3.0000; i=i+1
258c4762a1bSJed Brown      y(i) =     5.4400;  t(i) =   6.0000; i=i+1
259c4762a1bSJed Brown      y(i) =    76.8000;  t(i) =    .5000; i=i+1
260c4762a1bSJed Brown      y(i) =    60.0000;  t(i) =    .7500; i=i+1
261c4762a1bSJed Brown      y(i) =    47.8000;  t(i) =   1.0000; i=i+1
262c4762a1bSJed Brown      y(i) =    32.0000;  t(i) =   1.5000; i=i+1
263c4762a1bSJed Brown      y(i) =    22.2000;  t(i) =   2.0000; i=i+1
264c4762a1bSJed Brown      y(i) =    22.5700;  t(i) =   2.0000; i=i+1
265c4762a1bSJed Brown      y(i) =    18.8200;  t(i) =   2.5000; i=i+1
266c4762a1bSJed Brown      y(i) =    13.9500;  t(i) =   3.0000; i=i+1
267c4762a1bSJed Brown      y(i) =    11.2500;  t(i) =   4.0000; i=i+1
268c4762a1bSJed Brown      y(i) =     9.0000;  t(i) =   5.0000; i=i+1
269c4762a1bSJed Brown      y(i) =     6.6700;  t(i) =   6.0000; i=i+1
270c4762a1bSJed Brown      y(i) =    75.8000;  t(i) =    .5000; i=i+1
271c4762a1bSJed Brown      y(i) =    62.0000;  t(i) =    .7500; i=i+1
272c4762a1bSJed Brown      y(i) =    48.8000;  t(i) =   1.0000; i=i+1
273c4762a1bSJed Brown      y(i) =    35.2000;  t(i) =   1.5000; i=i+1
274c4762a1bSJed Brown      y(i) =    20.0000;  t(i) =   2.0000; i=i+1
275c4762a1bSJed Brown      y(i) =    20.3200;  t(i) =   2.0000; i=i+1
276c4762a1bSJed Brown      y(i) =    19.3100;  t(i) =   2.5000; i=i+1
277c4762a1bSJed Brown      y(i) =    12.7500;  t(i) =   3.0000; i=i+1
278c4762a1bSJed Brown      y(i) =    10.4200;  t(i) =   4.0000; i=i+1
279c4762a1bSJed Brown      y(i) =     7.3100;  t(i) =   5.0000; i=i+1
280c4762a1bSJed Brown      y(i) =     7.4200;  t(i) =   6.0000; i=i+1
281c4762a1bSJed Brown      y(i) =    70.5000;  t(i) =    .5000; i=i+1
282c4762a1bSJed Brown      y(i) =    59.5000;  t(i) =    .7500; i=i+1
283c4762a1bSJed Brown      y(i) =    48.5000;  t(i) =   1.0000; i=i+1
284c4762a1bSJed Brown      y(i) =    35.8000;  t(i) =   1.5000; i=i+1
285c4762a1bSJed Brown      y(i) =    21.0000;  t(i) =   2.0000; i=i+1
286c4762a1bSJed Brown      y(i) =    21.6700;  t(i) =   2.0000; i=i+1
287c4762a1bSJed Brown      y(i) =    21.0000;  t(i) =   2.5000; i=i+1
288c4762a1bSJed Brown      y(i) =    15.6400;  t(i) =   3.0000; i=i+1
289c4762a1bSJed Brown      y(i) =     8.1700;  t(i) =   4.0000; i=i+1
290c4762a1bSJed Brown      y(i) =     8.5500;  t(i) =   5.0000; i=i+1
291c4762a1bSJed Brown      y(i) =    10.1200;  t(i) =   6.0000; i=i+1
292c4762a1bSJed Brown      y(i) =    78.0000;  t(i) =    .5000; i=i+1
293c4762a1bSJed Brown      y(i) =    66.0000;  t(i) =    .6250; i=i+1
294c4762a1bSJed Brown      y(i) =    62.0000;  t(i) =    .7500; i=i+1
295c4762a1bSJed Brown      y(i) =    58.0000;  t(i) =    .8750; i=i+1
296c4762a1bSJed Brown      y(i) =    47.7000;  t(i) =   1.0000; i=i+1
297c4762a1bSJed Brown      y(i) =    37.8000;  t(i) =   1.2500; i=i+1
298c4762a1bSJed Brown      y(i) =    20.2000;  t(i) =   2.2500; i=i+1
299c4762a1bSJed Brown      y(i) =    21.0700;  t(i) =   2.2500; i=i+1
300c4762a1bSJed Brown      y(i) =    13.8700;  t(i) =   2.7500; i=i+1
301c4762a1bSJed Brown      y(i) =     9.6700;  t(i) =   3.2500; i=i+1
302c4762a1bSJed Brown      y(i) =     7.7600;  t(i) =   3.7500; i=i+1
303c4762a1bSJed Brown      y(i) =     5.4400;  t(i) =  4.2500; i=i+1
304c4762a1bSJed Brown      y(i) =     4.8700;  t(i) =  4.7500; i=i+1
305c4762a1bSJed Brown      y(i) =     4.0100;  t(i) =   5.2500; i=i+1
306c4762a1bSJed Brown      y(i) =     3.7500;  t(i) =   5.7500; i=i+1
307c4762a1bSJed Brown      y(i) =    24.1900;  t(i) =   3.0000; i=i+1
308c4762a1bSJed Brown      y(i) =    25.7600;  t(i) =   3.0000; i=i+1
309c4762a1bSJed Brown      y(i) =    18.0700;  t(i) =   3.0000; i=i+1
310c4762a1bSJed Brown      y(i) =    11.8100;  t(i) =   3.0000; i=i+1
311c4762a1bSJed Brown      y(i) =    12.0700;  t(i) =   3.0000; i=i+1
312c4762a1bSJed Brown      y(i) =    16.1200;  t(i) =   3.0000; i=i+1
313c4762a1bSJed Brown      y(i) =    70.8000;  t(i) =    .5000; i=i+1
314c4762a1bSJed Brown      y(i) =    54.7000;  t(i) =    .7500; i=i+1
315c4762a1bSJed Brown      y(i) =    48.0000;  t(i) =   1.0000; i=i+1
316c4762a1bSJed Brown      y(i) =    39.8000;  t(i) =   1.5000; i=i+1
317c4762a1bSJed Brown      y(i) =    29.8000;  t(i) =   2.0000; i=i+1
318c4762a1bSJed Brown      y(i) =    23.7000;  t(i) =   2.5000; i=i+1
319c4762a1bSJed Brown      y(i) =    29.6200;  t(i) =   2.0000; i=i+1
320c4762a1bSJed Brown      y(i) =    23.8100;  t(i) =   2.5000; i=i+1
321c4762a1bSJed Brown      y(i) =    17.7000;  t(i) =   3.0000; i=i+1
322c4762a1bSJed Brown      y(i) =    11.5500;  t(i) =   4.0000; i=i+1
323c4762a1bSJed Brown      y(i) =    12.0700;  t(i) =   5.0000; i=i+1
324c4762a1bSJed Brown      y(i) =     8.7400;  t(i) =   6.0000; i=i+1
325c4762a1bSJed Brown      y(i) =    80.7000;  t(i) =    .5000; i=i+1
326c4762a1bSJed Brown      y(i) =    61.3000;  t(i) =    .7500; i=i+1
327c4762a1bSJed Brown      y(i) =    47.5000;  t(i) =   1.0000; i=i+1
328c4762a1bSJed Brown      y(i) =    29.0000;  t(i) =   1.5000; i=i+1
329c4762a1bSJed Brown      y(i) =    24.0000;  t(i) =   2.0000; i=i+1
330c4762a1bSJed Brown      y(i) =    17.7000;  t(i) =   2.5000; i=i+1
331c4762a1bSJed Brown      y(i) =    24.5600;  t(i) =   2.0000; i=i+1
332c4762a1bSJed Brown      y(i) =    18.6700;  t(i) =   2.5000; i=i+1
333c4762a1bSJed Brown      y(i) =    16.2400;  t(i) =   3.0000; i=i+1
334c4762a1bSJed Brown      y(i) =     8.7400;  t(i) =   4.0000; i=i+1
335c4762a1bSJed Brown      y(i) =     7.8700;  t(i) =   5.0000; i=i+1
336c4762a1bSJed Brown      y(i) =     8.5100;  t(i) =   6.0000; i=i+1
337c4762a1bSJed Brown      y(i) =    66.7000;  t(i) =    .5000; i=i+1
338c4762a1bSJed Brown      y(i) =    59.2000;  t(i) =    .7500; i=i+1
339c4762a1bSJed Brown      y(i) =    40.8000;  t(i) =   1.0000; i=i+1
340c4762a1bSJed Brown      y(i) =    30.7000;  t(i) =   1.5000; i=i+1
341c4762a1bSJed Brown      y(i) =    25.7000;  t(i) =   2.0000; i=i+1
342c4762a1bSJed Brown      y(i) =    16.3000;  t(i) =   2.5000; i=i+1
343c4762a1bSJed Brown      y(i) =    25.9900;  t(i) =   2.0000; i=i+1
344c4762a1bSJed Brown      y(i) =    16.9500;  t(i) =   2.5000; i=i+1
345c4762a1bSJed Brown      y(i) =    13.3500;  t(i) =   3.0000; i=i+1
346c4762a1bSJed Brown      y(i) =     8.6200;  t(i) =   4.0000; i=i+1
347c4762a1bSJed Brown      y(i) =     7.2000;  t(i) =   5.0000; i=i+1
348c4762a1bSJed Brown      y(i) =     6.6400;  t(i) =   6.0000; i=i+1
349c4762a1bSJed Brown      y(i) =    13.6900;  t(i) =   3.0000; i=i+1
350c4762a1bSJed Brown      y(i) =    81.0000;  t(i) =    .5000; i=i+1
351c4762a1bSJed Brown      y(i) =    64.5000;  t(i) =    .7500; i=i+1
352c4762a1bSJed Brown      y(i) =    35.5000;  t(i) =   1.5000; i=i+1
353c4762a1bSJed Brown      y(i) =    13.3100;  t(i) =   3.0000; i=i+1
354c4762a1bSJed Brown      y(i) =     4.8700;  t(i) =   6.0000; i=i+1
355c4762a1bSJed Brown      y(i) =    12.9400;  t(i) =   3.0000; i=i+1
356c4762a1bSJed Brown      y(i) =     5.0600;  t(i) =   6.0000; i=i+1
357c4762a1bSJed Brown      y(i) =    15.1900;  t(i) =   3.0000; i=i+1
358c4762a1bSJed Brown      y(i) =    14.6200;  t(i) =   3.0000; i=i+1
359c4762a1bSJed Brown      y(i) =    15.6400;  t(i) =   3.0000; i=i+1
360c4762a1bSJed Brown      y(i) =    25.5000;  t(i) =   1.7500; i=i+1
361c4762a1bSJed Brown      y(i) =    25.9500;  t(i) =   1.7500; i=i+1
362c4762a1bSJed Brown      y(i) =    81.7000;  t(i) =    .5000; i=i+1
363c4762a1bSJed Brown      y(i) =    61.6000;  t(i) =    .7500; i=i+1
364c4762a1bSJed Brown      y(i) =    29.8000;  t(i) =   1.7500; i=i+1
365c4762a1bSJed Brown      y(i) =    29.8100;  t(i) =   1.7500; i=i+1
366c4762a1bSJed Brown      y(i) =    17.1700;  t(i) =   2.7500; i=i+1
367c4762a1bSJed Brown      y(i) =    10.3900;  t(i) =   3.7500; i=i+1
368c4762a1bSJed Brown      y(i) =    28.4000;  t(i) =   1.7500; i=i+1
369c4762a1bSJed Brown      y(i) =    28.6900;  t(i) =   1.7500; i=i+1
370c4762a1bSJed Brown      y(i) =    81.3000;  t(i) =    .5000; i=i+1
371c4762a1bSJed Brown      y(i) =    60.9000;  t(i) =    .7500; i=i+1
372c4762a1bSJed Brown      y(i) =    16.6500;  t(i) =   2.7500; i=i+1
373c4762a1bSJed Brown      y(i) =    10.0500;  t(i) =   3.7500; i=i+1
374c4762a1bSJed Brown      y(i) =    28.9000;  t(i) =   1.7500; i=i+1
375c4762a1bSJed Brown      y(i) =    28.9500;  t(i) =   1.7500; i=i+1
376c4762a1bSJed Brown
377c4762a1bSJed Brown      return
378c4762a1bSJed Brown      end
379c4762a1bSJed Brown
380c4762a1bSJed Brown!/*TEST
381c4762a1bSJed Brown!
382c4762a1bSJed Brown!   build:
383c4762a1bSJed Brown!      requires: !complex
384c4762a1bSJed Brown!
385c4762a1bSJed Brown!   test:
386c4762a1bSJed Brown!      args: -tao_smonitor -tao_max_it 100 -tao_type pounders -tao_gatol 1.e-5
387c4762a1bSJed Brown!      requires: !single
388c4762a1bSJed Brown!
389c4762a1bSJed Brown!TEST*/
390