Actual source code: chwirut2f.F

petsc-3.7.3 2016-08-01
Report Typos and Errors
  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: !/*T
 11: !  Concepts: TAO^Solving an unconstrained minimization problem
 12: !  Routines: TaoCreate();
 13: !  Routines: TaoSetType();
 14: !  Routines: TaoSetSeparableObjectiveRoutine();
 15: !  Routines: TaoSetInitialVector();
 16: !  Routines: TaoSetFromOptions();
 17: !  Routines: TaoSolve();
 18: !  Routines: TaoDestroy();
 19: !  Processors: n
 20: !T*/
 21: !
 22: ! ----------------------------------------------------------------------
 23: !
 24:       implicit none

 26: #include "chwirut2f.h"

 28: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 29: !                   Variable declarations
 30: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 31: !
 32: !  See additional variable declarations in the file chwirut2f.h

 34:       PetscErrorCode   ierr    ! used to check for functions returning nonzeros
 35:       Vec              x       ! solution vector
 36:       Vec              f       ! vector of functions
 37:       Tao        tao     ! Tao context

 39: !  Note: Any user-defined Fortran routines (such as FormGradient)
 40: !  MUST be declared as external.

 42:       external FormFunction

 44: !  Initialize TAO and PETSc
 45:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)

 47:       call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
 48:       CHKERRQ(ierr)
 49:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
 50:       CHKERRQ(ierr)

 52: !  Initialize problem parameters
 53:       call InitializeData()

 55:       if (rank .eq. 0) then
 56: !  Allocate vectors for the solution and gradient
 57:          call VecCreateSeq(PETSC_COMM_SELF,n,x,ierr)
 58:          CHKERRQ(ierr)
 59:          call VecCreateSeq(PETSC_COMM_SELF,m,f,ierr)
 60:          CHKERRQ(ierr)


 63: !     The TAO code begins here

 65: !     Create TAO solver
 66:          call TaoCreate(PETSC_COMM_SELF,tao,ierr)
 67:          CHKERRQ(ierr)
 68:          call TaoSetType(tao,TAOPOUNDERS,ierr)
 69:          CHKERRQ(ierr)

 71: !     Set routines for function, gradient, and hessian evaluation
 72:          call TaoSetSeparableObjectiveRoutine(tao,f,                    &
 73:      &        FormFunction,PETSC_NULL_OBJECT,ierr)
 74:          CHKERRQ(ierr)

 76: !     Optional: Set initial guess
 77:          call FormStartingPoint(x)
 78:          call TaoSetInitialVector(tao, x, ierr)
 79:          CHKERRQ(ierr)


 82: !     Check for TAO command line options
 83:          call TaoSetFromOptions(tao,ierr)
 84:          CHKERRQ(ierr)
 85: !     SOLVE THE APPLICATION
 86:          call TaoSolve(tao,ierr)
 87:          CHKERRQ(ierr)

 89: !     Free TAO data structures
 90:          call TaoDestroy(tao,ierr)
 91:          CHKERRQ(ierr)

 93: !     Free PETSc data structures
 94:          call VecDestroy(x,ierr)
 95:          CHKERRQ(ierr)
 96:          call VecDestroy(f,ierr)
 97:          CHKERRQ(ierr)
 98:          call StopWorkers(ierr)
 99:          CHKERRQ(ierr)

101:       else
102:          call TaskWorker(ierr)
103:          CHKERRQ(ierr)
104:       endif

106:       call PetscFinalize(ierr)
107:       end


110: ! --------------------------------------------------------------------
111: !  FormFunction - Evaluates the function f(X) and gradient G(X)
112: !
113: !  Input Parameters:
114: !  tao - the Tao context
115: !  X   - input vector
116: !  dummy - not used
117: !
118: !  Output Parameters:
119: !  f - function vector

121:       subroutine FormFunction(tao, x, f, dummy, ierr)
122:       implicit none

124: ! n,alpha defined in chwirut2f.h
125: #include "chwirut2f.h"

127:       Tao        tao
128:       Vec              x,f
129:       PetscErrorCode   ierr
130:       PetscInt         dummy

132:       PetscInt         i,checkedin
133:       PetscInt         finished_tasks
134:       integer          next_task,status(MPI_STATUS_SIZE),tag,source

136: ! PETSc's VecGetArray acts differently in Fortran than it does in C.
137: ! Calling VecGetArray((Vec) X, (PetscReal) x_array(0:1), (PetscOffset) x_index, ierr)
138: ! will return an array of doubles referenced by x_array offset by x_index.
139: !  i.e.,  to reference the kth element of X, use x_array(k + x_index).
140: ! Notice that by declaring the arrays with range (0:1), we are using the C 0-indexing practice.
141:       PetscReal        f_v(0:1),x_v(0:1),fval
142:       PetscOffset      f_i,x_i

144:       0

146: !     Get pointers to vector data
147:       call VecGetArray(x,x_v,x_i,ierr)
148:       CHKERRQ(ierr)
149:       call VecGetArray(f,f_v,f_i,ierr)
150:       CHKERRQ(ierr)


153: !     Compute F(X)
154:       if (size .eq. 1) then
155:          ! Single processor
156:          do i=0,m-1
157:             call RunSimulation(x_v(x_i),i,f_v(i+f_i),ierr)
158:          enddo
159:       else
160:          ! Multiprocessor master
161:          next_task = 0
162:          finished_tasks = 0
163:          checkedin = 0

165:          do while (finished_tasks .lt. m .or. checkedin .lt. size-1)
166:             call MPI_Recv(fval,1,MPIU_SCALAR,MPI_ANY_SOURCE,               &
167:      &           MPI_ANY_TAG,PETSC_COMM_WORLD,status,ierr)
168:             tag = status(MPI_TAG)
169:             source = status(MPI_SOURCE)
170:             if (tag .eq. IDLE_TAG) then
171:                checkedin = checkedin + 1
172:             else
173:                f_v(f_i+tag) = fval
174:                finished_tasks = finished_tasks + 1
175:             endif
176:             if (next_task .lt. m) then
177:                ! Send task to worker
178:                call MPI_Send(x_v(x_i),n,MPIU_SCALAR,source,next_task,             &
179:      &              PETSC_COMM_WORLD,ierr)
180:                next_task = next_task + 1
181:             else
182:                ! Send idle message to worker
183:                call MPI_Send(x_v(x_i),n,MPIU_SCALAR,source,IDLE_TAG,              &
184:      &              PETSC_COMM_WORLD,ierr)
185:             end if
186:          enddo
187:       endif

189: !     Restore vectors
190:       call VecRestoreArray(x,x_v,x_i,ierr)
191:       CHKERRQ(ierr)
192:       call VecRestoreArray(F,f_v,f_i,ierr)
193:       CHKERRQ(ierr)
194:       return
195:       end





201:       subroutine FormStartingPoint(x)
202:       implicit none

204: ! n,alpha defined in chwirut2f.h
205: #include "chwirut2f.h"
206:       Vec             x
207:       PetscReal       x_v(0:1)
208:       PetscOffset     x_i
209:       PetscErrorCode  ierr

211:       call VecGetArray(x,x_v,x_i,ierr)
212:       CHKERRQ(ierr)
213:       x_v(x_i) = 0.15
214:       x_v(x_i+1) = 0.008
215:       x_v(x_i+2) = 0.01
216:       call VecRestoreArray(x,x_v,x_i,ierr)
217:       CHKERRQ(ierr)
218:       return
219:       end


222:       subroutine InitializeData()
223:       implicit none

225: ! n,alpha defined in chwirut2f.h
226: #include "chwirut2f.h"
227:       PetscInt i
228:       i=0
229:       y(i) =    92.9000;  t(i) =  0.5000; i=i+1
230:       y(i) =    78.7000;  t(i) =   0.6250; i=i+1
231:       y(i) =    64.2000;  t(i) =   0.7500; i=i+1
232:       y(i) =    64.9000;  t(i) =   0.8750; i=i+1
233:       y(i) =    57.1000;  t(i) =   1.0000; i=i+1
234:       y(i) =    43.3000;  t(i) =   1.2500; i=i+1
235:       y(i) =    31.1000;  t(i) =  1.7500; i=i+1
236:       y(i) =    23.6000;  t(i) =  2.2500; i=i+1
237:       y(i) =    31.0500;  t(i) =  1.7500; i=i+1
238:       y(i) =    23.7750;  t(i) =  2.2500; i=i+1
239:       y(i) =    17.7375;  t(i) =  2.7500; i=i+1
240:       y(i) =    13.8000;  t(i) =  3.2500; i=i+1
241:       y(i) =    11.5875;  t(i) =  3.7500; i=i+1
242:       y(i) =     9.4125;  t(i) =  4.2500; i=i+1
243:       y(i) =     7.7250;  t(i) =  4.7500; i=i+1
244:       y(i) =     7.3500;  t(i) =  5.2500; i=i+1
245:       y(i) =     8.0250;  t(i) =  5.7500; i=i+1
246:       y(i) =    90.6000;  t(i) =  0.5000; i=i+1
247:       y(i) =    76.9000;  t(i) =  0.6250; i=i+1
248:       y(i) =    71.6000;  t(i) = 0.7500; i=i+1
249:       y(i) =    63.6000;  t(i) =  0.8750; i=i+1
250:       y(i) =    54.0000;  t(i) =  1.0000; i=i+1
251:       y(i) =    39.2000;  t(i) =  1.2500; i=i+1
252:       y(i) =    29.3000;  t(i) = 1.7500; i=i+1
253:       y(i) =    21.4000;  t(i) =  2.2500; i=i+1
254:       y(i) =    29.1750;  t(i) =  1.7500; i=i+1
255:       y(i) =    22.1250;  t(i) =  2.2500; i=i+1
256:       y(i) =    17.5125;  t(i) =  2.7500; i=i+1
257:       y(i) =    14.2500;  t(i) =  3.2500; i=i+1
258:       y(i) =     9.4500;  t(i) =  3.7500; i=i+1
259:       y(i) =     9.1500;  t(i) =  4.2500; i=i+1
260:       y(i) =     7.9125;  t(i) =  4.7500; i=i+1
261:       y(i) =     8.4750;  t(i) =  5.2500; i=i+1
262:       y(i) =     6.1125;  t(i) =  5.7500; i=i+1
263:       y(i) =    80.0000;  t(i) =  0.5000; i=i+1
264:       y(i) =    79.0000;  t(i) =  0.6250; i=i+1
265:       y(i) =    63.8000;  t(i) =  0.7500; i=i+1
266:       y(i) =    57.2000;  t(i) =  0.8750; i=i+1
267:       y(i) =    53.2000;  t(i) =  1.0000; i=i+1
268:       y(i) =    42.5000;  t(i) =  1.2500; i=i+1
269:       y(i) =    26.8000;  t(i) =  1.7500; i=i+1
270:       y(i) =    20.4000;  t(i) =  2.2500; i=i+1
271:       y(i) =    26.8500;  t(i) =   1.7500; i=i+1
272:       y(i) =    21.0000;  t(i) =   2.2500; i=i+1
273:       y(i) =    16.4625;  t(i) =   2.7500; i=i+1
274:       y(i) =    12.5250;  t(i) =   3.2500; i=i+1
275:       y(i) =    10.5375;  t(i) =   3.7500; i=i+1
276:       y(i) =     8.5875;  t(i) =   4.2500; i=i+1
277:       y(i) =     7.1250;  t(i) =   4.7500; i=i+1
278:       y(i) =     6.1125;  t(i) =   5.2500; i=i+1
279:       y(i) =     5.9625;  t(i) =   5.7500; i=i+1
280:       y(i) =    74.1000;  t(i) =   0.5000; i=i+1
281:       y(i) =    67.3000;  t(i) =   0.6250; i=i+1
282:       y(i) =    60.8000;  t(i) =   0.7500; i=i+1
283:       y(i) =    55.5000;  t(i) =   0.8750; i=i+1
284:       y(i) =    50.3000;  t(i) =   1.0000; i=i+1
285:       y(i) =    41.0000;  t(i) =   1.2500; i=i+1
286:       y(i) =    29.4000;  t(i) =   1.7500; i=i+1
287:       y(i) =    20.4000;  t(i) =   2.2500; i=i+1
288:       y(i) =    29.3625;  t(i) =   1.7500; i=i+1
289:       y(i) =    21.1500;  t(i) =   2.2500; i=i+1
290:       y(i) =    16.7625;  t(i) =   2.7500; i=i+1
291:       y(i) =    13.2000;  t(i) =   3.2500; i=i+1
292:       y(i) =    10.8750;  t(i) =   3.7500; i=i+1
293:       y(i) =     8.1750;  t(i) =   4.2500; i=i+1
294:       y(i) =     7.3500;  t(i) =   4.7500; i=i+1
295:       y(i) =     5.9625;  t(i) =  5.2500; i=i+1
296:       y(i) =     5.6250;  t(i) =   5.7500; i=i+1
297:       y(i) =    81.5000;  t(i) =    .5000; i=i+1
298:       y(i) =    62.4000;  t(i) =    .7500; i=i+1
299:       y(i) =    32.5000;  t(i) =   1.5000; i=i+1
300:       y(i) =    12.4100;  t(i) =   3.0000; i=i+1
301:       y(i) =    13.1200;  t(i) =   3.0000; i=i+1
302:       y(i) =    15.5600;  t(i) =   3.0000; i=i+1
303:       y(i) =     5.6300;  t(i) =   6.0000; i=i+1
304:       y(i) =    78.0000;  t(i) =   .5000; i=i+1
305:       y(i) =    59.9000;  t(i) =    .7500; i=i+1
306:       y(i) =    33.2000;  t(i) =   1.5000; i=i+1
307:       y(i) =    13.8400;  t(i) =   3.0000; i=i+1
308:       y(i) =    12.7500;  t(i) =   3.0000; i=i+1
309:       y(i) =    14.6200;  t(i) =   3.0000; i=i+1
310:       y(i) =     3.9400;  t(i) =   6.0000; i=i+1
311:       y(i) =    76.8000;  t(i) =    .5000; i=i+1
312:       y(i) =    61.0000;  t(i) =    .7500; i=i+1
313:       y(i) =    32.9000;  t(i) =   1.5000; i=i+1
314:       y(i) =    13.8700;  t(i) = 3.0000; i=i+1
315:       y(i) =    11.8100;  t(i) =   3.0000; i=i+1
316:       y(i) =    13.3100;  t(i) =   3.0000; i=i+1
317:       y(i) =     5.4400;  t(i) =   6.0000; i=i+1
318:       y(i) =    78.0000;  t(i) =    .5000; i=i+1
319:       y(i) =    63.5000;  t(i) =    .7500; i=i+1
320:       y(i) =    33.8000;  t(i) =   1.5000; i=i+1
321:       y(i) =    12.5600;  t(i) =   3.0000; i=i+1
322:       y(i) =     5.6300;  t(i) =   6.0000; i=i+1
323:       y(i) =    12.7500;  t(i) =   3.0000; i=i+1
324:       y(i) =    13.1200;  t(i) =   3.0000; i=i+1
325:       y(i) =     5.4400;  t(i) =   6.0000; i=i+1
326:       y(i) =    76.8000;  t(i) =    .5000; i=i+1
327:       y(i) =    60.0000;  t(i) =    .7500; i=i+1
328:       y(i) =    47.8000;  t(i) =   1.0000; i=i+1
329:       y(i) =    32.0000;  t(i) =   1.5000; i=i+1
330:       y(i) =    22.2000;  t(i) =   2.0000; i=i+1
331:       y(i) =    22.5700;  t(i) =   2.0000; i=i+1
332:       y(i) =    18.8200;  t(i) =   2.5000; i=i+1
333:       y(i) =    13.9500;  t(i) =   3.0000; i=i+1
334:       y(i) =    11.2500;  t(i) =   4.0000; i=i+1
335:       y(i) =     9.0000;  t(i) =   5.0000; i=i+1
336:       y(i) =     6.6700;  t(i) =   6.0000; i=i+1
337:       y(i) =    75.8000;  t(i) =    .5000; i=i+1
338:       y(i) =    62.0000;  t(i) =    .7500; i=i+1
339:       y(i) =    48.8000;  t(i) =   1.0000; i=i+1
340:       y(i) =    35.2000;  t(i) =   1.5000; i=i+1
341:       y(i) =    20.0000;  t(i) =   2.0000; i=i+1
342:       y(i) =    20.3200;  t(i) =   2.0000; i=i+1
343:       y(i) =    19.3100;  t(i) =   2.5000; i=i+1
344:       y(i) =    12.7500;  t(i) =   3.0000; i=i+1
345:       y(i) =    10.4200;  t(i) =   4.0000; i=i+1
346:       y(i) =     7.3100;  t(i) =   5.0000; i=i+1
347:       y(i) =     7.4200;  t(i) =   6.0000; i=i+1
348:       y(i) =    70.5000;  t(i) =    .5000; i=i+1
349:       y(i) =    59.5000;  t(i) =    .7500; i=i+1
350:       y(i) =    48.5000;  t(i) =   1.0000; i=i+1
351:       y(i) =    35.8000;  t(i) =   1.5000; i=i+1
352:       y(i) =    21.0000;  t(i) =   2.0000; i=i+1
353:       y(i) =    21.6700;  t(i) =   2.0000; i=i+1
354:       y(i) =    21.0000;  t(i) =   2.5000; i=i+1
355:       y(i) =    15.6400;  t(i) =   3.0000; i=i+1
356:       y(i) =     8.1700;  t(i) =   4.0000; i=i+1
357:       y(i) =     8.5500;  t(i) =   5.0000; i=i+1
358:       y(i) =    10.1200;  t(i) =   6.0000; i=i+1
359:       y(i) =    78.0000;  t(i) =    .5000; i=i+1
360:       y(i) =    66.0000;  t(i) =    .6250; i=i+1
361:       y(i) =    62.0000;  t(i) =    .7500; i=i+1
362:       y(i) =    58.0000;  t(i) =    .8750; i=i+1
363:       y(i) =    47.7000;  t(i) =   1.0000; i=i+1
364:       y(i) =    37.8000;  t(i) =   1.2500; i=i+1
365:       y(i) =    20.2000;  t(i) =   2.2500; i=i+1
366:       y(i) =    21.0700;  t(i) =   2.2500; i=i+1
367:       y(i) =    13.8700;  t(i) =   2.7500; i=i+1
368:       y(i) =     9.6700;  t(i) =   3.2500; i=i+1
369:       y(i) =     7.7600;  t(i) =   3.7500; i=i+1
370:       y(i) =     5.4400;  t(i) =  4.2500; i=i+1
371:       y(i) =     4.8700;  t(i) =  4.7500; i=i+1
372:       y(i) =     4.0100;  t(i) =   5.2500; i=i+1
373:       y(i) =     3.7500;  t(i) =   5.7500; i=i+1
374:       y(i) =    24.1900;  t(i) =   3.0000; i=i+1
375:       y(i) =    25.7600;  t(i) =   3.0000; i=i+1
376:       y(i) =    18.0700;  t(i) =   3.0000; i=i+1
377:       y(i) =    11.8100;  t(i) =   3.0000; i=i+1
378:       y(i) =    12.0700;  t(i) =   3.0000; i=i+1
379:       y(i) =    16.1200;  t(i) =   3.0000; i=i+1
380:       y(i) =    70.8000;  t(i) =    .5000; i=i+1
381:       y(i) =    54.7000;  t(i) =    .7500; i=i+1
382:       y(i) =    48.0000;  t(i) =   1.0000; i=i+1
383:       y(i) =    39.8000;  t(i) =   1.5000; i=i+1
384:       y(i) =    29.8000;  t(i) =   2.0000; i=i+1
385:       y(i) =    23.7000;  t(i) =   2.5000; i=i+1
386:       y(i) =    29.6200;  t(i) =   2.0000; i=i+1
387:       y(i) =    23.8100;  t(i) =   2.5000; i=i+1
388:       y(i) =    17.7000;  t(i) =   3.0000; i=i+1
389:       y(i) =    11.5500;  t(i) =   4.0000; i=i+1
390:       y(i) =    12.0700;  t(i) =   5.0000; i=i+1
391:       y(i) =     8.7400;  t(i) =   6.0000; i=i+1
392:       y(i) =    80.7000;  t(i) =    .5000; i=i+1
393:       y(i) =    61.3000;  t(i) =    .7500; i=i+1
394:       y(i) =    47.5000;  t(i) =   1.0000; i=i+1
395:       y(i) =    29.0000;  t(i) =   1.5000; i=i+1
396:       y(i) =    24.0000;  t(i) =   2.0000; i=i+1
397:       y(i) =    17.7000;  t(i) =   2.5000; i=i+1
398:       y(i) =    24.5600;  t(i) =   2.0000; i=i+1
399:       y(i) =    18.6700;  t(i) =   2.5000; i=i+1
400:       y(i) =    16.2400;  t(i) =   3.0000; i=i+1
401:       y(i) =     8.7400;  t(i) =   4.0000; i=i+1
402:       y(i) =     7.8700;  t(i) =   5.0000; i=i+1
403:       y(i) =     8.5100;  t(i) =   6.0000; i=i+1
404:       y(i) =    66.7000;  t(i) =    .5000; i=i+1
405:       y(i) =    59.2000;  t(i) =    .7500; i=i+1
406:       y(i) =    40.8000;  t(i) =   1.0000; i=i+1
407:       y(i) =    30.7000;  t(i) =   1.5000; i=i+1
408:       y(i) =    25.7000;  t(i) =   2.0000; i=i+1
409:       y(i) =    16.3000;  t(i) =   2.5000; i=i+1
410:       y(i) =    25.9900;  t(i) =   2.0000; i=i+1
411:       y(i) =    16.9500;  t(i) =   2.5000; i=i+1
412:       y(i) =    13.3500;  t(i) =   3.0000; i=i+1
413:       y(i) =     8.6200;  t(i) =   4.0000; i=i+1
414:       y(i) =     7.2000;  t(i) =   5.0000; i=i+1
415:       y(i) =     6.6400;  t(i) =   6.0000; i=i+1
416:       y(i) =    13.6900;  t(i) =   3.0000; i=i+1
417:       y(i) =    81.0000;  t(i) =    .5000; i=i+1
418:       y(i) =    64.5000;  t(i) =    .7500; i=i+1
419:       y(i) =    35.5000;  t(i) =   1.5000; i=i+1
420:       y(i) =    13.3100;  t(i) =   3.0000; i=i+1
421:       y(i) =     4.8700;  t(i) =   6.0000; i=i+1
422:       y(i) =    12.9400;  t(i) =   3.0000; i=i+1
423:       y(i) =     5.0600;  t(i) =   6.0000; i=i+1
424:       y(i) =    15.1900;  t(i) =   3.0000; i=i+1
425:       y(i) =    14.6200;  t(i) =   3.0000; i=i+1
426:       y(i) =    15.6400;  t(i) =   3.0000; i=i+1
427:       y(i) =    25.5000;  t(i) =   1.7500; i=i+1
428:       y(i) =    25.9500;  t(i) =   1.7500; i=i+1
429:       y(i) =    81.7000;  t(i) =    .5000; i=i+1
430:       y(i) =    61.6000;  t(i) =    .7500; i=i+1
431:       y(i) =    29.8000;  t(i) =   1.7500; i=i+1
432:       y(i) =    29.8100;  t(i) =   1.7500; i=i+1
433:       y(i) =    17.1700;  t(i) =   2.7500; i=i+1
434:       y(i) =    10.3900;  t(i) =   3.7500; i=i+1
435:       y(i) =    28.4000;  t(i) =   1.7500; i=i+1
436:       y(i) =    28.6900;  t(i) =   1.7500; i=i+1
437:       y(i) =    81.3000;  t(i) =    .5000; i=i+1
438:       y(i) =    60.9000;  t(i) =    .7500; i=i+1
439:       y(i) =    16.6500;  t(i) =   2.7500; i=i+1
440:       y(i) =    10.0500;  t(i) =   3.7500; i=i+1
441:       y(i) =    28.9000;  t(i) =   1.7500; i=i+1
442:       y(i) =    28.9500;  t(i) =   1.7500; i=i+1

444:       return
445:       end



449:       subroutine TaskWorker(ierr)
450:       implicit none
451: #include "chwirut2f.h"
452:       PetscErrorCode ierr
453:       PetscReal x(n),f
454:       integer tag
455:       PetscInt index
456:       integer status(MPI_STATUS_SIZE)

458:       tag = IDLE_TAG
459:       f   = 0.0
460:       ! Send check-in message to master
461:       call MPI_Send(f,1,MPIU_SCALAR,0,IDLE_TAG,PETSC_COMM_WORLD,ierr)
462:       CHKERRQ(ierr)
463:       do while (tag .ne. DIE_TAG)
464:          call MPI_Recv(x,n,MPIU_SCALAR,0,MPI_ANY_TAG,PETSC_COMM_WORLD,     &
465:      &        status,ierr)
466:          CHKERRQ(ierr)
467:          tag = status(MPI_TAG)
468:          if (tag .eq. IDLE_TAG) then
469:             call MPI_Send(f,1,MPIU_SCALAR,0,IDLE_TAG,PETSC_COMM_WORLD,     &
470:      &           ierr)
471:             CHKERRQ(ierr)
472:          else if (tag .ne. DIE_TAG) then
473:             index = tag
474:             ! Compute local part of residual
475:             call RunSimulation(x,index,f,ierr)
476:             CHKERRQ(ierr)

478:             ! Return residual to master
479:             call MPI_Send(f,1,MPIU_SCALAR,0,tag,PETSC_COMM_WORLD,ierr)
480:             CHKERRQ(ierr)
481:          end if
482:       enddo
483:       0
484:       return
485:       end



489:       subroutine RunSimulation(x,i,f,ierr)
490:       implicit none
491: #include "chwirut2f.h"
492:       PetscReal x(n),f
493:       PetscInt i
494:       PetscErrorCode ierr
495:       f = y(i) - exp(-x(1)*t(i))/(x(2)+x(3)*t(i))
496:       0
497:       return
498:       end

500:       subroutine StopWorkers(ierr)
501:       implicit none
502: #include "chwirut2f.h"
503:       integer checkedin
504:       integer status(MPI_STATUS_SIZE)
505:       integer source
506:       PetscReal f,x(n)
507:       PetscErrorCode ierr
508:       PetscInt i

510:       checkedin=0
511:       do while (checkedin .lt. size-1)
512:          call MPI_Recv(f,1,MPIU_SCALAR,MPI_ANY_SOURCE,MPI_ANY_TAG,         &
513:      &        PETSC_COMM_WORLD,status,ierr)
514:          CHKERRQ(ierr)
515:          checkedin=checkedin+1
516:          source = status(MPI_SOURCE)
517:          do i=1,n
518:            x(i) = 0.0
519:          enddo
520:          call MPI_Send(x,n,MPIU_SCALAR,source,DIE_TAG,PETSC_COMM_WORLD,    &
521:      &        ierr)
522:          CHKERRQ(ierr)
523:       enddo
524:       ierr=0
525:       return
526:       end