Actual source code: chwirut1.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: */

 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: */

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

 26: /*T
 27:    Concepts: TAO^Solving a system of nonlinear equations, nonlinear least squares
 28:    Routines: TaoCreate();
 29:    Routines: TaoSetType();
 30:    Routines: TaoSetSeparableObjectiveRoutine();
 31:    Routines: TaoSetJacobianRoutine();
 32:    Routines: TaoSetSolution();
 33:    Routines: TaoSetFromOptions();
 34:    Routines: TaoSetConvergenceHistory(); TaoGetConvergenceHistory();
 35:    Routines: TaoSolve();
 36:    Routines: TaoView(); TaoDestroy();
 37:    Processors: 1
 38: T*/

 40: #define NOBSERVATIONS 214
 41: #define NPARAMETERS 3

 43: /* User-defined application context */
 44: typedef struct {
 45:   /* Working space */
 46:   PetscReal t[NOBSERVATIONS];   /* array of independent variables of observation */
 47:   PetscReal y[NOBSERVATIONS];   /* array of dependent variables */
 48:   PetscReal j[NOBSERVATIONS][NPARAMETERS]; /* dense jacobian matrix array*/
 49:   PetscInt idm[NOBSERVATIONS];  /* Matrix indices for jacobian */
 50:   PetscInt idn[NPARAMETERS];
 51: } AppCtx;

 53: /* User provided Routines */
 54: PetscErrorCode InitializeData(AppCtx *user);
 55: PetscErrorCode FormStartingPoint(Vec);
 56: PetscErrorCode EvaluateFunction(Tao, Vec, Vec, void *);
 57: PetscErrorCode EvaluateJacobian(Tao, Vec, Mat, Mat, void *);

 59: /*--------------------------------------------------------------------*/
 60: int main(int argc,char **argv)
 61: {
 62:   Vec            x, f;               /* solution, function */
 63:   Mat            J;                  /* Jacobian matrix */
 64:   Tao            tao;                /* Tao solver context */
 65:   PetscInt       i;               /* iteration information */
 66:   PetscReal      hist[100],resid[100];
 67:   PetscInt       lits[100];
 68:   AppCtx         user;               /* user-defined work context */

 70:   PetscInitialize(&argc,&argv,(char *)0,help);
 71:   /* Allocate vectors */
 72:   VecCreateSeq(MPI_COMM_SELF,NPARAMETERS,&x);
 73:   VecCreateSeq(MPI_COMM_SELF,NOBSERVATIONS,&f);

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

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

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

 82:   /* Create TAO solver and set desired solution method */
 83:   TaoCreate(PETSC_COMM_SELF,&tao);
 84:   TaoSetType(tao,TAOPOUNDERS);

 86:  /* Set the function and Jacobian routines. */
 87:   InitializeData(&user);
 88:   FormStartingPoint(x);
 89:   TaoSetSolution(tao,x);
 90:   TaoSetResidualRoutine(tao,f,EvaluateFunction,(void*)&user);
 91:   TaoSetJacobianResidualRoutine(tao, J, J, EvaluateJacobian, (void*)&user);

 93:   /* Check for any TAO command line arguments */
 94:   TaoSetFromOptions(tao);

 96:   TaoSetConvergenceHistory(tao,hist,resid,0,lits,100,PETSC_TRUE);
 97:   /* Perform the Solve */
 98:   TaoSolve(tao);

100:   /* View the vector; then destroy it.  */
101:   VecView(x,PETSC_VIEWER_STDOUT_SELF);

103:   /* Free TAO data structures */
104:   TaoDestroy(&tao);

106:    /* Free PETSc data structures */
107:   VecDestroy(&x);
108:   VecDestroy(&f);
109:   MatDestroy(&J);

111:   PetscFinalize();
112:   return 0;
113: }

115: /*--------------------------------------------------------------------*/
116: PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr)
117: {
118:   AppCtx         *user = (AppCtx *)ptr;
119:   PetscInt       i;
120:   const PetscReal *x;
121:   PetscReal      *y=user->y,*f,*t=user->t;

123:   VecGetArrayRead(X,&x);
124:   VecGetArray(F,&f);

126:   for (i=0;i<NOBSERVATIONS;i++) {
127:     f[i] = y[i] - PetscExpScalar(-x[0]*t[i])/(x[1] + x[2]*t[i]);
128:   }
129:   VecRestoreArrayRead(X,&x);
130:   VecRestoreArray(F,&f);
131:   PetscLogFlops(6*NOBSERVATIONS);
132:   return 0;
133: }

135: /*------------------------------------------------------------*/
136: /* J[i][j] = df[i]/dt[j] */
137: PetscErrorCode EvaluateJacobian(Tao tao, Vec X, Mat J, Mat Jpre, void *ptr)
138: {
139:   AppCtx         *user = (AppCtx *)ptr;
140:   PetscInt       i;
141:   const PetscReal *x;
142:   PetscReal      *t=user->t;
143:   PetscReal      base;

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

149:     user->j[i][0] = t[i]*base;
150:     user->j[i][1] = base/(x[1] + x[2]*t[i]);
151:     user->j[i][2] = base*t[i]/(x[1] + x[2]*t[i]);
152:   }

154:   /* Assemble the matrix */
155:   MatSetValues(J,NOBSERVATIONS,user->idm, NPARAMETERS, user->idn,(PetscReal *)user->j,INSERT_VALUES);
156:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
157:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);

159:   VecRestoreArrayRead(X,&x);
160:   PetscLogFlops(NOBSERVATIONS * 13);
161:   return 0;
162: }

164: /* ------------------------------------------------------------ */
165: PetscErrorCode FormStartingPoint(Vec X)
166: {
167:   PetscReal      *x;

169:   VecGetArray(X,&x);
170:   x[0] = 0.15;
171:   x[1] = 0.008;
172:   x[2] = 0.010;
173:   VecRestoreArray(X,&x);
174:   return 0;
175: }

177: /* ---------------------------------------------------------------------- */
178: PetscErrorCode InitializeData(AppCtx *user)
179: {
180:   PetscReal *t=user->t,*y=user->y;
181:   PetscInt  i=0;

183:   y[i] =   92.9000;   t[i++] =  0.5000;
184:   y[i] =    78.7000;  t[i++] =   0.6250;
185:   y[i] =    64.2000;  t[i++] =   0.7500;
186:   y[i] =    64.9000;  t[i++] =   0.8750;
187:   y[i] =    57.1000;  t[i++] =   1.0000;
188:   y[i] =    43.3000;  t[i++] =   1.2500;
189:   y[i] =    31.1000;   t[i++] =  1.7500;
190:   y[i] =    23.6000;   t[i++] =  2.2500;
191:   y[i] =    31.0500;   t[i++] =  1.7500;
192:   y[i] =    23.7750;   t[i++] =  2.2500;
193:   y[i] =    17.7375;   t[i++] =  2.7500;
194:   y[i] =    13.8000;   t[i++] =  3.2500;
195:   y[i] =    11.5875;   t[i++] =  3.7500;
196:   y[i] =     9.4125;   t[i++] =  4.2500;
197:   y[i] =     7.7250;   t[i++] =  4.7500;
198:   y[i] =     7.3500;   t[i++] =  5.2500;
199:   y[i] =     8.0250;   t[i++] =  5.7500;
200:   y[i] =    90.6000;   t[i++] =  0.5000;
201:   y[i] =    76.9000;   t[i++] =  0.6250;
202:   y[i] =    71.6000;   t[i++] = 0.7500;
203:   y[i] =    63.6000;   t[i++] =  0.8750;
204:   y[i] =    54.0000;   t[i++] =  1.0000;
205:   y[i] =    39.2000;   t[i++] =  1.2500;
206:   y[i] =    29.3000;   t[i++] = 1.7500;
207:   y[i] =    21.4000;   t[i++] =  2.2500;
208:   y[i] =    29.1750;   t[i++] =  1.7500;
209:   y[i] =    22.1250;   t[i++] =  2.2500;
210:   y[i] =    17.5125;   t[i++] =  2.7500;
211:   y[i] =    14.2500;   t[i++] =  3.2500;
212:   y[i] =     9.4500;   t[i++] =  3.7500;
213:   y[i] =     9.1500;   t[i++] =  4.2500;
214:   y[i] =     7.9125;   t[i++] =  4.7500;
215:   y[i] =     8.4750;   t[i++] =  5.2500;
216:   y[i] =     6.1125;   t[i++] =  5.7500;
217:   y[i] =    80.0000;   t[i++] =  0.5000;
218:   y[i] =    79.0000;   t[i++] =  0.6250;
219:   y[i] =    63.8000;   t[i++] =  0.7500;
220:   y[i] =    57.2000;   t[i++] =  0.8750;
221:   y[i] =    53.2000;   t[i++] =  1.0000;
222:   y[i] =   42.5000;   t[i++] =  1.2500;
223:   y[i] =   26.8000;   t[i++] =  1.7500;
224:   y[i] =    20.4000;   t[i++] =  2.2500;
225:   y[i] =    26.8500;  t[i++] =   1.7500;
226:   y[i] =    21.0000;  t[i++] =   2.2500;
227:   y[i] =    16.4625;  t[i++] =   2.7500;
228:   y[i] =    12.5250;  t[i++] =   3.2500;
229:   y[i] =    10.5375;  t[i++] =   3.7500;
230:   y[i] =     8.5875;  t[i++] =   4.2500;
231:   y[i] =     7.1250;  t[i++] =   4.7500;
232:   y[i] =     6.1125;  t[i++] =   5.2500;
233:   y[i] =     5.9625;  t[i++] =   5.7500;
234:   y[i] =    74.1000;  t[i++] =   0.5000;
235:   y[i] =    67.3000;  t[i++] =   0.6250;
236:   y[i] =    60.8000;  t[i++] =   0.7500;
237:   y[i] =    55.5000;  t[i++] =   0.8750;
238:   y[i] =    50.3000;  t[i++] =   1.0000;
239:   y[i] =    41.0000;  t[i++] =   1.2500;
240:   y[i] =    29.4000;  t[i++] =   1.7500;
241:   y[i] =    20.4000;  t[i++] =   2.2500;
242:   y[i] =    29.3625;  t[i++] =   1.7500;
243:   y[i] =    21.1500;  t[i++] =   2.2500;
244:   y[i] =    16.7625;  t[i++] =   2.7500;
245:   y[i] =    13.2000;  t[i++] =   3.2500;
246:   y[i] =    10.8750;  t[i++] =   3.7500;
247:   y[i] =     8.1750;  t[i++] =   4.2500;
248:   y[i] =     7.3500;  t[i++] =   4.7500;
249:   y[i] =     5.9625;  t[i++] =  5.2500;
250:   y[i] =     5.6250;  t[i++] =   5.7500;
251:   y[i] =    81.5000;  t[i++] =    .5000;
252:   y[i] =    62.4000;  t[i++] =    .7500;
253:   y[i] =    32.5000;  t[i++] =   1.5000;
254:   y[i] =    12.4100;  t[i++] =   3.0000;
255:   y[i] =    13.1200;  t[i++] =   3.0000;
256:   y[i] =    15.5600;  t[i++] =   3.0000;
257:   y[i] =     5.6300;  t[i++] =   6.0000;
258:   y[i] =    78.0000;   t[i++] =   .5000;
259:   y[i] =    59.9000;  t[i++] =    .7500;
260:   y[i] =    33.2000;  t[i++] =   1.5000;
261:   y[i] =    13.8400;  t[i++] =   3.0000;
262:   y[i] =    12.7500;  t[i++] =   3.0000;
263:   y[i] =    14.6200;  t[i++] =   3.0000;
264:   y[i] =     3.9400;  t[i++] =   6.0000;
265:   y[i] =    76.8000;  t[i++] =    .5000;
266:   y[i] =    61.0000;  t[i++] =    .7500;
267:   y[i] =    32.9000;  t[i++] =   1.5000;
268:   y[i] =   13.8700;   t[i++] = 3.0000;
269:   y[i] =    11.8100;  t[i++] =   3.0000;
270:   y[i] =    13.3100;  t[i++] =   3.0000;
271:   y[i] =     5.4400;  t[i++] =   6.0000;
272:   y[i] =    78.0000;  t[i++] =    .5000;
273:   y[i] =    63.5000;  t[i++] =    .7500;
274:   y[i] =    33.8000;  t[i++] =   1.5000;
275:   y[i] =    12.5600;  t[i++] =   3.0000;
276:   y[i] =     5.6300;  t[i++] =   6.0000;
277:   y[i] =    12.7500;  t[i++] =   3.0000;
278:   y[i] =    13.1200;  t[i++] =   3.0000;
279:   y[i] =     5.4400;  t[i++] =   6.0000;
280:   y[i] =    76.8000;  t[i++] =    .5000;
281:   y[i] =    60.0000;  t[i++] =    .7500;
282:   y[i] =    47.8000;  t[i++] =   1.0000;
283:   y[i] =    32.0000;  t[i++] =   1.5000;
284:   y[i] =    22.2000;  t[i++] =   2.0000;
285:   y[i] =    22.5700;  t[i++] =   2.0000;
286:   y[i] =    18.8200;  t[i++] =   2.5000;
287:   y[i] =    13.9500;  t[i++] =   3.0000;
288:   y[i] =    11.2500;  t[i++] =   4.0000;
289:   y[i] =     9.0000;  t[i++] =   5.0000;
290:   y[i] =     6.6700;  t[i++] =   6.0000;
291:   y[i] =    75.8000;  t[i++] =    .5000;
292:   y[i] =    62.0000;  t[i++] =    .7500;
293:   y[i] =    48.8000;  t[i++] =   1.0000;
294:   y[i] =    35.2000;  t[i++] =   1.5000;
295:   y[i] =    20.0000;  t[i++] =   2.0000;
296:   y[i] =    20.3200;  t[i++] =   2.0000;
297:   y[i] =    19.3100;  t[i++] =   2.5000;
298:   y[i] =    12.7500;  t[i++] =   3.0000;
299:   y[i] =    10.4200;  t[i++] =   4.0000;
300:   y[i] =     7.3100;  t[i++] =   5.0000;
301:   y[i] =     7.4200;  t[i++] =   6.0000;
302:   y[i] =    70.5000;  t[i++] =    .5000;
303:   y[i] =    59.5000;  t[i++] =    .7500;
304:   y[i] =    48.5000;  t[i++] =   1.0000;
305:   y[i] =    35.8000;  t[i++] =   1.5000;
306:   y[i] =    21.0000;  t[i++] =   2.0000;
307:   y[i] =    21.6700;  t[i++] =   2.0000;
308:   y[i] =    21.0000;  t[i++] =   2.5000;
309:   y[i] =    15.6400;  t[i++] =   3.0000;
310:   y[i] =     8.1700;  t[i++] =   4.0000;
311:   y[i] =     8.5500;  t[i++] =   5.0000;
312:   y[i] =    10.1200;  t[i++] =   6.0000;
313:   y[i] =    78.0000;  t[i++] =    .5000;
314:   y[i] =    66.0000;  t[i++] =    .6250;
315:   y[i] =    62.0000;  t[i++] =    .7500;
316:   y[i] =    58.0000;  t[i++] =    .8750;
317:   y[i] =    47.7000;  t[i++] =   1.0000;
318:   y[i] =    37.8000;  t[i++] =   1.2500;
319:   y[i] =    20.2000;  t[i++] =   2.2500;
320:   y[i] =    21.0700;  t[i++] =   2.2500;
321:   y[i] =    13.8700;  t[i++] =   2.7500;
322:   y[i] =     9.6700;  t[i++] =   3.2500;
323:   y[i] =     7.7600;  t[i++] =   3.7500;
324:   y[i] =    5.4400;   t[i++] =  4.2500;
325:   y[i] =    4.8700;   t[i++] =  4.7500;
326:   y[i] =     4.0100;  t[i++] =   5.2500;
327:   y[i] =     3.7500;  t[i++] =   5.7500;
328:   y[i] =    24.1900;  t[i++] =   3.0000;
329:   y[i] =    25.7600;  t[i++] =   3.0000;
330:   y[i] =    18.0700;  t[i++] =   3.0000;
331:   y[i] =    11.8100;  t[i++] =   3.0000;
332:   y[i] =    12.0700;  t[i++] =   3.0000;
333:   y[i] =    16.1200;  t[i++] =   3.0000;
334:   y[i] =    70.8000;  t[i++] =    .5000;
335:   y[i] =    54.7000;  t[i++] =    .7500;
336:   y[i] =    48.0000;  t[i++] =   1.0000;
337:   y[i] =    39.8000;  t[i++] =   1.5000;
338:   y[i] =    29.8000;  t[i++] =   2.0000;
339:   y[i] =    23.7000;  t[i++] =   2.5000;
340:   y[i] =    29.6200;  t[i++] =   2.0000;
341:   y[i] =    23.8100;  t[i++] =   2.5000;
342:   y[i] =    17.7000;  t[i++] =   3.0000;
343:   y[i] =    11.5500;  t[i++] =   4.0000;
344:   y[i] =    12.0700;  t[i++] =   5.0000;
345:   y[i] =     8.7400;  t[i++] =   6.0000;
346:   y[i] =    80.7000;  t[i++] =    .5000;
347:   y[i] =    61.3000;  t[i++] =    .7500;
348:   y[i] =    47.5000;  t[i++] =   1.0000;
349:    y[i] =   29.0000;  t[i++] =   1.5000;
350:    y[i] =   24.0000;  t[i++] =   2.0000;
351:   y[i] =    17.7000;  t[i++] =   2.5000;
352:   y[i] =    24.5600;  t[i++] =   2.0000;
353:   y[i] =    18.6700;  t[i++] =   2.5000;
354:    y[i] =   16.2400;  t[i++] =   3.0000;
355:   y[i] =     8.7400;  t[i++] =   4.0000;
356:   y[i] =     7.8700;  t[i++] =   5.0000;
357:   y[i] =     8.5100;  t[i++] =   6.0000;
358:   y[i] =    66.7000;  t[i++] =    .5000;
359:   y[i] =    59.2000;  t[i++] =    .7500;
360:   y[i] =    40.8000;  t[i++] =   1.0000;
361:   y[i] =    30.7000;  t[i++] =   1.5000;
362:   y[i] =    25.7000;  t[i++] =   2.0000;
363:   y[i] =    16.3000;  t[i++] =   2.5000;
364:   y[i] =    25.9900;  t[i++] =   2.0000;
365:   y[i] =    16.9500;  t[i++] =   2.5000;
366:   y[i] =    13.3500;  t[i++] =   3.0000;
367:   y[i] =     8.6200;  t[i++] =   4.0000;
368:   y[i] =     7.2000;  t[i++] =   5.0000;
369:   y[i] =     6.6400;  t[i++] =   6.0000;
370:   y[i] =    13.6900;  t[i++] =   3.0000;
371:   y[i] =    81.0000;  t[i++] =    .5000;
372:   y[i] =    64.5000;  t[i++] =    .7500;
373:   y[i] =    35.5000;  t[i++] =   1.5000;
374:    y[i] =   13.3100;  t[i++] =   3.0000;
375:   y[i] =     4.8700;  t[i++] =   6.0000;
376:   y[i] =    12.9400;  t[i++] =   3.0000;
377:   y[i] =     5.0600;  t[i++] =   6.0000;
378:   y[i] =    15.1900;  t[i++] =   3.0000;
379:   y[i] =    14.6200;  t[i++] =   3.0000;
380:   y[i] =    15.6400;  t[i++] =   3.0000;
381:   y[i] =    25.5000;  t[i++] =   1.7500;
382:   y[i] =    25.9500;  t[i++] =   1.7500;
383:   y[i] =    81.7000;  t[i++] =    .5000;
384:   y[i] =    61.6000;  t[i++] =    .7500;
385:   y[i] =    29.8000;  t[i++] =   1.7500;
386:   y[i] =    29.8100;  t[i++] =   1.7500;
387:   y[i] =    17.1700;  t[i++] =   2.7500;
388:   y[i] =    10.3900;  t[i++] =   3.7500;
389:   y[i] =    28.4000;  t[i++] =   1.7500;
390:   y[i] =    28.6900;  t[i++] =   1.7500;
391:   y[i] =    81.3000;  t[i++] =    .5000;
392:   y[i] =    60.9000;  t[i++] =    .7500;
393:   y[i] =    16.6500;  t[i++] =   2.7500;
394:   y[i] =    10.0500;  t[i++] =   3.7500;
395:   y[i] =    28.9000;  t[i++] =   1.7500;
396:   y[i] =    28.9500;  t[i++] =   1.7500;
397:   return 0;
398: }

400: /*TEST

402:    build:
403:       requires: !complex !single

405:    test:
406:       args: -tao_smonitor -tao_max_it 100 -tao_type pounders -tao_gatol 1.e-5

408:    test:
409:       suffix: 2
410:       args: -tao_smonitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type l2prox -tao_brgn_regularizer_weight 1e-4 -tao_gatol 1.e-5

412:    test:
413:       suffix: 3
414:       args: -tao_smonitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type l1dict -tao_brgn_regularizer_weight 1e-4 -tao_brgn_l1_smooth_epsilon 1e-6 -tao_gatol 1.e-5

416:    test:
417:       suffix: 4
418:       args: -tao_smonitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type lm -tao_gatol 1.e-5 -tao_brgn_subsolver_tao_type bnls

420: TEST*/