Actual source code: chwirut1.c
1: /*
2: Include "tao.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 - sysem routines petscmat.h - matrices
6: petscis.h - index sets petscksp.h - Krylov subspace methods
7: petscviewer.h - viewers petscpc.h - preconditioners
9: */
11: #include "taosolver.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: */
26: static char help[]="Finds the nonlinear least-squares solution to the model \n\
27: y = exp[-b1*x]/(b2+b3*x) + e \n";
31: /*T
32: Concepts: TAO - Solving a system of nonlinear equations, nonlinear ;east squares
33: Routines: TaoInitialize(); TaoFinalize();
34: Routines: TaoCreate();
35: Routines: TaoSetType();
36: Routines: TaoSetSeparableObjectiveRoutine();
37: Routines: TaoSetJacobianRoutine();
38: Routines: TaoSetInitialVector();
39: Routines: TaoSetFromOptions();
40: Routines: TaoSetHistory(); TaoGetHistory();
41: Routines: TaoSolve();
42: Routines: TaoView(); TaoDestroy();
43: Processors: 1
44: T*/
46: #define NOBSERVATIONS 214
47: #define NPARAMETERS 3
49: /* User-defined application context */
50: typedef struct {
51: /* Working space */
52: PetscReal t[NOBSERVATIONS]; /* array of independent variables of observation */
53: PetscReal y[NOBSERVATIONS]; /* array of dependent variables */
54: PetscReal j[NOBSERVATIONS][NPARAMETERS]; /* dense jacobian matrix array*/
55: PetscInt idm[NOBSERVATIONS]; /* Matrix indices for jacobian */
56: PetscInt idn[NPARAMETERS];
57: } AppCtx;
59: /* User provided Routines */
60: PetscErrorCode InitializeData(AppCtx *user);
61: PetscErrorCode FormStartingPoint(Vec);
62: PetscErrorCode EvaluateFunction(TaoSolver, Vec, Vec, void *);
63: PetscErrorCode EvaluateJacobian(TaoSolver, Vec, Mat*, Mat*, MatStructure*,void *);
66: /*--------------------------------------------------------------------*/
69: int main(int argc,char **argv)
70: {
72: Vec x, f; /* solution, function */
73: Mat J; /* Jacobian matrix */
74: TaoSolver tao; /* TaoSolver solver context */
75: PetscInt i; /* iteration information */
76: PetscReal hist[100],resid[100];
77: PetscInt nhist;
78: PetscBool printhistory;
79: AppCtx user; /* user-defined work context */
81: /* Initialize TAO and PETSc */
82: PetscInitialize(&argc,&argv,(char *)0,help);
83: TaoInitialize(&argc,&argv,(char *)0,help);
85: printhistory = PETSC_FALSE;
86: PetscOptionsGetBool(PETSC_NULL,"-printhistory",&printhistory,0);
87: /* Allocate vectors */
88: VecCreateSeq(MPI_COMM_SELF,NPARAMETERS,&x);
89: VecCreateSeq(MPI_COMM_SELF,NOBSERVATIONS,&f);
91: /* Create the Jacobian matrix. */
92: MatCreateSeqDense(MPI_COMM_SELF,NOBSERVATIONS,NPARAMETERS,PETSC_NULL,&J);
93:
95: for (i=0;i<NOBSERVATIONS;i++)
96: user.idm[i] = i;
98: for (i=0;i<NPARAMETERS;i++)
99: user.idn[i] = i;
103: /* TAO code begins here */
105: /* Create TAO solver and set desired solution method */
106: TaoCreate(PETSC_COMM_SELF,&tao);
107: TaoSetType(tao,"tao_pounders");
109: /* Set the function and Jacobian routines. */
110: InitializeData(&user);
111: FormStartingPoint(x);
112: TaoSetInitialVector(tao,x);
113: TaoSetSeparableObjectiveRoutine(tao,f,EvaluateFunction,(void*)&user);
114: TaoSetJacobianRoutine(tao, J, J, EvaluateJacobian, (void*)&user);
117: /* Check for any TAO command line arguments */
118: TaoSetFromOptions(tao);
119:
120: TaoSetHistory(tao,hist,resid,0,100,PETSC_TRUE);
121: /* Perform the Solve */
122: TaoSolve(tao);
123: if (printhistory) {
124: TaoGetHistory(tao,0,0,0,&nhist);
125: for (i=0;i<nhist;i++) {
126: PetscPrintf(PETSC_COMM_WORLD,"%G\t%G\n",hist[i],resid[i]);
127: }
128: }
129: TaoView(tao,PETSC_VIEWER_STDOUT_SELF);
131: /* Free TAO data structures */
132: TaoDestroy(&tao);
134: /* Free PETSc data structures */
135: VecDestroy(&x);
136: VecDestroy(&f);
137: MatDestroy(&J);
140: /* Finalize TAO */
141: TaoFinalize();
142: PetscFinalize();
143:
144: return 0;
145: }
150: /*--------------------------------------------------------------------*/
153: PetscErrorCode EvaluateFunction(TaoSolver tao, Vec X, Vec F, void *ptr)
154: {
155: AppCtx *user = (AppCtx *)ptr;
156: PetscInt i;
157: PetscReal *y=user->y,*x,*f,*t=user->t;
161: VecGetArray(X,&x);
162: VecGetArray(F,&f);
164: for (i=0;i<NOBSERVATIONS;i++) {
165: f[i] = y[i] - PetscExpScalar(-x[0]*t[i])/(x[1] + x[2]*t[i]);
166: }
168:
170: VecRestoreArray(X,&x);
171: VecRestoreArray(F,&f);
173: PetscLogFlops(6*NOBSERVATIONS);
174: return(0);
175: }
177: /*------------------------------------------------------------*/
178: /* J[i][j] = df[i]/dt[j] */
181: PetscErrorCode EvaluateJacobian(TaoSolver tao, Vec X, Mat *J, Mat *Jpre, MatStructure*matstruct,void *ptr)
182: {
183: AppCtx *user = (AppCtx *)ptr;
184: PetscInt i;
185: PetscReal *x,*t=user->t;
186: PetscReal base;
192: /* Get handles to the Vectors */
193: VecGetArray(X,&x);
197: for (i=0;i<NOBSERVATIONS;i++) {
198: base = PetscExpScalar(-x[0]*t[i])/(x[1] + x[2]*t[i]);
200: user->j[i][0] = t[i]*base;
201: user->j[i][1] = base/(x[1] + x[2]*t[i]);
202: user->j[i][2] = base*t[i]/(x[1] + x[2]*t[i]);
204: }
206: /* Assemble the matrix */
207: MatSetValues(*J,NOBSERVATIONS,user->idm, NPARAMETERS, user->idn,(PetscReal *)user->j,
208: INSERT_VALUES);
209: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
210: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
211:
212: /* Return the handles */
213: VecRestoreArray(X,&x);
215: PetscLogFlops(NOBSERVATIONS * 13);
217: return(0);
218: }
220: /* ------------------------------------------------------------ */
223: PetscErrorCode FormStartingPoint(Vec X)
224: {
225: PetscReal *x;
227:
230: VecGetArray(X,&x);
232: x[0] = 0.15;
233: x[1] = 0.008;
234: x[2] = 0.010;
236: VecRestoreArray(X,&x);
238: return(0);
239: }
242: /* ---------------------------------------------------------------------- */
245: PetscErrorCode InitializeData(AppCtx *user)
246: {
247: PetscReal *t=user->t,*y=user->y;
248: PetscInt i=0;
253: y[i] = 92.9000; t[i++] = 0.5000;
254: y[i] = 78.7000; t[i++] = 0.6250;
255: y[i] = 64.2000; t[i++] = 0.7500;
256: y[i] = 64.9000; t[i++] = 0.8750;
257: y[i] = 57.1000; t[i++] = 1.0000;
258: y[i] = 43.3000; t[i++] = 1.2500;
259: y[i] = 31.1000; t[i++] = 1.7500;
260: y[i] = 23.6000; t[i++] = 2.2500;
261: y[i] = 31.0500; t[i++] = 1.7500;
262: y[i] = 23.7750; t[i++] = 2.2500;
263: y[i] = 17.7375; t[i++] = 2.7500;
264: y[i] = 13.8000; t[i++] = 3.2500;
265: y[i] = 11.5875; t[i++] = 3.7500;
266: y[i] = 9.4125; t[i++] = 4.2500;
267: y[i] = 7.7250; t[i++] = 4.7500;
268: y[i] = 7.3500; t[i++] = 5.2500;
269: y[i] = 8.0250; t[i++] = 5.7500;
270: y[i] = 90.6000; t[i++] = 0.5000;
271: y[i] = 76.9000; t[i++] = 0.6250;
272: y[i] = 71.6000; t[i++] = 0.7500;
273: y[i] = 63.6000; t[i++] = 0.8750;
274: y[i] = 54.0000; t[i++] = 1.0000;
275: y[i] = 39.2000; t[i++] = 1.2500;
276: y[i] = 29.3000; t[i++] = 1.7500;
277: y[i] = 21.4000; t[i++] = 2.2500;
278: y[i] = 29.1750; t[i++] = 1.7500;
279: y[i] = 22.1250; t[i++] = 2.2500;
280: y[i] = 17.5125; t[i++] = 2.7500;
281: y[i] = 14.2500; t[i++] = 3.2500;
282: y[i] = 9.4500; t[i++] = 3.7500;
283: y[i] = 9.1500; t[i++] = 4.2500;
284: y[i] = 7.9125; t[i++] = 4.7500;
285: y[i] = 8.4750; t[i++] = 5.2500;
286: y[i] = 6.1125; t[i++] = 5.7500;
287: y[i] = 80.0000; t[i++] = 0.5000;
288: y[i] = 79.0000; t[i++] = 0.6250;
289: y[i] = 63.8000; t[i++] = 0.7500;
290: y[i] = 57.2000; t[i++] = 0.8750;
291: y[i] = 53.2000; t[i++] = 1.0000;
292: y[i] = 42.5000; t[i++] = 1.2500;
293: y[i] = 26.8000; t[i++] = 1.7500;
294: y[i] = 20.4000; t[i++] = 2.2500;
295: y[i] = 26.8500; t[i++] = 1.7500;
296: y[i] = 21.0000; t[i++] = 2.2500;
297: y[i] = 16.4625; t[i++] = 2.7500;
298: y[i] = 12.5250; t[i++] = 3.2500;
299: y[i] = 10.5375; t[i++] = 3.7500;
300: y[i] = 8.5875; t[i++] = 4.2500;
301: y[i] = 7.1250; t[i++] = 4.7500;
302: y[i] = 6.1125; t[i++] = 5.2500;
303: y[i] = 5.9625; t[i++] = 5.7500;
304: y[i] = 74.1000; t[i++] = 0.5000;
305: y[i] = 67.3000; t[i++] = 0.6250;
306: y[i] = 60.8000; t[i++] = 0.7500;
307: y[i] = 55.5000; t[i++] = 0.8750;
308: y[i] = 50.3000; t[i++] = 1.0000;
309: y[i] = 41.0000; t[i++] = 1.2500;
310: y[i] = 29.4000; t[i++] = 1.7500;
311: y[i] = 20.4000; t[i++] = 2.2500;
312: y[i] = 29.3625; t[i++] = 1.7500;
313: y[i] = 21.1500; t[i++] = 2.2500;
314: y[i] = 16.7625; t[i++] = 2.7500;
315: y[i] = 13.2000; t[i++] = 3.2500;
316: y[i] = 10.8750; t[i++] = 3.7500;
317: y[i] = 8.1750; t[i++] = 4.2500;
318: y[i] = 7.3500; t[i++] = 4.7500;
319: y[i] = 5.9625; t[i++] = 5.2500;
320: y[i] = 5.6250; t[i++] = 5.7500;
321: y[i] = 81.5000; t[i++] = .5000;
322: y[i] = 62.4000; t[i++] = .7500;
323: y[i] = 32.5000; t[i++] = 1.5000;
324: y[i] = 12.4100; t[i++] = 3.0000;
325: y[i] = 13.1200; t[i++] = 3.0000;
326: y[i] = 15.5600; t[i++] = 3.0000;
327: y[i] = 5.6300; t[i++] = 6.0000;
328: y[i] = 78.0000; t[i++] = .5000;
329: y[i] = 59.9000; t[i++] = .7500;
330: y[i] = 33.2000; t[i++] = 1.5000;
331: y[i] = 13.8400; t[i++] = 3.0000;
332: y[i] = 12.7500; t[i++] = 3.0000;
333: y[i] = 14.6200; t[i++] = 3.0000;
334: y[i] = 3.9400; t[i++] = 6.0000;
335: y[i] = 76.8000; t[i++] = .5000;
336: y[i] = 61.0000; t[i++] = .7500;
337: y[i] = 32.9000; t[i++] = 1.5000;
338: y[i] = 13.8700; t[i++] = 3.0000;
339: y[i] = 11.8100; t[i++] = 3.0000;
340: y[i] = 13.3100; t[i++] = 3.0000;
341: y[i] = 5.4400; t[i++] = 6.0000;
342: y[i] = 78.0000; t[i++] = .5000;
343: y[i] = 63.5000; t[i++] = .7500;
344: y[i] = 33.8000; t[i++] = 1.5000;
345: y[i] = 12.5600; t[i++] = 3.0000;
346: y[i] = 5.6300; t[i++] = 6.0000;
347: y[i] = 12.7500; t[i++] = 3.0000;
348: y[i] = 13.1200; t[i++] = 3.0000;
349: y[i] = 5.4400; t[i++] = 6.0000;
350: y[i] = 76.8000; t[i++] = .5000;
351: y[i] = 60.0000; t[i++] = .7500;
352: y[i] = 47.8000; t[i++] = 1.0000;
353: y[i] = 32.0000; t[i++] = 1.5000;
354: y[i] = 22.2000; t[i++] = 2.0000;
355: y[i] = 22.5700; t[i++] = 2.0000;
356: y[i] = 18.8200; t[i++] = 2.5000;
357: y[i] = 13.9500; t[i++] = 3.0000;
358: y[i] = 11.2500; t[i++] = 4.0000;
359: y[i] = 9.0000; t[i++] = 5.0000;
360: y[i] = 6.6700; t[i++] = 6.0000;
361: y[i] = 75.8000; t[i++] = .5000;
362: y[i] = 62.0000; t[i++] = .7500;
363: y[i] = 48.8000; t[i++] = 1.0000;
364: y[i] = 35.2000; t[i++] = 1.5000;
365: y[i] = 20.0000; t[i++] = 2.0000;
366: y[i] = 20.3200; t[i++] = 2.0000;
367: y[i] = 19.3100; t[i++] = 2.5000;
368: y[i] = 12.7500; t[i++] = 3.0000;
369: y[i] = 10.4200; t[i++] = 4.0000;
370: y[i] = 7.3100; t[i++] = 5.0000;
371: y[i] = 7.4200; t[i++] = 6.0000;
372: y[i] = 70.5000; t[i++] = .5000;
373: y[i] = 59.5000; t[i++] = .7500;
374: y[i] = 48.5000; t[i++] = 1.0000;
375: y[i] = 35.8000; t[i++] = 1.5000;
376: y[i] = 21.0000; t[i++] = 2.0000;
377: y[i] = 21.6700; t[i++] = 2.0000;
378: y[i] = 21.0000; t[i++] = 2.5000;
379: y[i] = 15.6400; t[i++] = 3.0000;
380: y[i] = 8.1700; t[i++] = 4.0000;
381: y[i] = 8.5500; t[i++] = 5.0000;
382: y[i] = 10.1200; t[i++] = 6.0000;
383: y[i] = 78.0000; t[i++] = .5000;
384: y[i] = 66.0000; t[i++] = .6250;
385: y[i] = 62.0000; t[i++] = .7500;
386: y[i] = 58.0000; t[i++] = .8750;
387: y[i] = 47.7000; t[i++] = 1.0000;
388: y[i] = 37.8000; t[i++] = 1.2500;
389: y[i] = 20.2000; t[i++] = 2.2500;
390: y[i] = 21.0700; t[i++] = 2.2500;
391: y[i] = 13.8700; t[i++] = 2.7500;
392: y[i] = 9.6700; t[i++] = 3.2500;
393: y[i] = 7.7600; t[i++] = 3.7500;
394: y[i] = 5.4400; t[i++] = 4.2500;
395: y[i] = 4.8700; t[i++] = 4.7500;
396: y[i] = 4.0100; t[i++] = 5.2500;
397: y[i] = 3.7500; t[i++] = 5.7500;
398: y[i] = 24.1900; t[i++] = 3.0000;
399: y[i] = 25.7600; t[i++] = 3.0000;
400: y[i] = 18.0700; t[i++] = 3.0000;
401: y[i] = 11.8100; t[i++] = 3.0000;
402: y[i] = 12.0700; t[i++] = 3.0000;
403: y[i] = 16.1200; t[i++] = 3.0000;
404: y[i] = 70.8000; t[i++] = .5000;
405: y[i] = 54.7000; t[i++] = .7500;
406: y[i] = 48.0000; t[i++] = 1.0000;
407: y[i] = 39.8000; t[i++] = 1.5000;
408: y[i] = 29.8000; t[i++] = 2.0000;
409: y[i] = 23.7000; t[i++] = 2.5000;
410: y[i] = 29.6200; t[i++] = 2.0000;
411: y[i] = 23.8100; t[i++] = 2.5000;
412: y[i] = 17.7000; t[i++] = 3.0000;
413: y[i] = 11.5500; t[i++] = 4.0000;
414: y[i] = 12.0700; t[i++] = 5.0000;
415: y[i] = 8.7400; t[i++] = 6.0000;
416: y[i] = 80.7000; t[i++] = .5000;
417: y[i] = 61.3000; t[i++] = .7500;
418: y[i] = 47.5000; t[i++] = 1.0000;
419: y[i] = 29.0000; t[i++] = 1.5000;
420: y[i] = 24.0000; t[i++] = 2.0000;
421: y[i] = 17.7000; t[i++] = 2.5000;
422: y[i] = 24.5600; t[i++] = 2.0000;
423: y[i] = 18.6700; t[i++] = 2.5000;
424: y[i] = 16.2400; t[i++] = 3.0000;
425: y[i] = 8.7400; t[i++] = 4.0000;
426: y[i] = 7.8700; t[i++] = 5.0000;
427: y[i] = 8.5100; t[i++] = 6.0000;
428: y[i] = 66.7000; t[i++] = .5000;
429: y[i] = 59.2000; t[i++] = .7500;
430: y[i] = 40.8000; t[i++] = 1.0000;
431: y[i] = 30.7000; t[i++] = 1.5000;
432: y[i] = 25.7000; t[i++] = 2.0000;
433: y[i] = 16.3000; t[i++] = 2.5000;
434: y[i] = 25.9900; t[i++] = 2.0000;
435: y[i] = 16.9500; t[i++] = 2.5000;
436: y[i] = 13.3500; t[i++] = 3.0000;
437: y[i] = 8.6200; t[i++] = 4.0000;
438: y[i] = 7.2000; t[i++] = 5.0000;
439: y[i] = 6.6400; t[i++] = 6.0000;
440: y[i] = 13.6900; t[i++] = 3.0000;
441: y[i] = 81.0000; t[i++] = .5000;
442: y[i] = 64.5000; t[i++] = .7500;
443: y[i] = 35.5000; t[i++] = 1.5000;
444: y[i] = 13.3100; t[i++] = 3.0000;
445: y[i] = 4.8700; t[i++] = 6.0000;
446: y[i] = 12.9400; t[i++] = 3.0000;
447: y[i] = 5.0600; t[i++] = 6.0000;
448: y[i] = 15.1900; t[i++] = 3.0000;
449: y[i] = 14.6200; t[i++] = 3.0000;
450: y[i] = 15.6400; t[i++] = 3.0000;
451: y[i] = 25.5000; t[i++] = 1.7500;
452: y[i] = 25.9500; t[i++] = 1.7500;
453: y[i] = 81.7000; t[i++] = .5000;
454: y[i] = 61.6000; t[i++] = .7500;
455: y[i] = 29.8000; t[i++] = 1.7500;
456: y[i] = 29.8100; t[i++] = 1.7500;
457: y[i] = 17.1700; t[i++] = 2.7500;
458: y[i] = 10.3900; t[i++] = 3.7500;
459: y[i] = 28.4000; t[i++] = 1.7500;
460: y[i] = 28.6900; t[i++] = 1.7500;
461: y[i] = 81.3000; t[i++] = .5000;
462: y[i] = 60.9000; t[i++] = .7500;
463: y[i] = 16.6500; t[i++] = 2.7500;
464: y[i] = 10.0500; t[i++] = 3.7500;
465: y[i] = 28.9000; t[i++] = 1.7500;
466: y[i] = 28.9500; t[i++] = 1.7500;
468: return(0);
469: }