Actual source code: chwirut2.c

  1: /*
  2:    Include "petsctao.h" so that we can use TAO solvers.  Note that this
  3:    file automatically includes libraries such as:
  4:      petsc.h       - base PETSc routines   petscvec.h - vectors
  5:      petscsys.h    - system routines        petscmat.h - matrices
  6:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
  7:      petscviewer.h - viewers               petscpc.h  - preconditioners

  9:  This version tests correlated terms using both vector and listed forms
 10: */

 12: #include <petsctao.h>

 14: /*
 15: Description:   These data are the result of a NIST study involving
 16:                ultrasonic calibration.  The response variable is
 17:                ultrasonic response, and the predictor variable is
 18:                metal distance.

 20: Reference:     Chwirut, D., NIST (197?).
 21:                Ultrasonic Reference Block Study.
 22: */

 24: static char help[] = "Finds the nonlinear least-squares solution to the model \n\
 25:             y = exp[-b1*x]/(b2+b3*x)  +  e \n";

 27: #define NOBSERVATIONS 214
 28: #define NPARAMETERS   3

 30: /* User-defined application context */
 31: typedef struct {
 32:   /* Working space */
 33:   PetscReal t[NOBSERVATIONS];              /* array of independent variables of observation */
 34:   PetscReal y[NOBSERVATIONS];              /* array of dependent variables */
 35:   PetscReal j[NOBSERVATIONS][NPARAMETERS]; /* dense jacobian matrix array*/
 36:   PetscInt  idm[NOBSERVATIONS];            /* Matrix indices for jacobian */
 37:   PetscInt  idn[NPARAMETERS];
 38: } AppCtx;

 40: /* User provided Routines */
 41: PetscErrorCode InitializeData(AppCtx *user);
 42: PetscErrorCode FormStartingPoint(Vec);
 43: PetscErrorCode EvaluateFunction(Tao, Vec, Vec, void *);
 44: PetscErrorCode EvaluateJacobian(Tao, Vec, Mat, Mat, void *);

 46: /*--------------------------------------------------------------------*/
 47: int main(int argc, char **argv)
 48: {
 49:   PetscInt  wtype = 0;
 50:   Vec       x, f; /* solution, function */
 51:   Vec       w;    /* weights */
 52:   Mat       J;    /* Jacobian matrix */
 53:   Tao       tao;  /* Tao solver context */
 54:   PetscInt  i;    /* iteration information */
 55:   PetscReal hist[100], resid[100];
 56:   PetscInt  lits[100];
 57:   PetscInt  w_row[NOBSERVATIONS]; /* explicit weights */
 58:   PetscInt  w_col[NOBSERVATIONS];
 59:   PetscReal w_vals[NOBSERVATIONS];
 60:   PetscBool flg;
 61:   AppCtx    user; /* user-defined work context */

 64:   PetscInitialize(&argc, &argv, (char *)0, help);
 65:   PetscOptionsGetInt(NULL, NULL, "-wtype", &wtype, &flg);
 66:   PetscPrintf(PETSC_COMM_WORLD, "wtype=%" PetscInt_FMT "\n", wtype);
 67:   /* Allocate vectors */
 68:   VecCreateSeq(MPI_COMM_SELF, NPARAMETERS, &x);
 69:   VecCreateSeq(MPI_COMM_SELF, NOBSERVATIONS, &f);

 71:   VecDuplicate(f, &w);

 73:   /* no correlation, but set in different ways */
 74:   VecSet(w, 1.0);
 75:   for (i = 0; i < NOBSERVATIONS; i++) {
 76:     w_row[i]  = i;
 77:     w_col[i]  = i;
 78:     w_vals[i] = 1.0;
 79:   }

 81:   /* Create the Jacobian matrix. */
 82:   MatCreateSeqDense(MPI_COMM_SELF, NOBSERVATIONS, NPARAMETERS, NULL, &J);

 84:   for (i = 0; i < NOBSERVATIONS; i++) user.idm[i] = i;

 86:   for (i = 0; i < NPARAMETERS; i++) user.idn[i] = i;

 88:   /* Create TAO solver and set desired solution method */
 89:   TaoCreate(PETSC_COMM_SELF, &tao);
 90:   TaoSetType(tao, TAOPOUNDERS);

 92:   /* Set the function and Jacobian routines. */
 93:   InitializeData(&user);
 94:   FormStartingPoint(x);
 95:   TaoSetSolution(tao, x);
 96:   TaoSetResidualRoutine(tao, f, EvaluateFunction, (void *)&user);
 97:   if (wtype == 1) {
 98:     TaoSetResidualWeights(tao, w, 0, NULL, NULL, NULL);
 99:   } else if (wtype == 2) {
100:     TaoSetResidualWeights(tao, NULL, NOBSERVATIONS, w_row, w_col, w_vals);
101:   }
102:   TaoSetJacobianResidualRoutine(tao, J, J, EvaluateJacobian, (void *)&user);
103:   TaoSetTolerances(tao, 1e-5, 0.0, PETSC_DEFAULT);

105:   /* Check for any TAO command line arguments */
106:   TaoSetFromOptions(tao);

108:   TaoSetConvergenceHistory(tao, hist, resid, 0, lits, 100, PETSC_TRUE);
109:   /* Perform the Solve */
110:   TaoSolve(tao);

112:   /* Free TAO data structures */
113:   TaoDestroy(&tao);

115:   /* Free PETSc data structures */
116:   VecDestroy(&x);
117:   VecDestroy(&w);
118:   VecDestroy(&f);
119:   MatDestroy(&J);

121:   PetscFinalize();
122:   return 0;
123: }

125: /*--------------------------------------------------------------------*/
126: PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr)
127: {
128:   AppCtx          *user = (AppCtx *)ptr;
129:   PetscInt         i;
130:   PetscReal       *y = user->y, *f, *t = user->t;
131:   const PetscReal *x;

133:   VecGetArrayRead(X, &x);
134:   VecGetArray(F, &f);

136:   for (i = 0; i < NOBSERVATIONS; i++) f[i] = y[i] - PetscExpScalar(-x[0] * t[i]) / (x[1] + x[2] * t[i]);
137:   VecRestoreArrayRead(X, &x);
138:   VecRestoreArray(F, &f);
139:   PetscLogFlops(6 * NOBSERVATIONS);
140:   return 0;
141: }

143: /*------------------------------------------------------------*/
144: /* J[i][j] = df[i]/dt[j] */
145: PetscErrorCode EvaluateJacobian(Tao tao, Vec X, Mat J, Mat Jpre, void *ptr)
146: {
147:   AppCtx          *user = (AppCtx *)ptr;
148:   PetscInt         i;
149:   PetscReal       *t = user->t;
150:   const PetscReal *x;
151:   PetscReal        base;

153:   VecGetArrayRead(X, &x);
154:   for (i = 0; i < NOBSERVATIONS; i++) {
155:     base = PetscExpScalar(-x[0] * t[i]) / (x[1] + x[2] * t[i]);

157:     user->j[i][0] = t[i] * base;
158:     user->j[i][1] = base / (x[1] + x[2] * t[i]);
159:     user->j[i][2] = base * t[i] / (x[1] + x[2] * t[i]);
160:   }

162:   /* Assemble the matrix */
163:   MatSetValues(J, NOBSERVATIONS, user->idm, NPARAMETERS, user->idn, (PetscReal *)user->j, INSERT_VALUES);
164:   MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);
165:   MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);

167:   VecRestoreArrayRead(X, &x);
168:   PetscLogFlops(NOBSERVATIONS * 13);
169:   return 0;
170: }

172: /* ------------------------------------------------------------ */
173: PetscErrorCode FormStartingPoint(Vec X)
174: {
175:   PetscReal *x;

177:   VecGetArray(X, &x);
178:   x[0] = 1.19;
179:   x[1] = -1.86;
180:   x[2] = 1.08;
181:   VecRestoreArray(X, &x);
182:   return 0;
183: }

185: /* ---------------------------------------------------------------------- */
186: PetscErrorCode InitializeData(AppCtx *user)
187: {
188:   PetscReal *t = user->t, *y = user->y;
189:   PetscInt   i = 0;

191:   y[i]   = 92.9000;
192:   t[i++] = 0.5000;
193:   y[i]   = 78.7000;
194:   t[i++] = 0.6250;
195:   y[i]   = 64.2000;
196:   t[i++] = 0.7500;
197:   y[i]   = 64.9000;
198:   t[i++] = 0.8750;
199:   y[i]   = 57.1000;
200:   t[i++] = 1.0000;
201:   y[i]   = 43.3000;
202:   t[i++] = 1.2500;
203:   y[i]   = 31.1000;
204:   t[i++] = 1.7500;
205:   y[i]   = 23.6000;
206:   t[i++] = 2.2500;
207:   y[i]   = 31.0500;
208:   t[i++] = 1.7500;
209:   y[i]   = 23.7750;
210:   t[i++] = 2.2500;
211:   y[i]   = 17.7375;
212:   t[i++] = 2.7500;
213:   y[i]   = 13.8000;
214:   t[i++] = 3.2500;
215:   y[i]   = 11.5875;
216:   t[i++] = 3.7500;
217:   y[i]   = 9.4125;
218:   t[i++] = 4.2500;
219:   y[i]   = 7.7250;
220:   t[i++] = 4.7500;
221:   y[i]   = 7.3500;
222:   t[i++] = 5.2500;
223:   y[i]   = 8.0250;
224:   t[i++] = 5.7500;
225:   y[i]   = 90.6000;
226:   t[i++] = 0.5000;
227:   y[i]   = 76.9000;
228:   t[i++] = 0.6250;
229:   y[i]   = 71.6000;
230:   t[i++] = 0.7500;
231:   y[i]   = 63.6000;
232:   t[i++] = 0.8750;
233:   y[i]   = 54.0000;
234:   t[i++] = 1.0000;
235:   y[i]   = 39.2000;
236:   t[i++] = 1.2500;
237:   y[i]   = 29.3000;
238:   t[i++] = 1.7500;
239:   y[i]   = 21.4000;
240:   t[i++] = 2.2500;
241:   y[i]   = 29.1750;
242:   t[i++] = 1.7500;
243:   y[i]   = 22.1250;
244:   t[i++] = 2.2500;
245:   y[i]   = 17.5125;
246:   t[i++] = 2.7500;
247:   y[i]   = 14.2500;
248:   t[i++] = 3.2500;
249:   y[i]   = 9.4500;
250:   t[i++] = 3.7500;
251:   y[i]   = 9.1500;
252:   t[i++] = 4.2500;
253:   y[i]   = 7.9125;
254:   t[i++] = 4.7500;
255:   y[i]   = 8.4750;
256:   t[i++] = 5.2500;
257:   y[i]   = 6.1125;
258:   t[i++] = 5.7500;
259:   y[i]   = 80.0000;
260:   t[i++] = 0.5000;
261:   y[i]   = 79.0000;
262:   t[i++] = 0.6250;
263:   y[i]   = 63.8000;
264:   t[i++] = 0.7500;
265:   y[i]   = 57.2000;
266:   t[i++] = 0.8750;
267:   y[i]   = 53.2000;
268:   t[i++] = 1.0000;
269:   y[i]   = 42.5000;
270:   t[i++] = 1.2500;
271:   y[i]   = 26.8000;
272:   t[i++] = 1.7500;
273:   y[i]   = 20.4000;
274:   t[i++] = 2.2500;
275:   y[i]   = 26.8500;
276:   t[i++] = 1.7500;
277:   y[i]   = 21.0000;
278:   t[i++] = 2.2500;
279:   y[i]   = 16.4625;
280:   t[i++] = 2.7500;
281:   y[i]   = 12.5250;
282:   t[i++] = 3.2500;
283:   y[i]   = 10.5375;
284:   t[i++] = 3.7500;
285:   y[i]   = 8.5875;
286:   t[i++] = 4.2500;
287:   y[i]   = 7.1250;
288:   t[i++] = 4.7500;
289:   y[i]   = 6.1125;
290:   t[i++] = 5.2500;
291:   y[i]   = 5.9625;
292:   t[i++] = 5.7500;
293:   y[i]   = 74.1000;
294:   t[i++] = 0.5000;
295:   y[i]   = 67.3000;
296:   t[i++] = 0.6250;
297:   y[i]   = 60.8000;
298:   t[i++] = 0.7500;
299:   y[i]   = 55.5000;
300:   t[i++] = 0.8750;
301:   y[i]   = 50.3000;
302:   t[i++] = 1.0000;
303:   y[i]   = 41.0000;
304:   t[i++] = 1.2500;
305:   y[i]   = 29.4000;
306:   t[i++] = 1.7500;
307:   y[i]   = 20.4000;
308:   t[i++] = 2.2500;
309:   y[i]   = 29.3625;
310:   t[i++] = 1.7500;
311:   y[i]   = 21.1500;
312:   t[i++] = 2.2500;
313:   y[i]   = 16.7625;
314:   t[i++] = 2.7500;
315:   y[i]   = 13.2000;
316:   t[i++] = 3.2500;
317:   y[i]   = 10.8750;
318:   t[i++] = 3.7500;
319:   y[i]   = 8.1750;
320:   t[i++] = 4.2500;
321:   y[i]   = 7.3500;
322:   t[i++] = 4.7500;
323:   y[i]   = 5.9625;
324:   t[i++] = 5.2500;
325:   y[i]   = 5.6250;
326:   t[i++] = 5.7500;
327:   y[i]   = 81.5000;
328:   t[i++] = .5000;
329:   y[i]   = 62.4000;
330:   t[i++] = .7500;
331:   y[i]   = 32.5000;
332:   t[i++] = 1.5000;
333:   y[i]   = 12.4100;
334:   t[i++] = 3.0000;
335:   y[i]   = 13.1200;
336:   t[i++] = 3.0000;
337:   y[i]   = 15.5600;
338:   t[i++] = 3.0000;
339:   y[i]   = 5.6300;
340:   t[i++] = 6.0000;
341:   y[i]   = 78.0000;
342:   t[i++] = .5000;
343:   y[i]   = 59.9000;
344:   t[i++] = .7500;
345:   y[i]   = 33.2000;
346:   t[i++] = 1.5000;
347:   y[i]   = 13.8400;
348:   t[i++] = 3.0000;
349:   y[i]   = 12.7500;
350:   t[i++] = 3.0000;
351:   y[i]   = 14.6200;
352:   t[i++] = 3.0000;
353:   y[i]   = 3.9400;
354:   t[i++] = 6.0000;
355:   y[i]   = 76.8000;
356:   t[i++] = .5000;
357:   y[i]   = 61.0000;
358:   t[i++] = .7500;
359:   y[i]   = 32.9000;
360:   t[i++] = 1.5000;
361:   y[i]   = 13.8700;
362:   t[i++] = 3.0000;
363:   y[i]   = 11.8100;
364:   t[i++] = 3.0000;
365:   y[i]   = 13.3100;
366:   t[i++] = 3.0000;
367:   y[i]   = 5.4400;
368:   t[i++] = 6.0000;
369:   y[i]   = 78.0000;
370:   t[i++] = .5000;
371:   y[i]   = 63.5000;
372:   t[i++] = .7500;
373:   y[i]   = 33.8000;
374:   t[i++] = 1.5000;
375:   y[i]   = 12.5600;
376:   t[i++] = 3.0000;
377:   y[i]   = 5.6300;
378:   t[i++] = 6.0000;
379:   y[i]   = 12.7500;
380:   t[i++] = 3.0000;
381:   y[i]   = 13.1200;
382:   t[i++] = 3.0000;
383:   y[i]   = 5.4400;
384:   t[i++] = 6.0000;
385:   y[i]   = 76.8000;
386:   t[i++] = .5000;
387:   y[i]   = 60.0000;
388:   t[i++] = .7500;
389:   y[i]   = 47.8000;
390:   t[i++] = 1.0000;
391:   y[i]   = 32.0000;
392:   t[i++] = 1.5000;
393:   y[i]   = 22.2000;
394:   t[i++] = 2.0000;
395:   y[i]   = 22.5700;
396:   t[i++] = 2.0000;
397:   y[i]   = 18.8200;
398:   t[i++] = 2.5000;
399:   y[i]   = 13.9500;
400:   t[i++] = 3.0000;
401:   y[i]   = 11.2500;
402:   t[i++] = 4.0000;
403:   y[i]   = 9.0000;
404:   t[i++] = 5.0000;
405:   y[i]   = 6.6700;
406:   t[i++] = 6.0000;
407:   y[i]   = 75.8000;
408:   t[i++] = .5000;
409:   y[i]   = 62.0000;
410:   t[i++] = .7500;
411:   y[i]   = 48.8000;
412:   t[i++] = 1.0000;
413:   y[i]   = 35.2000;
414:   t[i++] = 1.5000;
415:   y[i]   = 20.0000;
416:   t[i++] = 2.0000;
417:   y[i]   = 20.3200;
418:   t[i++] = 2.0000;
419:   y[i]   = 19.3100;
420:   t[i++] = 2.5000;
421:   y[i]   = 12.7500;
422:   t[i++] = 3.0000;
423:   y[i]   = 10.4200;
424:   t[i++] = 4.0000;
425:   y[i]   = 7.3100;
426:   t[i++] = 5.0000;
427:   y[i]   = 7.4200;
428:   t[i++] = 6.0000;
429:   y[i]   = 70.5000;
430:   t[i++] = .5000;
431:   y[i]   = 59.5000;
432:   t[i++] = .7500;
433:   y[i]   = 48.5000;
434:   t[i++] = 1.0000;
435:   y[i]   = 35.8000;
436:   t[i++] = 1.5000;
437:   y[i]   = 21.0000;
438:   t[i++] = 2.0000;
439:   y[i]   = 21.6700;
440:   t[i++] = 2.0000;
441:   y[i]   = 21.0000;
442:   t[i++] = 2.5000;
443:   y[i]   = 15.6400;
444:   t[i++] = 3.0000;
445:   y[i]   = 8.1700;
446:   t[i++] = 4.0000;
447:   y[i]   = 8.5500;
448:   t[i++] = 5.0000;
449:   y[i]   = 10.1200;
450:   t[i++] = 6.0000;
451:   y[i]   = 78.0000;
452:   t[i++] = .5000;
453:   y[i]   = 66.0000;
454:   t[i++] = .6250;
455:   y[i]   = 62.0000;
456:   t[i++] = .7500;
457:   y[i]   = 58.0000;
458:   t[i++] = .8750;
459:   y[i]   = 47.7000;
460:   t[i++] = 1.0000;
461:   y[i]   = 37.8000;
462:   t[i++] = 1.2500;
463:   y[i]   = 20.2000;
464:   t[i++] = 2.2500;
465:   y[i]   = 21.0700;
466:   t[i++] = 2.2500;
467:   y[i]   = 13.8700;
468:   t[i++] = 2.7500;
469:   y[i]   = 9.6700;
470:   t[i++] = 3.2500;
471:   y[i]   = 7.7600;
472:   t[i++] = 3.7500;
473:   y[i]   = 5.4400;
474:   t[i++] = 4.2500;
475:   y[i]   = 4.8700;
476:   t[i++] = 4.7500;
477:   y[i]   = 4.0100;
478:   t[i++] = 5.2500;
479:   y[i]   = 3.7500;
480:   t[i++] = 5.7500;
481:   y[i]   = 24.1900;
482:   t[i++] = 3.0000;
483:   y[i]   = 25.7600;
484:   t[i++] = 3.0000;
485:   y[i]   = 18.0700;
486:   t[i++] = 3.0000;
487:   y[i]   = 11.8100;
488:   t[i++] = 3.0000;
489:   y[i]   = 12.0700;
490:   t[i++] = 3.0000;
491:   y[i]   = 16.1200;
492:   t[i++] = 3.0000;
493:   y[i]   = 70.8000;
494:   t[i++] = .5000;
495:   y[i]   = 54.7000;
496:   t[i++] = .7500;
497:   y[i]   = 48.0000;
498:   t[i++] = 1.0000;
499:   y[i]   = 39.8000;
500:   t[i++] = 1.5000;
501:   y[i]   = 29.8000;
502:   t[i++] = 2.0000;
503:   y[i]   = 23.7000;
504:   t[i++] = 2.5000;
505:   y[i]   = 29.6200;
506:   t[i++] = 2.0000;
507:   y[i]   = 23.8100;
508:   t[i++] = 2.5000;
509:   y[i]   = 17.7000;
510:   t[i++] = 3.0000;
511:   y[i]   = 11.5500;
512:   t[i++] = 4.0000;
513:   y[i]   = 12.0700;
514:   t[i++] = 5.0000;
515:   y[i]   = 8.7400;
516:   t[i++] = 6.0000;
517:   y[i]   = 80.7000;
518:   t[i++] = .5000;
519:   y[i]   = 61.3000;
520:   t[i++] = .7500;
521:   y[i]   = 47.5000;
522:   t[i++] = 1.0000;
523:   y[i]   = 29.0000;
524:   t[i++] = 1.5000;
525:   y[i]   = 24.0000;
526:   t[i++] = 2.0000;
527:   y[i]   = 17.7000;
528:   t[i++] = 2.5000;
529:   y[i]   = 24.5600;
530:   t[i++] = 2.0000;
531:   y[i]   = 18.6700;
532:   t[i++] = 2.5000;
533:   y[i]   = 16.2400;
534:   t[i++] = 3.0000;
535:   y[i]   = 8.7400;
536:   t[i++] = 4.0000;
537:   y[i]   = 7.8700;
538:   t[i++] = 5.0000;
539:   y[i]   = 8.5100;
540:   t[i++] = 6.0000;
541:   y[i]   = 66.7000;
542:   t[i++] = .5000;
543:   y[i]   = 59.2000;
544:   t[i++] = .7500;
545:   y[i]   = 40.8000;
546:   t[i++] = 1.0000;
547:   y[i]   = 30.7000;
548:   t[i++] = 1.5000;
549:   y[i]   = 25.7000;
550:   t[i++] = 2.0000;
551:   y[i]   = 16.3000;
552:   t[i++] = 2.5000;
553:   y[i]   = 25.9900;
554:   t[i++] = 2.0000;
555:   y[i]   = 16.9500;
556:   t[i++] = 2.5000;
557:   y[i]   = 13.3500;
558:   t[i++] = 3.0000;
559:   y[i]   = 8.6200;
560:   t[i++] = 4.0000;
561:   y[i]   = 7.2000;
562:   t[i++] = 5.0000;
563:   y[i]   = 6.6400;
564:   t[i++] = 6.0000;
565:   y[i]   = 13.6900;
566:   t[i++] = 3.0000;
567:   y[i]   = 81.0000;
568:   t[i++] = .5000;
569:   y[i]   = 64.5000;
570:   t[i++] = .7500;
571:   y[i]   = 35.5000;
572:   t[i++] = 1.5000;
573:   y[i]   = 13.3100;
574:   t[i++] = 3.0000;
575:   y[i]   = 4.8700;
576:   t[i++] = 6.0000;
577:   y[i]   = 12.9400;
578:   t[i++] = 3.0000;
579:   y[i]   = 5.0600;
580:   t[i++] = 6.0000;
581:   y[i]   = 15.1900;
582:   t[i++] = 3.0000;
583:   y[i]   = 14.6200;
584:   t[i++] = 3.0000;
585:   y[i]   = 15.6400;
586:   t[i++] = 3.0000;
587:   y[i]   = 25.5000;
588:   t[i++] = 1.7500;
589:   y[i]   = 25.9500;
590:   t[i++] = 1.7500;
591:   y[i]   = 81.7000;
592:   t[i++] = .5000;
593:   y[i]   = 61.6000;
594:   t[i++] = .7500;
595:   y[i]   = 29.8000;
596:   t[i++] = 1.7500;
597:   y[i]   = 29.8100;
598:   t[i++] = 1.7500;
599:   y[i]   = 17.1700;
600:   t[i++] = 2.7500;
601:   y[i]   = 10.3900;
602:   t[i++] = 3.7500;
603:   y[i]   = 28.4000;
604:   t[i++] = 1.7500;
605:   y[i]   = 28.6900;
606:   t[i++] = 1.7500;
607:   y[i]   = 81.3000;
608:   t[i++] = .5000;
609:   y[i]   = 60.9000;
610:   t[i++] = .7500;
611:   y[i]   = 16.6500;
612:   t[i++] = 2.7500;
613:   y[i]   = 10.0500;
614:   t[i++] = 3.7500;
615:   y[i]   = 28.9000;
616:   t[i++] = 1.7500;
617:   y[i]   = 28.9500;
618:   t[i++] = 1.7500;
619:   return 0;
620: }

622: /*TEST

624:      build:
625:        requires: !complex

627:      test:
628:        args:  -tao_smonitor -tao_max_it 100 -tao_type pounders -tao_pounders_delta 0.05 -tao_gatol 1.e-5
629:        requires: !single
630:        TODO: produces different output for many different systems

632:      test:
633:        suffix: 2
634:        args: -tao_smonitor -tao_max_it 100 -wtype 1 -tao_type pounders -tao_pounders_delta 0.05 -tao_gatol 1.e-5
635:        requires: !single
636:        TODO: produces different output for many different systems

638:      test:
639:        suffix: 3
640:        args: -tao_smonitor -tao_max_it 100 -wtype 2 -tao_type pounders -tao_pounders_delta 0.05 -tao_gatol 1.e-5
641:        requires: !single
642:        TODO: produces different output for many different systems

644:      test:
645:        suffix: 4
646:        args: -tao_smonitor -tao_max_it 100 -tao_type pounders -tao_pounders_delta 0.05 -pounders_subsolver_tao_type blmvm -tao_gatol 1.e-5
647:        requires: !single
648:        TODO: produces different output for many different systems

650:  TEST*/