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 */
63: PetscFunctionBeginUser;
64: PetscCall(PetscInitialize(&argc, &argv, NULL, help));
65: PetscCall(PetscOptionsGetInt(NULL, NULL, "-wtype", &wtype, &flg));
66: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "wtype=%" PetscInt_FMT "\n", wtype));
67: /* Allocate vectors */
68: PetscCall(VecCreateSeq(MPI_COMM_SELF, NPARAMETERS, &x));
69: PetscCall(VecCreateSeq(MPI_COMM_SELF, NOBSERVATIONS, &f));
71: PetscCall(VecDuplicate(f, &w));
73: /* no correlation, but set in different ways */
74: PetscCall(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: PetscCall(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: PetscCall(TaoCreate(PETSC_COMM_SELF, &tao));
90: PetscCall(TaoSetType(tao, TAOPOUNDERS));
92: /* Set the function and Jacobian routines. */
93: PetscCall(InitializeData(&user));
94: PetscCall(FormStartingPoint(x));
95: PetscCall(TaoSetSolution(tao, x));
96: PetscCall(TaoSetResidualRoutine(tao, f, EvaluateFunction, (void *)&user));
97: if (wtype == 1) PetscCall(TaoSetResidualWeights(tao, w, 0, NULL, NULL, NULL));
98: else if (wtype == 2) PetscCall(TaoSetResidualWeights(tao, NULL, NOBSERVATIONS, w_row, w_col, w_vals));
99: PetscCall(TaoSetJacobianResidualRoutine(tao, J, J, EvaluateJacobian, (void *)&user));
100: PetscCall(TaoSetTolerances(tao, 1e-5, 0.0, PETSC_CURRENT));
102: /* Check for any TAO command line arguments */
103: PetscCall(TaoSetFromOptions(tao));
105: PetscCall(TaoSetConvergenceHistory(tao, hist, resid, 0, lits, 100, PETSC_TRUE));
106: /* Perform the Solve */
107: PetscCall(TaoSolve(tao));
109: /* Free TAO data structures */
110: PetscCall(TaoDestroy(&tao));
112: /* Free PETSc data structures */
113: PetscCall(VecDestroy(&x));
114: PetscCall(VecDestroy(&w));
115: PetscCall(VecDestroy(&f));
116: PetscCall(MatDestroy(&J));
118: PetscCall(PetscFinalize());
119: return 0;
120: }
122: /*--------------------------------------------------------------------*/
123: PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr)
124: {
125: AppCtx *user = (AppCtx *)ptr;
126: PetscInt i;
127: PetscReal *y = user->y, *f, *t = user->t;
128: const PetscReal *x;
130: PetscFunctionBegin;
131: PetscCall(VecGetArrayRead(X, &x));
132: PetscCall(VecGetArray(F, &f));
134: for (i = 0; i < NOBSERVATIONS; i++) f[i] = y[i] - PetscExpScalar(-x[0] * t[i]) / (x[1] + x[2] * t[i]);
135: PetscCall(VecRestoreArrayRead(X, &x));
136: PetscCall(VecRestoreArray(F, &f));
137: PetscCall(PetscLogFlops(6 * NOBSERVATIONS));
138: PetscFunctionReturn(PETSC_SUCCESS);
139: }
141: /*------------------------------------------------------------*/
142: /* J[i][j] = df[i]/dt[j] */
143: PetscErrorCode EvaluateJacobian(Tao tao, Vec X, Mat J, Mat Jpre, void *ptr)
144: {
145: AppCtx *user = (AppCtx *)ptr;
146: PetscInt i;
147: PetscReal *t = user->t;
148: const PetscReal *x;
149: PetscReal base;
151: PetscFunctionBegin;
152: PetscCall(VecGetArrayRead(X, &x));
153: for (i = 0; i < NOBSERVATIONS; i++) {
154: base = PetscExpScalar(-x[0] * t[i]) / (x[1] + x[2] * t[i]);
156: user->j[i][0] = t[i] * base;
157: user->j[i][1] = base / (x[1] + x[2] * t[i]);
158: user->j[i][2] = base * t[i] / (x[1] + x[2] * t[i]);
159: }
161: /* Assemble the matrix */
162: PetscCall(MatSetValues(J, NOBSERVATIONS, user->idm, NPARAMETERS, user->idn, (PetscReal *)user->j, INSERT_VALUES));
163: PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
164: PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
166: PetscCall(VecRestoreArrayRead(X, &x));
167: PetscCall(PetscLogFlops(NOBSERVATIONS * 13));
168: PetscFunctionReturn(PETSC_SUCCESS);
169: }
171: /* ------------------------------------------------------------ */
172: PetscErrorCode FormStartingPoint(Vec X)
173: {
174: PetscReal *x;
176: PetscFunctionBegin;
177: PetscCall(VecGetArray(X, &x));
178: x[0] = 1.19;
179: x[1] = -1.86;
180: x[2] = 1.08;
181: PetscCall(VecRestoreArray(X, &x));
182: PetscFunctionReturn(PETSC_SUCCESS);
183: }
185: /* ---------------------------------------------------------------------- */
186: PetscErrorCode InitializeData(AppCtx *user)
187: {
188: PetscReal *t = user->t, *y = user->y;
189: PetscInt i = 0;
191: PetscFunctionBegin;
192: y[i] = 92.9000;
193: t[i++] = 0.5000;
194: y[i] = 78.7000;
195: t[i++] = 0.6250;
196: y[i] = 64.2000;
197: t[i++] = 0.7500;
198: y[i] = 64.9000;
199: t[i++] = 0.8750;
200: y[i] = 57.1000;
201: t[i++] = 1.0000;
202: y[i] = 43.3000;
203: t[i++] = 1.2500;
204: y[i] = 31.1000;
205: t[i++] = 1.7500;
206: y[i] = 23.6000;
207: t[i++] = 2.2500;
208: y[i] = 31.0500;
209: t[i++] = 1.7500;
210: y[i] = 23.7750;
211: t[i++] = 2.2500;
212: y[i] = 17.7375;
213: t[i++] = 2.7500;
214: y[i] = 13.8000;
215: t[i++] = 3.2500;
216: y[i] = 11.5875;
217: t[i++] = 3.7500;
218: y[i] = 9.4125;
219: t[i++] = 4.2500;
220: y[i] = 7.7250;
221: t[i++] = 4.7500;
222: y[i] = 7.3500;
223: t[i++] = 5.2500;
224: y[i] = 8.0250;
225: t[i++] = 5.7500;
226: y[i] = 90.6000;
227: t[i++] = 0.5000;
228: y[i] = 76.9000;
229: t[i++] = 0.6250;
230: y[i] = 71.6000;
231: t[i++] = 0.7500;
232: y[i] = 63.6000;
233: t[i++] = 0.8750;
234: y[i] = 54.0000;
235: t[i++] = 1.0000;
236: y[i] = 39.2000;
237: t[i++] = 1.2500;
238: y[i] = 29.3000;
239: t[i++] = 1.7500;
240: y[i] = 21.4000;
241: t[i++] = 2.2500;
242: y[i] = 29.1750;
243: t[i++] = 1.7500;
244: y[i] = 22.1250;
245: t[i++] = 2.2500;
246: y[i] = 17.5125;
247: t[i++] = 2.7500;
248: y[i] = 14.2500;
249: t[i++] = 3.2500;
250: y[i] = 9.4500;
251: t[i++] = 3.7500;
252: y[i] = 9.1500;
253: t[i++] = 4.2500;
254: y[i] = 7.9125;
255: t[i++] = 4.7500;
256: y[i] = 8.4750;
257: t[i++] = 5.2500;
258: y[i] = 6.1125;
259: t[i++] = 5.7500;
260: y[i] = 80.0000;
261: t[i++] = 0.5000;
262: y[i] = 79.0000;
263: t[i++] = 0.6250;
264: y[i] = 63.8000;
265: t[i++] = 0.7500;
266: y[i] = 57.2000;
267: t[i++] = 0.8750;
268: y[i] = 53.2000;
269: t[i++] = 1.0000;
270: y[i] = 42.5000;
271: t[i++] = 1.2500;
272: y[i] = 26.8000;
273: t[i++] = 1.7500;
274: y[i] = 20.4000;
275: t[i++] = 2.2500;
276: y[i] = 26.8500;
277: t[i++] = 1.7500;
278: y[i] = 21.0000;
279: t[i++] = 2.2500;
280: y[i] = 16.4625;
281: t[i++] = 2.7500;
282: y[i] = 12.5250;
283: t[i++] = 3.2500;
284: y[i] = 10.5375;
285: t[i++] = 3.7500;
286: y[i] = 8.5875;
287: t[i++] = 4.2500;
288: y[i] = 7.1250;
289: t[i++] = 4.7500;
290: y[i] = 6.1125;
291: t[i++] = 5.2500;
292: y[i] = 5.9625;
293: t[i++] = 5.7500;
294: y[i] = 74.1000;
295: t[i++] = 0.5000;
296: y[i] = 67.3000;
297: t[i++] = 0.6250;
298: y[i] = 60.8000;
299: t[i++] = 0.7500;
300: y[i] = 55.5000;
301: t[i++] = 0.8750;
302: y[i] = 50.3000;
303: t[i++] = 1.0000;
304: y[i] = 41.0000;
305: t[i++] = 1.2500;
306: y[i] = 29.4000;
307: t[i++] = 1.7500;
308: y[i] = 20.4000;
309: t[i++] = 2.2500;
310: y[i] = 29.3625;
311: t[i++] = 1.7500;
312: y[i] = 21.1500;
313: t[i++] = 2.2500;
314: y[i] = 16.7625;
315: t[i++] = 2.7500;
316: y[i] = 13.2000;
317: t[i++] = 3.2500;
318: y[i] = 10.8750;
319: t[i++] = 3.7500;
320: y[i] = 8.1750;
321: t[i++] = 4.2500;
322: y[i] = 7.3500;
323: t[i++] = 4.7500;
324: y[i] = 5.9625;
325: t[i++] = 5.2500;
326: y[i] = 5.6250;
327: t[i++] = 5.7500;
328: y[i] = 81.5000;
329: t[i++] = .5000;
330: y[i] = 62.4000;
331: t[i++] = .7500;
332: y[i] = 32.5000;
333: t[i++] = 1.5000;
334: y[i] = 12.4100;
335: t[i++] = 3.0000;
336: y[i] = 13.1200;
337: t[i++] = 3.0000;
338: y[i] = 15.5600;
339: t[i++] = 3.0000;
340: y[i] = 5.6300;
341: t[i++] = 6.0000;
342: y[i] = 78.0000;
343: t[i++] = .5000;
344: y[i] = 59.9000;
345: t[i++] = .7500;
346: y[i] = 33.2000;
347: t[i++] = 1.5000;
348: y[i] = 13.8400;
349: t[i++] = 3.0000;
350: y[i] = 12.7500;
351: t[i++] = 3.0000;
352: y[i] = 14.6200;
353: t[i++] = 3.0000;
354: y[i] = 3.9400;
355: t[i++] = 6.0000;
356: y[i] = 76.8000;
357: t[i++] = .5000;
358: y[i] = 61.0000;
359: t[i++] = .7500;
360: y[i] = 32.9000;
361: t[i++] = 1.5000;
362: y[i] = 13.8700;
363: t[i++] = 3.0000;
364: y[i] = 11.8100;
365: t[i++] = 3.0000;
366: y[i] = 13.3100;
367: t[i++] = 3.0000;
368: y[i] = 5.4400;
369: t[i++] = 6.0000;
370: y[i] = 78.0000;
371: t[i++] = .5000;
372: y[i] = 63.5000;
373: t[i++] = .7500;
374: y[i] = 33.8000;
375: t[i++] = 1.5000;
376: y[i] = 12.5600;
377: t[i++] = 3.0000;
378: y[i] = 5.6300;
379: t[i++] = 6.0000;
380: y[i] = 12.7500;
381: t[i++] = 3.0000;
382: y[i] = 13.1200;
383: t[i++] = 3.0000;
384: y[i] = 5.4400;
385: t[i++] = 6.0000;
386: y[i] = 76.8000;
387: t[i++] = .5000;
388: y[i] = 60.0000;
389: t[i++] = .7500;
390: y[i] = 47.8000;
391: t[i++] = 1.0000;
392: y[i] = 32.0000;
393: t[i++] = 1.5000;
394: y[i] = 22.2000;
395: t[i++] = 2.0000;
396: y[i] = 22.5700;
397: t[i++] = 2.0000;
398: y[i] = 18.8200;
399: t[i++] = 2.5000;
400: y[i] = 13.9500;
401: t[i++] = 3.0000;
402: y[i] = 11.2500;
403: t[i++] = 4.0000;
404: y[i] = 9.0000;
405: t[i++] = 5.0000;
406: y[i] = 6.6700;
407: t[i++] = 6.0000;
408: y[i] = 75.8000;
409: t[i++] = .5000;
410: y[i] = 62.0000;
411: t[i++] = .7500;
412: y[i] = 48.8000;
413: t[i++] = 1.0000;
414: y[i] = 35.2000;
415: t[i++] = 1.5000;
416: y[i] = 20.0000;
417: t[i++] = 2.0000;
418: y[i] = 20.3200;
419: t[i++] = 2.0000;
420: y[i] = 19.3100;
421: t[i++] = 2.5000;
422: y[i] = 12.7500;
423: t[i++] = 3.0000;
424: y[i] = 10.4200;
425: t[i++] = 4.0000;
426: y[i] = 7.3100;
427: t[i++] = 5.0000;
428: y[i] = 7.4200;
429: t[i++] = 6.0000;
430: y[i] = 70.5000;
431: t[i++] = .5000;
432: y[i] = 59.5000;
433: t[i++] = .7500;
434: y[i] = 48.5000;
435: t[i++] = 1.0000;
436: y[i] = 35.8000;
437: t[i++] = 1.5000;
438: y[i] = 21.0000;
439: t[i++] = 2.0000;
440: y[i] = 21.6700;
441: t[i++] = 2.0000;
442: y[i] = 21.0000;
443: t[i++] = 2.5000;
444: y[i] = 15.6400;
445: t[i++] = 3.0000;
446: y[i] = 8.1700;
447: t[i++] = 4.0000;
448: y[i] = 8.5500;
449: t[i++] = 5.0000;
450: y[i] = 10.1200;
451: t[i++] = 6.0000;
452: y[i] = 78.0000;
453: t[i++] = .5000;
454: y[i] = 66.0000;
455: t[i++] = .6250;
456: y[i] = 62.0000;
457: t[i++] = .7500;
458: y[i] = 58.0000;
459: t[i++] = .8750;
460: y[i] = 47.7000;
461: t[i++] = 1.0000;
462: y[i] = 37.8000;
463: t[i++] = 1.2500;
464: y[i] = 20.2000;
465: t[i++] = 2.2500;
466: y[i] = 21.0700;
467: t[i++] = 2.2500;
468: y[i] = 13.8700;
469: t[i++] = 2.7500;
470: y[i] = 9.6700;
471: t[i++] = 3.2500;
472: y[i] = 7.7600;
473: t[i++] = 3.7500;
474: y[i] = 5.4400;
475: t[i++] = 4.2500;
476: y[i] = 4.8700;
477: t[i++] = 4.7500;
478: y[i] = 4.0100;
479: t[i++] = 5.2500;
480: y[i] = 3.7500;
481: t[i++] = 5.7500;
482: y[i] = 24.1900;
483: t[i++] = 3.0000;
484: y[i] = 25.7600;
485: t[i++] = 3.0000;
486: y[i] = 18.0700;
487: t[i++] = 3.0000;
488: y[i] = 11.8100;
489: t[i++] = 3.0000;
490: y[i] = 12.0700;
491: t[i++] = 3.0000;
492: y[i] = 16.1200;
493: t[i++] = 3.0000;
494: y[i] = 70.8000;
495: t[i++] = .5000;
496: y[i] = 54.7000;
497: t[i++] = .7500;
498: y[i] = 48.0000;
499: t[i++] = 1.0000;
500: y[i] = 39.8000;
501: t[i++] = 1.5000;
502: y[i] = 29.8000;
503: t[i++] = 2.0000;
504: y[i] = 23.7000;
505: t[i++] = 2.5000;
506: y[i] = 29.6200;
507: t[i++] = 2.0000;
508: y[i] = 23.8100;
509: t[i++] = 2.5000;
510: y[i] = 17.7000;
511: t[i++] = 3.0000;
512: y[i] = 11.5500;
513: t[i++] = 4.0000;
514: y[i] = 12.0700;
515: t[i++] = 5.0000;
516: y[i] = 8.7400;
517: t[i++] = 6.0000;
518: y[i] = 80.7000;
519: t[i++] = .5000;
520: y[i] = 61.3000;
521: t[i++] = .7500;
522: y[i] = 47.5000;
523: t[i++] = 1.0000;
524: y[i] = 29.0000;
525: t[i++] = 1.5000;
526: y[i] = 24.0000;
527: t[i++] = 2.0000;
528: y[i] = 17.7000;
529: t[i++] = 2.5000;
530: y[i] = 24.5600;
531: t[i++] = 2.0000;
532: y[i] = 18.6700;
533: t[i++] = 2.5000;
534: y[i] = 16.2400;
535: t[i++] = 3.0000;
536: y[i] = 8.7400;
537: t[i++] = 4.0000;
538: y[i] = 7.8700;
539: t[i++] = 5.0000;
540: y[i] = 8.5100;
541: t[i++] = 6.0000;
542: y[i] = 66.7000;
543: t[i++] = .5000;
544: y[i] = 59.2000;
545: t[i++] = .7500;
546: y[i] = 40.8000;
547: t[i++] = 1.0000;
548: y[i] = 30.7000;
549: t[i++] = 1.5000;
550: y[i] = 25.7000;
551: t[i++] = 2.0000;
552: y[i] = 16.3000;
553: t[i++] = 2.5000;
554: y[i] = 25.9900;
555: t[i++] = 2.0000;
556: y[i] = 16.9500;
557: t[i++] = 2.5000;
558: y[i] = 13.3500;
559: t[i++] = 3.0000;
560: y[i] = 8.6200;
561: t[i++] = 4.0000;
562: y[i] = 7.2000;
563: t[i++] = 5.0000;
564: y[i] = 6.6400;
565: t[i++] = 6.0000;
566: y[i] = 13.6900;
567: t[i++] = 3.0000;
568: y[i] = 81.0000;
569: t[i++] = .5000;
570: y[i] = 64.5000;
571: t[i++] = .7500;
572: y[i] = 35.5000;
573: t[i++] = 1.5000;
574: y[i] = 13.3100;
575: t[i++] = 3.0000;
576: y[i] = 4.8700;
577: t[i++] = 6.0000;
578: y[i] = 12.9400;
579: t[i++] = 3.0000;
580: y[i] = 5.0600;
581: t[i++] = 6.0000;
582: y[i] = 15.1900;
583: t[i++] = 3.0000;
584: y[i] = 14.6200;
585: t[i++] = 3.0000;
586: y[i] = 15.6400;
587: t[i++] = 3.0000;
588: y[i] = 25.5000;
589: t[i++] = 1.7500;
590: y[i] = 25.9500;
591: t[i++] = 1.7500;
592: y[i] = 81.7000;
593: t[i++] = .5000;
594: y[i] = 61.6000;
595: t[i++] = .7500;
596: y[i] = 29.8000;
597: t[i++] = 1.7500;
598: y[i] = 29.8100;
599: t[i++] = 1.7500;
600: y[i] = 17.1700;
601: t[i++] = 2.7500;
602: y[i] = 10.3900;
603: t[i++] = 3.7500;
604: y[i] = 28.4000;
605: t[i++] = 1.7500;
606: y[i] = 28.6900;
607: t[i++] = 1.7500;
608: y[i] = 81.3000;
609: t[i++] = .5000;
610: y[i] = 60.9000;
611: t[i++] = .7500;
612: y[i] = 16.6500;
613: t[i++] = 2.7500;
614: y[i] = 10.0500;
615: t[i++] = 3.7500;
616: y[i] = 28.9000;
617: t[i++] = 1.7500;
618: y[i] = 28.9500;
619: t[i++] = 1.7500;
620: PetscFunctionReturn(PETSC_SUCCESS);
621: }
623: /*TEST
625: build:
626: requires: !complex
628: test:
629: args: -tao_monitor_short -tao_max_it 100 -tao_type pounders -tao_pounders_delta 0.05 -tao_gatol 1.e-5
630: requires: !single
631: TODO: produces different output for many different systems
633: test:
634: suffix: 2
635: args: -tao_monitor_short -tao_max_it 100 -wtype 1 -tao_type pounders -tao_pounders_delta 0.05 -tao_gatol 1.e-5
636: requires: !single
637: TODO: produces different output for many different systems
639: test:
640: suffix: 3
641: args: -tao_monitor_short -tao_max_it 100 -wtype 2 -tao_type pounders -tao_pounders_delta 0.05 -tao_gatol 1.e-5
642: requires: !single
643: TODO: produces different output for many different systems
645: test:
646: suffix: 4
647: args: -tao_monitor_short -tao_max_it 100 -tao_type pounders -tao_pounders_delta 0.05 -pounders_subsolver_tao_type blmvm -tao_gatol 1.e-5
648: requires: !single
649: TODO: produces different output for many different systems
651: TEST*/