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