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