Actual source code: chwirut1.c
petsc-3.6.4 2016-04-12
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 - 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 <petsctao.h>
13: /*
14: Description: These data are the result of a NIST study involving
15: ultrasonic calibration. The response variable is
16: ultrasonic response, and the predictor variable is
17: metal distance.
19: Reference: Chwirut, D., NIST (197?).
20: Ultrasonic Reference Block Study.
21: */
25: static char help[]="Finds the nonlinear least-squares solution to the model \n\
26: y = exp[-b1*x]/(b2+b3*x) + e \n";
28: /*T
29: Concepts: TAO^Solving a system of nonlinear equations, nonlinear least squares
30: Routines: TaoCreate();
31: Routines: TaoSetType();
32: Routines: TaoSetSeparableObjectiveRoutine();
33: Routines: TaoSetJacobianRoutine();
34: Routines: TaoSetInitialVector();
35: Routines: TaoSetFromOptions();
36: Routines: TaoSetConvergenceHistory(); TaoGetConvergenceHistory();
37: Routines: TaoSolve();
38: Routines: TaoView(); TaoDestroy();
39: Processors: 1
40: T*/
42: #define NOBSERVATIONS 214
43: #define NPARAMETERS 3
45: /* User-defined application context */
46: typedef struct {
47: /* Working space */
48: PetscReal t[NOBSERVATIONS]; /* array of independent variables of observation */
49: PetscReal y[NOBSERVATIONS]; /* array of dependent variables */
50: PetscReal j[NOBSERVATIONS][NPARAMETERS]; /* dense jacobian matrix array*/
51: PetscInt idm[NOBSERVATIONS]; /* Matrix indices for jacobian */
52: PetscInt idn[NPARAMETERS];
53: } AppCtx;
55: /* User provided Routines */
56: PetscErrorCode InitializeData(AppCtx *user);
57: PetscErrorCode FormStartingPoint(Vec);
58: PetscErrorCode EvaluateFunction(Tao, Vec, Vec, void *);
59: PetscErrorCode EvaluateJacobian(Tao, Vec, Mat, Mat, void *);
62: /*--------------------------------------------------------------------*/
65: int main(int argc,char **argv)
66: {
68: Vec x, f; /* solution, function */
69: Mat J; /* Jacobian matrix */
70: Tao tao; /* Tao solver context */
71: PetscInt i; /* iteration information */
72: PetscReal hist[100],resid[100];
73: PetscInt nhist,lits[100];
74: PetscBool printhistory;
75: AppCtx user; /* user-defined work context */
77: PetscInitialize(&argc,&argv,(char *)0,help);
79: printhistory = PETSC_FALSE;
80: PetscOptionsGetBool(NULL,"-printhistory",&printhistory,0);
81: /* Allocate vectors */
82: VecCreateSeq(MPI_COMM_SELF,NPARAMETERS,&x);
83: VecCreateSeq(MPI_COMM_SELF,NOBSERVATIONS,&f);
85: /* Create the Jacobian matrix. */
86: MatCreateSeqDense(MPI_COMM_SELF,NOBSERVATIONS,NPARAMETERS,NULL,&J);
88: for (i=0;i<NOBSERVATIONS;i++) user.idm[i] = i;
90: for (i=0;i<NPARAMETERS;i++) user.idn[i] = i;
92: /* Create TAO solver and set desired solution method */
93: TaoCreate(PETSC_COMM_SELF,&tao);
94: TaoSetType(tao,TAOPOUNDERS);
96: /* Set the function and Jacobian routines. */
97: InitializeData(&user);
98: FormStartingPoint(x);
99: TaoSetInitialVector(tao,x);
100: TaoSetSeparableObjectiveRoutine(tao,f,EvaluateFunction,(void*)&user);
101: TaoSetJacobianRoutine(tao, J, J, EvaluateJacobian, (void*)&user);
103: /* Check for any TAO command line arguments */
104: TaoSetFromOptions(tao);
106: TaoSetConvergenceHistory(tao,hist,resid,0,lits,100,PETSC_TRUE);
107: /* Perform the Solve */
108: TaoSolve(tao);
109: if (printhistory) {
110: TaoGetConvergenceHistory(tao,0,0,0,0,&nhist);
111: for (i=0;i<nhist;i++) {
112: PetscPrintf(PETSC_COMM_WORLD,"%g\t%g\n",(double)hist[i],(double)resid[i]);
113: }
114: }
115: TaoView(tao,PETSC_VIEWER_STDOUT_SELF);
117: /* Free TAO data structures */
118: TaoDestroy(&tao);
120: /* Free PETSc data structures */
121: VecDestroy(&x);
122: VecDestroy(&f);
123: MatDestroy(&J);
125: PetscFinalize();
126: return 0;
127: }
129: /*--------------------------------------------------------------------*/
132: PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr)
133: {
134: AppCtx *user = (AppCtx *)ptr;
135: PetscInt i;
136: PetscReal *y=user->y,*x,*f,*t=user->t;
140: VecGetArray(X,&x);
141: VecGetArray(F,&f);
143: for (i=0;i<NOBSERVATIONS;i++) {
144: f[i] = y[i] - PetscExpScalar(-x[0]*t[i])/(x[1] + x[2]*t[i]);
145: }
146: VecRestoreArray(X,&x);
147: VecRestoreArray(F,&f);
148: PetscLogFlops(6*NOBSERVATIONS);
149: return(0);
150: }
152: /*------------------------------------------------------------*/
153: /* J[i][j] = df[i]/dt[j] */
156: PetscErrorCode EvaluateJacobian(Tao tao, Vec X, Mat J, Mat Jpre, void *ptr)
157: {
158: AppCtx *user = (AppCtx *)ptr;
159: PetscInt i;
160: PetscReal *x,*t=user->t;
161: PetscReal base;
165: VecGetArray(X,&x);
166: for (i=0;i<NOBSERVATIONS;i++) {
167: base = PetscExpScalar(-x[0]*t[i])/(x[1] + x[2]*t[i]);
169: user->j[i][0] = t[i]*base;
170: user->j[i][1] = base/(x[1] + x[2]*t[i]);
171: user->j[i][2] = base*t[i]/(x[1] + x[2]*t[i]);
172: }
174: /* Assemble the matrix */
175: MatSetValues(J,NOBSERVATIONS,user->idm, NPARAMETERS, user->idn,(PetscReal *)user->j,INSERT_VALUES);
176: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
177: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
179: VecRestoreArray(X,&x);
180: PetscLogFlops(NOBSERVATIONS * 13);
181: return(0);
182: }
184: /* ------------------------------------------------------------ */
187: PetscErrorCode FormStartingPoint(Vec X)
188: {
189: PetscReal *x;
193: VecGetArray(X,&x);
194: x[0] = 0.15;
195: x[1] = 0.008;
196: x[2] = 0.010;
197: VecRestoreArray(X,&x);
198: return(0);
199: }
201: /* ---------------------------------------------------------------------- */
204: PetscErrorCode InitializeData(AppCtx *user)
205: {
206: PetscReal *t=user->t,*y=user->y;
207: PetscInt i=0;
210: y[i] = 92.9000; t[i++] = 0.5000;
211: y[i] = 78.7000; t[i++] = 0.6250;
212: y[i] = 64.2000; t[i++] = 0.7500;
213: y[i] = 64.9000; t[i++] = 0.8750;
214: y[i] = 57.1000; t[i++] = 1.0000;
215: y[i] = 43.3000; t[i++] = 1.2500;
216: y[i] = 31.1000; t[i++] = 1.7500;
217: y[i] = 23.6000; t[i++] = 2.2500;
218: y[i] = 31.0500; t[i++] = 1.7500;
219: y[i] = 23.7750; t[i++] = 2.2500;
220: y[i] = 17.7375; t[i++] = 2.7500;
221: y[i] = 13.8000; t[i++] = 3.2500;
222: y[i] = 11.5875; t[i++] = 3.7500;
223: y[i] = 9.4125; t[i++] = 4.2500;
224: y[i] = 7.7250; t[i++] = 4.7500;
225: y[i] = 7.3500; t[i++] = 5.2500;
226: y[i] = 8.0250; t[i++] = 5.7500;
227: y[i] = 90.6000; t[i++] = 0.5000;
228: y[i] = 76.9000; t[i++] = 0.6250;
229: y[i] = 71.6000; t[i++] = 0.7500;
230: y[i] = 63.6000; t[i++] = 0.8750;
231: y[i] = 54.0000; t[i++] = 1.0000;
232: y[i] = 39.2000; t[i++] = 1.2500;
233: y[i] = 29.3000; t[i++] = 1.7500;
234: y[i] = 21.4000; t[i++] = 2.2500;
235: y[i] = 29.1750; t[i++] = 1.7500;
236: y[i] = 22.1250; t[i++] = 2.2500;
237: y[i] = 17.5125; t[i++] = 2.7500;
238: y[i] = 14.2500; t[i++] = 3.2500;
239: y[i] = 9.4500; t[i++] = 3.7500;
240: y[i] = 9.1500; t[i++] = 4.2500;
241: y[i] = 7.9125; t[i++] = 4.7500;
242: y[i] = 8.4750; t[i++] = 5.2500;
243: y[i] = 6.1125; t[i++] = 5.7500;
244: y[i] = 80.0000; t[i++] = 0.5000;
245: y[i] = 79.0000; t[i++] = 0.6250;
246: y[i] = 63.8000; t[i++] = 0.7500;
247: y[i] = 57.2000; t[i++] = 0.8750;
248: y[i] = 53.2000; t[i++] = 1.0000;
249: y[i] = 42.5000; t[i++] = 1.2500;
250: y[i] = 26.8000; t[i++] = 1.7500;
251: y[i] = 20.4000; t[i++] = 2.2500;
252: y[i] = 26.8500; t[i++] = 1.7500;
253: y[i] = 21.0000; t[i++] = 2.2500;
254: y[i] = 16.4625; t[i++] = 2.7500;
255: y[i] = 12.5250; t[i++] = 3.2500;
256: y[i] = 10.5375; t[i++] = 3.7500;
257: y[i] = 8.5875; t[i++] = 4.2500;
258: y[i] = 7.1250; t[i++] = 4.7500;
259: y[i] = 6.1125; t[i++] = 5.2500;
260: y[i] = 5.9625; t[i++] = 5.7500;
261: y[i] = 74.1000; t[i++] = 0.5000;
262: y[i] = 67.3000; t[i++] = 0.6250;
263: y[i] = 60.8000; t[i++] = 0.7500;
264: y[i] = 55.5000; t[i++] = 0.8750;
265: y[i] = 50.3000; t[i++] = 1.0000;
266: y[i] = 41.0000; t[i++] = 1.2500;
267: y[i] = 29.4000; t[i++] = 1.7500;
268: y[i] = 20.4000; t[i++] = 2.2500;
269: y[i] = 29.3625; t[i++] = 1.7500;
270: y[i] = 21.1500; t[i++] = 2.2500;
271: y[i] = 16.7625; t[i++] = 2.7500;
272: y[i] = 13.2000; t[i++] = 3.2500;
273: y[i] = 10.8750; t[i++] = 3.7500;
274: y[i] = 8.1750; t[i++] = 4.2500;
275: y[i] = 7.3500; t[i++] = 4.7500;
276: y[i] = 5.9625; t[i++] = 5.2500;
277: y[i] = 5.6250; t[i++] = 5.7500;
278: y[i] = 81.5000; t[i++] = .5000;
279: y[i] = 62.4000; t[i++] = .7500;
280: y[i] = 32.5000; t[i++] = 1.5000;
281: y[i] = 12.4100; t[i++] = 3.0000;
282: y[i] = 13.1200; t[i++] = 3.0000;
283: y[i] = 15.5600; t[i++] = 3.0000;
284: y[i] = 5.6300; t[i++] = 6.0000;
285: y[i] = 78.0000; t[i++] = .5000;
286: y[i] = 59.9000; t[i++] = .7500;
287: y[i] = 33.2000; t[i++] = 1.5000;
288: y[i] = 13.8400; t[i++] = 3.0000;
289: y[i] = 12.7500; t[i++] = 3.0000;
290: y[i] = 14.6200; t[i++] = 3.0000;
291: y[i] = 3.9400; t[i++] = 6.0000;
292: y[i] = 76.8000; t[i++] = .5000;
293: y[i] = 61.0000; t[i++] = .7500;
294: y[i] = 32.9000; t[i++] = 1.5000;
295: y[i] = 13.8700; t[i++] = 3.0000;
296: y[i] = 11.8100; t[i++] = 3.0000;
297: y[i] = 13.3100; t[i++] = 3.0000;
298: y[i] = 5.4400; t[i++] = 6.0000;
299: y[i] = 78.0000; t[i++] = .5000;
300: y[i] = 63.5000; t[i++] = .7500;
301: y[i] = 33.8000; t[i++] = 1.5000;
302: y[i] = 12.5600; t[i++] = 3.0000;
303: y[i] = 5.6300; t[i++] = 6.0000;
304: y[i] = 12.7500; t[i++] = 3.0000;
305: y[i] = 13.1200; t[i++] = 3.0000;
306: y[i] = 5.4400; t[i++] = 6.0000;
307: y[i] = 76.8000; t[i++] = .5000;
308: y[i] = 60.0000; t[i++] = .7500;
309: y[i] = 47.8000; t[i++] = 1.0000;
310: y[i] = 32.0000; t[i++] = 1.5000;
311: y[i] = 22.2000; t[i++] = 2.0000;
312: y[i] = 22.5700; t[i++] = 2.0000;
313: y[i] = 18.8200; t[i++] = 2.5000;
314: y[i] = 13.9500; t[i++] = 3.0000;
315: y[i] = 11.2500; t[i++] = 4.0000;
316: y[i] = 9.0000; t[i++] = 5.0000;
317: y[i] = 6.6700; t[i++] = 6.0000;
318: y[i] = 75.8000; t[i++] = .5000;
319: y[i] = 62.0000; t[i++] = .7500;
320: y[i] = 48.8000; t[i++] = 1.0000;
321: y[i] = 35.2000; t[i++] = 1.5000;
322: y[i] = 20.0000; t[i++] = 2.0000;
323: y[i] = 20.3200; t[i++] = 2.0000;
324: y[i] = 19.3100; t[i++] = 2.5000;
325: y[i] = 12.7500; t[i++] = 3.0000;
326: y[i] = 10.4200; t[i++] = 4.0000;
327: y[i] = 7.3100; t[i++] = 5.0000;
328: y[i] = 7.4200; t[i++] = 6.0000;
329: y[i] = 70.5000; t[i++] = .5000;
330: y[i] = 59.5000; t[i++] = .7500;
331: y[i] = 48.5000; t[i++] = 1.0000;
332: y[i] = 35.8000; t[i++] = 1.5000;
333: y[i] = 21.0000; t[i++] = 2.0000;
334: y[i] = 21.6700; t[i++] = 2.0000;
335: y[i] = 21.0000; t[i++] = 2.5000;
336: y[i] = 15.6400; t[i++] = 3.0000;
337: y[i] = 8.1700; t[i++] = 4.0000;
338: y[i] = 8.5500; t[i++] = 5.0000;
339: y[i] = 10.1200; t[i++] = 6.0000;
340: y[i] = 78.0000; t[i++] = .5000;
341: y[i] = 66.0000; t[i++] = .6250;
342: y[i] = 62.0000; t[i++] = .7500;
343: y[i] = 58.0000; t[i++] = .8750;
344: y[i] = 47.7000; t[i++] = 1.0000;
345: y[i] = 37.8000; t[i++] = 1.2500;
346: y[i] = 20.2000; t[i++] = 2.2500;
347: y[i] = 21.0700; t[i++] = 2.2500;
348: y[i] = 13.8700; t[i++] = 2.7500;
349: y[i] = 9.6700; t[i++] = 3.2500;
350: y[i] = 7.7600; t[i++] = 3.7500;
351: y[i] = 5.4400; t[i++] = 4.2500;
352: y[i] = 4.8700; t[i++] = 4.7500;
353: y[i] = 4.0100; t[i++] = 5.2500;
354: y[i] = 3.7500; t[i++] = 5.7500;
355: y[i] = 24.1900; t[i++] = 3.0000;
356: y[i] = 25.7600; t[i++] = 3.0000;
357: y[i] = 18.0700; t[i++] = 3.0000;
358: y[i] = 11.8100; t[i++] = 3.0000;
359: y[i] = 12.0700; t[i++] = 3.0000;
360: y[i] = 16.1200; t[i++] = 3.0000;
361: y[i] = 70.8000; t[i++] = .5000;
362: y[i] = 54.7000; t[i++] = .7500;
363: y[i] = 48.0000; t[i++] = 1.0000;
364: y[i] = 39.8000; t[i++] = 1.5000;
365: y[i] = 29.8000; t[i++] = 2.0000;
366: y[i] = 23.7000; t[i++] = 2.5000;
367: y[i] = 29.6200; t[i++] = 2.0000;
368: y[i] = 23.8100; t[i++] = 2.5000;
369: y[i] = 17.7000; t[i++] = 3.0000;
370: y[i] = 11.5500; t[i++] = 4.0000;
371: y[i] = 12.0700; t[i++] = 5.0000;
372: y[i] = 8.7400; t[i++] = 6.0000;
373: y[i] = 80.7000; t[i++] = .5000;
374: y[i] = 61.3000; t[i++] = .7500;
375: y[i] = 47.5000; t[i++] = 1.0000;
376: y[i] = 29.0000; t[i++] = 1.5000;
377: y[i] = 24.0000; t[i++] = 2.0000;
378: y[i] = 17.7000; t[i++] = 2.5000;
379: y[i] = 24.5600; t[i++] = 2.0000;
380: y[i] = 18.6700; t[i++] = 2.5000;
381: y[i] = 16.2400; t[i++] = 3.0000;
382: y[i] = 8.7400; t[i++] = 4.0000;
383: y[i] = 7.8700; t[i++] = 5.0000;
384: y[i] = 8.5100; t[i++] = 6.0000;
385: y[i] = 66.7000; t[i++] = .5000;
386: y[i] = 59.2000; t[i++] = .7500;
387: y[i] = 40.8000; t[i++] = 1.0000;
388: y[i] = 30.7000; t[i++] = 1.5000;
389: y[i] = 25.7000; t[i++] = 2.0000;
390: y[i] = 16.3000; t[i++] = 2.5000;
391: y[i] = 25.9900; t[i++] = 2.0000;
392: y[i] = 16.9500; t[i++] = 2.5000;
393: y[i] = 13.3500; t[i++] = 3.0000;
394: y[i] = 8.6200; t[i++] = 4.0000;
395: y[i] = 7.2000; t[i++] = 5.0000;
396: y[i] = 6.6400; t[i++] = 6.0000;
397: y[i] = 13.6900; t[i++] = 3.0000;
398: y[i] = 81.0000; t[i++] = .5000;
399: y[i] = 64.5000; t[i++] = .7500;
400: y[i] = 35.5000; t[i++] = 1.5000;
401: y[i] = 13.3100; t[i++] = 3.0000;
402: y[i] = 4.8700; t[i++] = 6.0000;
403: y[i] = 12.9400; t[i++] = 3.0000;
404: y[i] = 5.0600; t[i++] = 6.0000;
405: y[i] = 15.1900; t[i++] = 3.0000;
406: y[i] = 14.6200; t[i++] = 3.0000;
407: y[i] = 15.6400; t[i++] = 3.0000;
408: y[i] = 25.5000; t[i++] = 1.7500;
409: y[i] = 25.9500; t[i++] = 1.7500;
410: y[i] = 81.7000; t[i++] = .5000;
411: y[i] = 61.6000; t[i++] = .7500;
412: y[i] = 29.8000; t[i++] = 1.7500;
413: y[i] = 29.8100; t[i++] = 1.7500;
414: y[i] = 17.1700; t[i++] = 2.7500;
415: y[i] = 10.3900; t[i++] = 3.7500;
416: y[i] = 28.4000; t[i++] = 1.7500;
417: y[i] = 28.6900; t[i++] = 1.7500;
418: y[i] = 81.3000; t[i++] = .5000;
419: y[i] = 60.9000; t[i++] = .7500;
420: y[i] = 16.6500; t[i++] = 2.7500;
421: y[i] = 10.0500; t[i++] = 3.7500;
422: y[i] = 28.9000; t[i++] = 1.7500;
423: y[i] = 28.9500; t[i++] = 1.7500;
424: return(0);
425: }