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