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