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