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