xref: /petsc/src/tao/leastsquares/tutorials/chwirut2f.F90 (revision 749c190bad46ba447444c173d8c7a4080c70750e)
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 chwirut1.c
9c4762a1bSJed Brown!
10c5e229c2SMartin Diehl#include <petsc/finclude/petsctao.h>
11dfbbaf82SBarry Smithmodule chwirut2fmodule
12dfbbaf82SBarry Smith  use petsctao
13e7a95102SMartin Diehl  implicit none
14dfbbaf82SBarry Smith  PetscReal t(0:213)
15dfbbaf82SBarry Smith  PetscReal y(0:213)
16e7a95102SMartin Diehl  PetscInt, parameter :: m = 214, n = 3
17e7a95102SMartin Diehl  PetscMPIInt, parameter :: nn = n
18dfbbaf82SBarry Smith  PetscMPIInt rank
19dfbbaf82SBarry Smith  PetscMPIInt size
20e7a95102SMartin Diehl  PetscMPIInt, parameter :: idle_tag = 2000, die_tag = 3000
21e7a95102SMartin Diehl  PetscMPIInt, parameter :: zero = 0, one = 1
22dfbbaf82SBarry Smith
23e7a95102SMartin Diehlcontains
24c4762a1bSJed Brown  subroutine InitializeData()
25c4762a1bSJed Brown
26c4762a1bSJed Brown    PetscInt i
27c4762a1bSJed Brown    i = 0
28c4762a1bSJed Brown    y(i) = 92.9000; t(i) = 0.5000; i = i + 1
29c4762a1bSJed Brown    y(i) = 78.7000; t(i) = 0.6250; i = i + 1
30c4762a1bSJed Brown    y(i) = 64.2000; t(i) = 0.7500; i = i + 1
31c4762a1bSJed Brown    y(i) = 64.9000; t(i) = 0.8750; i = i + 1
32c4762a1bSJed Brown    y(i) = 57.1000; t(i) = 1.0000; i = i + 1
33c4762a1bSJed Brown    y(i) = 43.3000; t(i) = 1.2500; i = i + 1
34c4762a1bSJed Brown    y(i) = 31.1000; t(i) = 1.7500; i = i + 1
35c4762a1bSJed Brown    y(i) = 23.6000; t(i) = 2.2500; i = i + 1
36c4762a1bSJed Brown    y(i) = 31.0500; t(i) = 1.7500; i = i + 1
37c4762a1bSJed Brown    y(i) = 23.7750; t(i) = 2.2500; i = i + 1
38c4762a1bSJed Brown    y(i) = 17.7375; t(i) = 2.7500; i = i + 1
39c4762a1bSJed Brown    y(i) = 13.8000; t(i) = 3.2500; i = i + 1
40c4762a1bSJed Brown    y(i) = 11.5875; t(i) = 3.7500; i = i + 1
41c4762a1bSJed Brown    y(i) = 9.4125; t(i) = 4.2500; i = i + 1
42c4762a1bSJed Brown    y(i) = 7.7250; t(i) = 4.7500; i = i + 1
43c4762a1bSJed Brown    y(i) = 7.3500; t(i) = 5.2500; i = i + 1
44c4762a1bSJed Brown    y(i) = 8.0250; t(i) = 5.7500; i = i + 1
45c4762a1bSJed Brown    y(i) = 90.6000; t(i) = 0.5000; i = i + 1
46c4762a1bSJed Brown    y(i) = 76.9000; t(i) = 0.6250; i = i + 1
47c4762a1bSJed Brown    y(i) = 71.6000; t(i) = 0.7500; i = i + 1
48c4762a1bSJed Brown    y(i) = 63.6000; t(i) = 0.8750; i = i + 1
49c4762a1bSJed Brown    y(i) = 54.0000; t(i) = 1.0000; i = i + 1
50c4762a1bSJed Brown    y(i) = 39.2000; t(i) = 1.2500; i = i + 1
51c4762a1bSJed Brown    y(i) = 29.3000; t(i) = 1.7500; i = i + 1
52c4762a1bSJed Brown    y(i) = 21.4000; t(i) = 2.2500; i = i + 1
53c4762a1bSJed Brown    y(i) = 29.1750; t(i) = 1.7500; i = i + 1
54c4762a1bSJed Brown    y(i) = 22.1250; t(i) = 2.2500; i = i + 1
55c4762a1bSJed Brown    y(i) = 17.5125; t(i) = 2.7500; i = i + 1
56c4762a1bSJed Brown    y(i) = 14.2500; t(i) = 3.2500; i = i + 1
57c4762a1bSJed Brown    y(i) = 9.4500; t(i) = 3.7500; i = i + 1
58c4762a1bSJed Brown    y(i) = 9.1500; t(i) = 4.2500; i = i + 1
59c4762a1bSJed Brown    y(i) = 7.9125; t(i) = 4.7500; i = i + 1
60c4762a1bSJed Brown    y(i) = 8.4750; t(i) = 5.2500; i = i + 1
61c4762a1bSJed Brown    y(i) = 6.1125; t(i) = 5.7500; i = i + 1
62c4762a1bSJed Brown    y(i) = 80.0000; t(i) = 0.5000; i = i + 1
63c4762a1bSJed Brown    y(i) = 79.0000; t(i) = 0.6250; i = i + 1
64c4762a1bSJed Brown    y(i) = 63.8000; t(i) = 0.7500; i = i + 1
65c4762a1bSJed Brown    y(i) = 57.2000; t(i) = 0.8750; i = i + 1
66c4762a1bSJed Brown    y(i) = 53.2000; t(i) = 1.0000; i = i + 1
67c4762a1bSJed Brown    y(i) = 42.5000; t(i) = 1.2500; i = i + 1
68c4762a1bSJed Brown    y(i) = 26.8000; t(i) = 1.7500; i = i + 1
69c4762a1bSJed Brown    y(i) = 20.4000; t(i) = 2.2500; i = i + 1
70c4762a1bSJed Brown    y(i) = 26.8500; t(i) = 1.7500; i = i + 1
71c4762a1bSJed Brown    y(i) = 21.0000; t(i) = 2.2500; i = i + 1
72c4762a1bSJed Brown    y(i) = 16.4625; t(i) = 2.7500; i = i + 1
73c4762a1bSJed Brown    y(i) = 12.5250; t(i) = 3.2500; i = i + 1
74c4762a1bSJed Brown    y(i) = 10.5375; t(i) = 3.7500; i = i + 1
75c4762a1bSJed Brown    y(i) = 8.5875; t(i) = 4.2500; i = i + 1
76c4762a1bSJed Brown    y(i) = 7.1250; t(i) = 4.7500; i = i + 1
77c4762a1bSJed Brown    y(i) = 6.1125; t(i) = 5.2500; i = i + 1
78c4762a1bSJed Brown    y(i) = 5.9625; t(i) = 5.7500; i = i + 1
79c4762a1bSJed Brown    y(i) = 74.1000; t(i) = 0.5000; i = i + 1
80c4762a1bSJed Brown    y(i) = 67.3000; t(i) = 0.6250; i = i + 1
81c4762a1bSJed Brown    y(i) = 60.8000; t(i) = 0.7500; i = i + 1
82c4762a1bSJed Brown    y(i) = 55.5000; t(i) = 0.8750; i = i + 1
83c4762a1bSJed Brown    y(i) = 50.3000; t(i) = 1.0000; i = i + 1
84c4762a1bSJed Brown    y(i) = 41.0000; t(i) = 1.2500; i = i + 1
85c4762a1bSJed Brown    y(i) = 29.4000; t(i) = 1.7500; i = i + 1
86c4762a1bSJed Brown    y(i) = 20.4000; t(i) = 2.2500; i = i + 1
87c4762a1bSJed Brown    y(i) = 29.3625; t(i) = 1.7500; i = i + 1
88c4762a1bSJed Brown    y(i) = 21.1500; t(i) = 2.2500; i = i + 1
89c4762a1bSJed Brown    y(i) = 16.7625; t(i) = 2.7500; i = i + 1
90c4762a1bSJed Brown    y(i) = 13.2000; t(i) = 3.2500; i = i + 1
91c4762a1bSJed Brown    y(i) = 10.8750; t(i) = 3.7500; i = i + 1
92c4762a1bSJed Brown    y(i) = 8.1750; t(i) = 4.2500; i = i + 1
93c4762a1bSJed Brown    y(i) = 7.3500; t(i) = 4.7500; i = i + 1
94c4762a1bSJed Brown    y(i) = 5.9625; t(i) = 5.2500; i = i + 1
95c4762a1bSJed Brown    y(i) = 5.6250; t(i) = 5.7500; i = i + 1
96c4762a1bSJed Brown    y(i) = 81.5000; t(i) = .5000; i = i + 1
97c4762a1bSJed Brown    y(i) = 62.4000; t(i) = .7500; i = i + 1
98c4762a1bSJed Brown    y(i) = 32.5000; t(i) = 1.5000; i = i + 1
99c4762a1bSJed Brown    y(i) = 12.4100; t(i) = 3.0000; i = i + 1
100c4762a1bSJed Brown    y(i) = 13.1200; t(i) = 3.0000; i = i + 1
101c4762a1bSJed Brown    y(i) = 15.5600; t(i) = 3.0000; i = i + 1
102c4762a1bSJed Brown    y(i) = 5.6300; t(i) = 6.0000; i = i + 1
103c4762a1bSJed Brown    y(i) = 78.0000; t(i) = .5000; i = i + 1
104c4762a1bSJed Brown    y(i) = 59.9000; t(i) = .7500; i = i + 1
105c4762a1bSJed Brown    y(i) = 33.2000; t(i) = 1.5000; i = i + 1
106c4762a1bSJed Brown    y(i) = 13.8400; t(i) = 3.0000; i = i + 1
107c4762a1bSJed Brown    y(i) = 12.7500; t(i) = 3.0000; i = i + 1
108c4762a1bSJed Brown    y(i) = 14.6200; t(i) = 3.0000; i = i + 1
109c4762a1bSJed Brown    y(i) = 3.9400; t(i) = 6.0000; i = i + 1
110c4762a1bSJed Brown    y(i) = 76.8000; t(i) = .5000; i = i + 1
111c4762a1bSJed Brown    y(i) = 61.0000; t(i) = .7500; i = i + 1
112c4762a1bSJed Brown    y(i) = 32.9000; t(i) = 1.5000; i = i + 1
113c4762a1bSJed Brown    y(i) = 13.8700; t(i) = 3.0000; i = i + 1
114c4762a1bSJed Brown    y(i) = 11.8100; t(i) = 3.0000; i = i + 1
115c4762a1bSJed Brown    y(i) = 13.3100; t(i) = 3.0000; i = i + 1
116c4762a1bSJed Brown    y(i) = 5.4400; t(i) = 6.0000; i = i + 1
117c4762a1bSJed Brown    y(i) = 78.0000; t(i) = .5000; i = i + 1
118c4762a1bSJed Brown    y(i) = 63.5000; t(i) = .7500; i = i + 1
119c4762a1bSJed Brown    y(i) = 33.8000; t(i) = 1.5000; i = i + 1
120c4762a1bSJed Brown    y(i) = 12.5600; t(i) = 3.0000; i = i + 1
121c4762a1bSJed Brown    y(i) = 5.6300; t(i) = 6.0000; i = i + 1
122c4762a1bSJed Brown    y(i) = 12.7500; t(i) = 3.0000; i = i + 1
123c4762a1bSJed Brown    y(i) = 13.1200; t(i) = 3.0000; i = i + 1
124c4762a1bSJed Brown    y(i) = 5.4400; t(i) = 6.0000; i = i + 1
125c4762a1bSJed Brown    y(i) = 76.8000; t(i) = .5000; i = i + 1
126c4762a1bSJed Brown    y(i) = 60.0000; t(i) = .7500; i = i + 1
127c4762a1bSJed Brown    y(i) = 47.8000; t(i) = 1.0000; i = i + 1
128c4762a1bSJed Brown    y(i) = 32.0000; t(i) = 1.5000; i = i + 1
129c4762a1bSJed Brown    y(i) = 22.2000; t(i) = 2.0000; i = i + 1
130c4762a1bSJed Brown    y(i) = 22.5700; t(i) = 2.0000; i = i + 1
131c4762a1bSJed Brown    y(i) = 18.8200; t(i) = 2.5000; i = i + 1
132c4762a1bSJed Brown    y(i) = 13.9500; t(i) = 3.0000; i = i + 1
133c4762a1bSJed Brown    y(i) = 11.2500; t(i) = 4.0000; i = i + 1
134c4762a1bSJed Brown    y(i) = 9.0000; t(i) = 5.0000; i = i + 1
135c4762a1bSJed Brown    y(i) = 6.6700; t(i) = 6.0000; i = i + 1
136c4762a1bSJed Brown    y(i) = 75.8000; t(i) = .5000; i = i + 1
137c4762a1bSJed Brown    y(i) = 62.0000; t(i) = .7500; i = i + 1
138c4762a1bSJed Brown    y(i) = 48.8000; t(i) = 1.0000; i = i + 1
139c4762a1bSJed Brown    y(i) = 35.2000; t(i) = 1.5000; i = i + 1
140c4762a1bSJed Brown    y(i) = 20.0000; t(i) = 2.0000; i = i + 1
141c4762a1bSJed Brown    y(i) = 20.3200; t(i) = 2.0000; i = i + 1
142c4762a1bSJed Brown    y(i) = 19.3100; t(i) = 2.5000; i = i + 1
143c4762a1bSJed Brown    y(i) = 12.7500; t(i) = 3.0000; i = i + 1
144c4762a1bSJed Brown    y(i) = 10.4200; t(i) = 4.0000; i = i + 1
145c4762a1bSJed Brown    y(i) = 7.3100; t(i) = 5.0000; i = i + 1
146c4762a1bSJed Brown    y(i) = 7.4200; t(i) = 6.0000; i = i + 1
147c4762a1bSJed Brown    y(i) = 70.5000; t(i) = .5000; i = i + 1
148c4762a1bSJed Brown    y(i) = 59.5000; t(i) = .7500; i = i + 1
149c4762a1bSJed Brown    y(i) = 48.5000; t(i) = 1.0000; i = i + 1
150c4762a1bSJed Brown    y(i) = 35.8000; t(i) = 1.5000; i = i + 1
151c4762a1bSJed Brown    y(i) = 21.0000; t(i) = 2.0000; i = i + 1
152c4762a1bSJed Brown    y(i) = 21.6700; t(i) = 2.0000; i = i + 1
153c4762a1bSJed Brown    y(i) = 21.0000; t(i) = 2.5000; i = i + 1
154c4762a1bSJed Brown    y(i) = 15.6400; t(i) = 3.0000; i = i + 1
155c4762a1bSJed Brown    y(i) = 8.1700; t(i) = 4.0000; i = i + 1
156c4762a1bSJed Brown    y(i) = 8.5500; t(i) = 5.0000; i = i + 1
157c4762a1bSJed Brown    y(i) = 10.1200; t(i) = 6.0000; i = i + 1
158c4762a1bSJed Brown    y(i) = 78.0000; t(i) = .5000; i = i + 1
159c4762a1bSJed Brown    y(i) = 66.0000; t(i) = .6250; i = i + 1
160c4762a1bSJed Brown    y(i) = 62.0000; t(i) = .7500; i = i + 1
161c4762a1bSJed Brown    y(i) = 58.0000; t(i) = .8750; i = i + 1
162c4762a1bSJed Brown    y(i) = 47.7000; t(i) = 1.0000; i = i + 1
163c4762a1bSJed Brown    y(i) = 37.8000; t(i) = 1.2500; i = i + 1
164c4762a1bSJed Brown    y(i) = 20.2000; t(i) = 2.2500; i = i + 1
165c4762a1bSJed Brown    y(i) = 21.0700; t(i) = 2.2500; i = i + 1
166c4762a1bSJed Brown    y(i) = 13.8700; t(i) = 2.7500; i = i + 1
167c4762a1bSJed Brown    y(i) = 9.6700; t(i) = 3.2500; i = i + 1
168c4762a1bSJed Brown    y(i) = 7.7600; t(i) = 3.7500; i = i + 1
169c4762a1bSJed Brown    y(i) = 5.4400; t(i) = 4.2500; i = i + 1
170c4762a1bSJed Brown    y(i) = 4.8700; t(i) = 4.7500; i = i + 1
171c4762a1bSJed Brown    y(i) = 4.0100; t(i) = 5.2500; i = i + 1
172c4762a1bSJed Brown    y(i) = 3.7500; t(i) = 5.7500; i = i + 1
173c4762a1bSJed Brown    y(i) = 24.1900; t(i) = 3.0000; i = i + 1
174c4762a1bSJed Brown    y(i) = 25.7600; t(i) = 3.0000; i = i + 1
175c4762a1bSJed Brown    y(i) = 18.0700; t(i) = 3.0000; i = i + 1
176c4762a1bSJed Brown    y(i) = 11.8100; t(i) = 3.0000; i = i + 1
177c4762a1bSJed Brown    y(i) = 12.0700; t(i) = 3.0000; i = i + 1
178c4762a1bSJed Brown    y(i) = 16.1200; t(i) = 3.0000; i = i + 1
179c4762a1bSJed Brown    y(i) = 70.8000; t(i) = .5000; i = i + 1
180c4762a1bSJed Brown    y(i) = 54.7000; t(i) = .7500; i = i + 1
181c4762a1bSJed Brown    y(i) = 48.0000; t(i) = 1.0000; i = i + 1
182c4762a1bSJed Brown    y(i) = 39.8000; t(i) = 1.5000; i = i + 1
183c4762a1bSJed Brown    y(i) = 29.8000; t(i) = 2.0000; i = i + 1
184c4762a1bSJed Brown    y(i) = 23.7000; t(i) = 2.5000; i = i + 1
185c4762a1bSJed Brown    y(i) = 29.6200; t(i) = 2.0000; i = i + 1
186c4762a1bSJed Brown    y(i) = 23.8100; t(i) = 2.5000; i = i + 1
187c4762a1bSJed Brown    y(i) = 17.7000; t(i) = 3.0000; i = i + 1
188c4762a1bSJed Brown    y(i) = 11.5500; t(i) = 4.0000; i = i + 1
189c4762a1bSJed Brown    y(i) = 12.0700; t(i) = 5.0000; i = i + 1
190c4762a1bSJed Brown    y(i) = 8.7400; t(i) = 6.0000; i = i + 1
191c4762a1bSJed Brown    y(i) = 80.7000; t(i) = .5000; i = i + 1
192c4762a1bSJed Brown    y(i) = 61.3000; t(i) = .7500; i = i + 1
193c4762a1bSJed Brown    y(i) = 47.5000; t(i) = 1.0000; i = i + 1
194c4762a1bSJed Brown    y(i) = 29.0000; t(i) = 1.5000; i = i + 1
195c4762a1bSJed Brown    y(i) = 24.0000; t(i) = 2.0000; i = i + 1
196c4762a1bSJed Brown    y(i) = 17.7000; t(i) = 2.5000; i = i + 1
197c4762a1bSJed Brown    y(i) = 24.5600; t(i) = 2.0000; i = i + 1
198c4762a1bSJed Brown    y(i) = 18.6700; t(i) = 2.5000; i = i + 1
199c4762a1bSJed Brown    y(i) = 16.2400; t(i) = 3.0000; i = i + 1
200c4762a1bSJed Brown    y(i) = 8.7400; t(i) = 4.0000; i = i + 1
201c4762a1bSJed Brown    y(i) = 7.8700; t(i) = 5.0000; i = i + 1
202c4762a1bSJed Brown    y(i) = 8.5100; t(i) = 6.0000; i = i + 1
203c4762a1bSJed Brown    y(i) = 66.7000; t(i) = .5000; i = i + 1
204c4762a1bSJed Brown    y(i) = 59.2000; t(i) = .7500; i = i + 1
205c4762a1bSJed Brown    y(i) = 40.8000; t(i) = 1.0000; i = i + 1
206c4762a1bSJed Brown    y(i) = 30.7000; t(i) = 1.5000; i = i + 1
207c4762a1bSJed Brown    y(i) = 25.7000; t(i) = 2.0000; i = i + 1
208c4762a1bSJed Brown    y(i) = 16.3000; t(i) = 2.5000; i = i + 1
209c4762a1bSJed Brown    y(i) = 25.9900; t(i) = 2.0000; i = i + 1
210c4762a1bSJed Brown    y(i) = 16.9500; t(i) = 2.5000; i = i + 1
211c4762a1bSJed Brown    y(i) = 13.3500; t(i) = 3.0000; i = i + 1
212c4762a1bSJed Brown    y(i) = 8.6200; t(i) = 4.0000; i = i + 1
213c4762a1bSJed Brown    y(i) = 7.2000; t(i) = 5.0000; i = i + 1
214c4762a1bSJed Brown    y(i) = 6.6400; t(i) = 6.0000; i = i + 1
215c4762a1bSJed Brown    y(i) = 13.6900; t(i) = 3.0000; i = i + 1
216c4762a1bSJed Brown    y(i) = 81.0000; t(i) = .5000; i = i + 1
217c4762a1bSJed Brown    y(i) = 64.5000; t(i) = .7500; i = i + 1
218c4762a1bSJed Brown    y(i) = 35.5000; t(i) = 1.5000; i = i + 1
219c4762a1bSJed Brown    y(i) = 13.3100; t(i) = 3.0000; i = i + 1
220c4762a1bSJed Brown    y(i) = 4.8700; t(i) = 6.0000; i = i + 1
221c4762a1bSJed Brown    y(i) = 12.9400; t(i) = 3.0000; i = i + 1
222c4762a1bSJed Brown    y(i) = 5.0600; t(i) = 6.0000; i = i + 1
223c4762a1bSJed Brown    y(i) = 15.1900; t(i) = 3.0000; i = i + 1
224c4762a1bSJed Brown    y(i) = 14.6200; t(i) = 3.0000; i = i + 1
225c4762a1bSJed Brown    y(i) = 15.6400; t(i) = 3.0000; i = i + 1
226c4762a1bSJed Brown    y(i) = 25.5000; t(i) = 1.7500; i = i + 1
227c4762a1bSJed Brown    y(i) = 25.9500; t(i) = 1.7500; i = i + 1
228c4762a1bSJed Brown    y(i) = 81.7000; t(i) = .5000; i = i + 1
229c4762a1bSJed Brown    y(i) = 61.6000; t(i) = .7500; i = i + 1
230c4762a1bSJed Brown    y(i) = 29.8000; t(i) = 1.7500; i = i + 1
231c4762a1bSJed Brown    y(i) = 29.8100; t(i) = 1.7500; i = i + 1
232c4762a1bSJed Brown    y(i) = 17.1700; t(i) = 2.7500; i = i + 1
233c4762a1bSJed Brown    y(i) = 10.3900; t(i) = 3.7500; i = i + 1
234c4762a1bSJed Brown    y(i) = 28.4000; t(i) = 1.7500; i = i + 1
235c4762a1bSJed Brown    y(i) = 28.6900; t(i) = 1.7500; i = i + 1
236c4762a1bSJed Brown    y(i) = 81.3000; t(i) = .5000; i = i + 1
237c4762a1bSJed Brown    y(i) = 60.9000; t(i) = .7500; i = i + 1
238c4762a1bSJed Brown    y(i) = 16.6500; t(i) = 2.7500; i = i + 1
239c4762a1bSJed Brown    y(i) = 10.0500; t(i) = 3.7500; i = i + 1
240c4762a1bSJed Brown    y(i) = 28.9000; t(i) = 1.7500; i = i + 1
241c4762a1bSJed Brown    y(i) = 28.9500; t(i) = 1.7500; i = i + 1
242c4762a1bSJed Brown
243c4762a1bSJed Brown  end
244c4762a1bSJed Brown
245c4762a1bSJed Brown  subroutine TaskWorker(ierr)
246c4762a1bSJed Brown
247c4762a1bSJed Brown    PetscErrorCode ierr
24817a42bb7SSatish Balay    PetscReal x(n), f(1)
249c4762a1bSJed Brown    PetscMPIInt tag
250c4762a1bSJed Brown    PetscInt index
251*b06eb4cdSBarry Smith#if defined(PETSC_USE_MPI_F08)
252*b06eb4cdSBarry Smith    MPIU_Status status
253*b06eb4cdSBarry Smith#else
254*b06eb4cdSBarry Smith    MPIU_Status status(MPI_STATUS_SIZE)
255*b06eb4cdSBarry Smith#endif
256c4762a1bSJed Brown    tag = IDLE_TAG
257c4762a1bSJed Brown    f = 0.0
2589dddd249SSatish Balay    ! Send check-in message to rank-0
259d8606c27SBarry Smith    PetscCallMPI(MPI_Send(f, one, MPIU_SCALAR, zero, IDLE_TAG, PETSC_COMM_WORLD, ierr))
2604820e4eaSBarry Smith    do while (tag /= DIE_TAG)
261d8606c27SBarry Smith      PetscCallMPI(MPI_Recv(x, nn, MPIU_SCALAR, zero, MPI_ANY_TAG, PETSC_COMM_WORLD, status, ierr))
262*b06eb4cdSBarry Smith#if defined(PETSC_USE_MPI_F08)
263*b06eb4cdSBarry Smith      tag = status%MPI_TAG
264*b06eb4cdSBarry Smith#else
265c4762a1bSJed Brown      tag = status(MPI_TAG)
266*b06eb4cdSBarry Smith#endif
2674820e4eaSBarry Smith      if (tag == IDLE_TAG) then
268d8606c27SBarry Smith        PetscCallMPI(MPI_Send(f, one, MPIU_SCALAR, zero, IDLE_TAG, PETSC_COMM_WORLD, ierr))
2694820e4eaSBarry Smith      else if (tag /= DIE_TAG) then
270c4762a1bSJed Brown        index = tag
271c4762a1bSJed Brown        ! Compute local part of residual
272d8606c27SBarry Smith        PetscCall(RunSimulation(x, index, f(1), ierr))
273c4762a1bSJed Brown
2749dddd249SSatish Balay        ! Return residual to rank-0
275d8606c27SBarry Smith        PetscCallMPI(MPI_Send(f, one, MPIU_SCALAR, zero, tag, PETSC_COMM_WORLD, ierr))
276c4762a1bSJed Brown      end if
277c4762a1bSJed Brown    end do
278c4762a1bSJed Brown    ierr = 0
279c4762a1bSJed Brown  end
280c4762a1bSJed Brown
281c4762a1bSJed Brown  subroutine RunSimulation(x, i, f, ierr)
282c4762a1bSJed Brown
283c4762a1bSJed Brown    PetscReal x(n), f
284c4762a1bSJed Brown    PetscInt i
285c4762a1bSJed Brown    PetscErrorCode ierr
286c4762a1bSJed Brown    f = y(i) - exp(-x(1)*t(i))/(x(2) + x(3)*t(i))
287c4762a1bSJed Brown    ierr = 0
288c4762a1bSJed Brown  end
289c4762a1bSJed Brown
290c4762a1bSJed Brown  subroutine StopWorkers(ierr)
291c4762a1bSJed Brown
292c4762a1bSJed Brown    integer checkedin
293*b06eb4cdSBarry Smith#if defined(PETSC_USE_MPI_F08)
294*b06eb4cdSBarry Smith    MPIU_Status status
295*b06eb4cdSBarry Smith#else
296*b06eb4cdSBarry Smith    MPIU_Status status(MPI_STATUS_SIZE)
297*b06eb4cdSBarry Smith#endif
298c4762a1bSJed Brown    PetscMPIInt source
29917a42bb7SSatish Balay    PetscReal f(1), x(n)
300c4762a1bSJed Brown    PetscErrorCode ierr
301c4762a1bSJed Brown    PetscInt i
302c4762a1bSJed Brown
303c4762a1bSJed Brown    checkedin = 0
3044820e4eaSBarry Smith    do while (checkedin < size - 1)
305d8606c27SBarry Smith      PetscCallMPI(MPI_Recv(f, one, MPIU_SCALAR, MPI_ANY_SOURCE, MPI_ANY_TAG, PETSC_COMM_WORLD, status, ierr))
306c4762a1bSJed Brown      checkedin = checkedin + 1
307*b06eb4cdSBarry Smith#if defined(PETSC_USE_MPI_F08)
308*b06eb4cdSBarry Smith      source = status%MPI_SOURCE
309*b06eb4cdSBarry Smith#else
310c4762a1bSJed Brown      source = status(MPI_SOURCE)
311*b06eb4cdSBarry Smith#endif
312c4762a1bSJed Brown      do i = 1, n
313c4762a1bSJed Brown        x(i) = 0.0
314c4762a1bSJed Brown      end do
315d8606c27SBarry Smith      PetscCallMPI(MPI_Send(x, nn, MPIU_SCALAR, source, DIE_TAG, PETSC_COMM_WORLD, ierr))
316c4762a1bSJed Brown    end do
317c4762a1bSJed Brown    ierr = 0
318c4762a1bSJed Brown  end
319c4762a1bSJed Brown
320e7a95102SMartin Diehl! --------------------------------------------------------------------
321e7a95102SMartin Diehl!  FormFunction - Evaluates the function f(X) and gradient G(X)
322e7a95102SMartin Diehl!
323e7a95102SMartin Diehl!  Input Parameters:
324e7a95102SMartin Diehl!  tao - the Tao context
325e7a95102SMartin Diehl!  X   - input vector
326e7a95102SMartin Diehl!  dummy - not used
327e7a95102SMartin Diehl!
328e7a95102SMartin Diehl!  Output Parameters:
329e7a95102SMartin Diehl!  f - function vector
330e7a95102SMartin Diehl
331e7a95102SMartin Diehl  subroutine FormFunction(ta, x, f, dummy, ierr)
332e7a95102SMartin Diehl
333e7a95102SMartin Diehl    Tao ta
334e7a95102SMartin Diehl    Vec x, f
335e7a95102SMartin Diehl    PetscErrorCode ierr
336e7a95102SMartin Diehl
337e7a95102SMartin Diehl    PetscInt i, checkedin
338e7a95102SMartin Diehl    PetscInt finished_tasks
339e7a95102SMartin Diehl    PetscMPIInt next_task
340*b06eb4cdSBarry Smith    PetscMPIInt tag, source
341*b06eb4cdSBarry Smith#if defined(PETSC_USE_MPI_F08)
342*b06eb4cdSBarry Smith    MPIU_Status status
343*b06eb4cdSBarry Smith#else
344*b06eb4cdSBarry Smith    MPIU_Status status(MPI_STATUS_SIZE)
345*b06eb4cdSBarry Smith#endif
346e7a95102SMartin Diehl    PetscInt dummy
347e7a95102SMartin Diehl
348e7a95102SMartin Diehl    PetscReal, pointer :: f_v(:), x_v(:)
349e7a95102SMartin Diehl    PetscReal fval(1)
350e7a95102SMartin Diehl
351e7a95102SMartin Diehl    ierr = 0
352e7a95102SMartin Diehl
353e7a95102SMartin Diehl!     Get pointers to vector data
354e7a95102SMartin Diehl    PetscCall(VecGetArrayRead(x, x_v, ierr))
355e7a95102SMartin Diehl    PetscCall(VecGetArray(f, f_v, ierr))
356e7a95102SMartin Diehl
357e7a95102SMartin Diehl!     Compute F(X)
358e7a95102SMartin Diehl    if (size == 1) then
359e7a95102SMartin Diehl      ! Single processor
360e7a95102SMartin Diehl      do i = 1, m
361e7a95102SMartin Diehl        PetscCall(RunSimulation(x_v, i, f_v(i), ierr))
362e7a95102SMartin Diehl      end do
363e7a95102SMartin Diehl    else
364e7a95102SMartin Diehl      ! Multiprocessor main
365e7a95102SMartin Diehl      next_task = zero
366e7a95102SMartin Diehl      finished_tasks = 0
367e7a95102SMartin Diehl      checkedin = 0
368e7a95102SMartin Diehl
369e7a95102SMartin Diehl      do while (finished_tasks < m .or. checkedin < size - 1)
370e7a95102SMartin Diehl        PetscCallMPI(MPI_Recv(fval, one, MPIU_SCALAR, MPI_ANY_SOURCE, MPI_ANY_TAG, PETSC_COMM_WORLD, status, ierr))
371*b06eb4cdSBarry Smith#if defined(PETSC_USE_MPI_F08)
372*b06eb4cdSBarry Smith        tag = status%MPI_TAG
373*b06eb4cdSBarry Smith        source = status%MPI_SOURCE
374*b06eb4cdSBarry Smith#else
375e7a95102SMartin Diehl        tag = status(MPI_TAG)
376e7a95102SMartin Diehl        source = status(MPI_SOURCE)
377*b06eb4cdSBarry Smith#endif
378e7a95102SMartin Diehl        if (tag == IDLE_TAG) then
379e7a95102SMartin Diehl          checkedin = checkedin + 1
380e7a95102SMartin Diehl        else
381e7a95102SMartin Diehl          f_v(tag + 1) = fval(1)
382e7a95102SMartin Diehl          finished_tasks = finished_tasks + 1
383e7a95102SMartin Diehl        end if
384e7a95102SMartin Diehl        if (next_task < m) then
385e7a95102SMartin Diehl          ! Send task to worker
386e7a95102SMartin Diehl          PetscCallMPI(MPI_Send(x_v, nn, MPIU_SCALAR, source, next_task, PETSC_COMM_WORLD, ierr))
387e7a95102SMartin Diehl          next_task = next_task + one
388e7a95102SMartin Diehl        else
389e7a95102SMartin Diehl          ! Send idle message to worker
390e7a95102SMartin Diehl          PetscCallMPI(MPI_Send(x_v, nn, MPIU_SCALAR, source, IDLE_TAG, PETSC_COMM_WORLD, ierr))
391e7a95102SMartin Diehl        end if
392e7a95102SMartin Diehl      end do
393e7a95102SMartin Diehl    end if
394e7a95102SMartin Diehl
395e7a95102SMartin Diehl!     Restore vectors
396e7a95102SMartin Diehl    PetscCall(VecRestoreArrayRead(x, x_v, ierr))
397e7a95102SMartin Diehl    PetscCall(VecRestoreArray(F, f_v, ierr))
398e7a95102SMartin Diehl  end
399e7a95102SMartin Diehl
400e7a95102SMartin Diehl  subroutine FormStartingPoint(x)
401e7a95102SMartin Diehl
402e7a95102SMartin Diehl    Vec x
403e7a95102SMartin Diehl    PetscReal, pointer :: x_v(:)
404e7a95102SMartin Diehl    PetscErrorCode ierr
405e7a95102SMartin Diehl
406e7a95102SMartin Diehl    PetscCall(VecGetArray(x, x_v, ierr))
407e7a95102SMartin Diehl    x_v(1) = 0.15
408e7a95102SMartin Diehl    x_v(2) = 0.008
409e7a95102SMartin Diehl    x_v(3) = 0.01
410e7a95102SMartin Diehl    PetscCall(VecRestoreArray(x, x_v, ierr))
411e7a95102SMartin Diehl  end
412e7a95102SMartin Diehlend module chwirut2fmodule
413e7a95102SMartin Diehl
414e7a95102SMartin Diehlprogram main
415e7a95102SMartin Diehl  use chwirut2fmodule
416e7a95102SMartin Diehl  implicit none
417e7a95102SMartin Diehl! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
418e7a95102SMartin Diehl!                   Variable declarations
419e7a95102SMartin Diehl! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
420e7a95102SMartin Diehl!
421e7a95102SMartin Diehl!  See additional variable declarations in the file chwirut2f.h
422e7a95102SMartin Diehl
423e7a95102SMartin Diehl  PetscErrorCode ierr    ! used to check for functions returning nonzeros
424e7a95102SMartin Diehl  Vec x       ! solution vector
425e7a95102SMartin Diehl  Vec f       ! vector of functions
426e7a95102SMartin Diehl  Tao ta     ! Tao context
427e7a95102SMartin Diehl
428e7a95102SMartin Diehl!  Initialize TAO and PETSc
429e7a95102SMartin Diehl  PetscCallA(PetscInitialize(ierr))
430e7a95102SMartin Diehl  PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, size, ierr))
431e7a95102SMartin Diehl  PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))
432e7a95102SMartin Diehl
433e7a95102SMartin Diehl!  Initialize problem parameters
434e7a95102SMartin Diehl  call InitializeData()
435e7a95102SMartin Diehl
436e7a95102SMartin Diehl  if (rank == 0) then
437e7a95102SMartin Diehl!  Allocate vectors for the solution and gradient
438e7a95102SMartin Diehl    PetscCallA(VecCreateSeq(PETSC_COMM_SELF, n, x, ierr))
439e7a95102SMartin Diehl    PetscCallA(VecCreateSeq(PETSC_COMM_SELF, m, f, ierr))
440e7a95102SMartin Diehl
441e7a95102SMartin Diehl!     The TAO code begins here
442e7a95102SMartin Diehl
443e7a95102SMartin Diehl!     Create TAO solver
444e7a95102SMartin Diehl    PetscCallA(TaoCreate(PETSC_COMM_SELF, ta, ierr))
445e7a95102SMartin Diehl    PetscCallA(TaoSetType(ta, TAOPOUNDERS, ierr))
446e7a95102SMartin Diehl
447e7a95102SMartin Diehl!     Set routines for function, gradient, and hessian evaluation
448e7a95102SMartin Diehl    PetscCallA(TaoSetResidualRoutine(ta, f, FormFunction, 0, ierr))
449e7a95102SMartin Diehl
450e7a95102SMartin Diehl!     Optional: Set initial guess
451e7a95102SMartin Diehl    call FormStartingPoint(x)
452e7a95102SMartin Diehl    PetscCallA(TaoSetSolution(ta, x, ierr))
453e7a95102SMartin Diehl
454e7a95102SMartin Diehl!     Check for TAO command line options
455e7a95102SMartin Diehl    PetscCallA(TaoSetFromOptions(ta, ierr))
456e7a95102SMartin Diehl!     SOLVE THE APPLICATION
457e7a95102SMartin Diehl    PetscCallA(TaoSolve(ta, ierr))
458e7a95102SMartin Diehl
459e7a95102SMartin Diehl!     Free TAO data structures
460e7a95102SMartin Diehl    PetscCallA(TaoDestroy(ta, ierr))
461e7a95102SMartin Diehl
462e7a95102SMartin Diehl!     Free PETSc data structures
463e7a95102SMartin Diehl    PetscCallA(VecDestroy(x, ierr))
464e7a95102SMartin Diehl    PetscCallA(VecDestroy(f, ierr))
465e7a95102SMartin Diehl    PetscCallA(StopWorkers(ierr))
466e7a95102SMartin Diehl
467e7a95102SMartin Diehl  else
468e7a95102SMartin Diehl    PetscCallA(TaskWorker(ierr))
469e7a95102SMartin Diehl  end if
470e7a95102SMartin Diehl
471e7a95102SMartin Diehl  PetscCallA(PetscFinalize(ierr))
472e7a95102SMartin Diehlend
473c4762a1bSJed Brown!/*TEST
474c4762a1bSJed Brown!
475c4762a1bSJed Brown!   build:
476c4762a1bSJed Brown!      requires: !complex
477c4762a1bSJed Brown!
478c4762a1bSJed Brown!   test:
479c4762a1bSJed Brown!      nsize: 3
48010978b7dSBarry Smith!      args: -tao_monitor_short -tao_max_it 100 -tao_type pounders -tao_gatol 1.e-5
481c4762a1bSJed Brown!      requires: !single
482c4762a1bSJed Brown!
483c4762a1bSJed Brown!
484c4762a1bSJed Brown!TEST*/
485