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

 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: TaoSetInitialVector();
 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: {
 63:   Vec            x, f;               /* solution, function */
 64:   Mat            J;                  /* Jacobian matrix */
 65:   Tao            tao;                /* Tao solver context */
 66:   PetscInt       i;               /* iteration information */
 67:   PetscReal      hist[100],resid[100];
 68:   PetscInt       lits[100];
 69:   AppCtx         user;               /* user-defined work context */

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

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

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

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

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

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

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

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

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

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

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

112:   PetscFinalize();
113:   return ierr;
114: }

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

126:   VecGetArrayRead(X,&x);
127:   VecGetArray(F,&f);

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

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

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

154:     user->j[i][0] = t[i]*base;
155:     user->j[i][1] = base/(x[1] + x[2]*t[i]);
156:     user->j[i][2] = base*t[i]/(x[1] + x[2]*t[i]);
157:   }

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

164:   VecRestoreArrayRead(X,&x);
165:   PetscLogFlops(NOBSERVATIONS * 13);
166:   return(0);
167: }

169: /* ------------------------------------------------------------ */
170: PetscErrorCode FormStartingPoint(Vec X)
171: {
172:   PetscReal      *x;

176:   VecGetArray(X,&x);
177:   x[0] = 0.15;
178:   x[1] = 0.008;
179:   x[2] = 0.010;
180:   VecRestoreArray(X,&x);
181:   return(0);
182: }

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

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

408: /*TEST

410:    build:
411:       requires: !complex !single

413:    test:
414:       args: -tao_smonitor -tao_max_it 100 -tao_type pounders -tao_gatol 1.e-5

416:    test:
417:       suffix: 2
418:       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

420:    test:
421:       suffix: 3
422:       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

424:    test:
425:       suffix: 4
426:       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

428: TEST*/