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