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 /* T 27 Concepts: TAO^Solving a system of nonlinear equations, nonlinear least squares 28 Routines: TaoCreate(); 29 Routines: TaoSetType(); 30 Routines: TaoSetResidualRoutine(); 31 Routines: TaoSetMonitor(); 32 Routines: TaoSetSolution(); 33 Routines: TaoSetFromOptions(); 34 Routines: TaoSolve(); 35 Routines: TaoDestroy(); 36 Processors: n 37 T*/ 38 39 #define NOBSERVATIONS 214 40 #define NPARAMETERS 3 41 42 #define DIE_TAG 2000 43 #define IDLE_TAG 1000 44 45 /* User-defined application context */ 46 typedef struct { 47 /* Working space */ 48 PetscReal t[NOBSERVATIONS]; /* array of independent variables of observation */ 49 PetscReal y[NOBSERVATIONS]; /* array of dependent variables */ 50 PetscMPIInt size,rank; 51 } AppCtx; 52 53 /* User provided Routines */ 54 PetscErrorCode InitializeData(AppCtx *user); 55 PetscErrorCode FormStartingPoint(Vec); 56 PetscErrorCode EvaluateFunction(Tao, Vec, Vec, void *); 57 PetscErrorCode TaskWorker(AppCtx *user); 58 PetscErrorCode StopWorkers(AppCtx *user); 59 PetscErrorCode RunSimulation(PetscReal *x, PetscInt i, PetscReal*f, AppCtx *user); 60 61 /*--------------------------------------------------------------------*/ 62 int main(int argc,char **argv) 63 { 64 Vec x, f; /* solution, function */ 65 Tao tao; /* Tao solver context */ 66 AppCtx user; /* user-defined work context */ 67 68 /* Initialize TAO and PETSc */ 69 PetscCall(PetscInitialize(&argc,&argv,(char *)0,help)); 70 MPI_Comm_size(MPI_COMM_WORLD,&user.size); 71 MPI_Comm_rank(MPI_COMM_WORLD,&user.rank); 72 PetscCall(InitializeData(&user)); 73 74 /* Run optimization on rank 0 */ 75 if (user.rank == 0) { 76 /* Allocate vectors */ 77 PetscCall(VecCreateSeq(PETSC_COMM_SELF,NPARAMETERS,&x)); 78 PetscCall(VecCreateSeq(PETSC_COMM_SELF,NOBSERVATIONS,&f)); 79 80 /* TAO code begins here */ 81 82 /* Create TAO solver and set desired solution method */ 83 PetscCall(TaoCreate(PETSC_COMM_SELF,&tao)); 84 PetscCall(TaoSetType(tao,TAOPOUNDERS)); 85 86 /* Set the function and Jacobian routines. */ 87 PetscCall(FormStartingPoint(x)); 88 PetscCall(TaoSetSolution(tao,x)); 89 PetscCall(TaoSetResidualRoutine(tao,f,EvaluateFunction,(void*)&user)); 90 91 /* Check for any TAO command line arguments */ 92 PetscCall(TaoSetFromOptions(tao)); 93 94 /* Perform the Solve */ 95 PetscCall(TaoSolve(tao)); 96 97 /* Free TAO data structures */ 98 PetscCall(TaoDestroy(&tao)); 99 100 /* Free PETSc data structures */ 101 PetscCall(VecDestroy(&x)); 102 PetscCall(VecDestroy(&f)); 103 StopWorkers(&user); 104 } else { 105 TaskWorker(&user); 106 } 107 PetscCall(PetscFinalize()); 108 return 0; 109 } 110 111 /*--------------------------------------------------------------------*/ 112 PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr) 113 { 114 AppCtx *user = (AppCtx *)ptr; 115 PetscInt i; 116 PetscReal *x,*f; 117 118 PetscFunctionBegin; 119 PetscCall(VecGetArray(X,&x)); 120 PetscCall(VecGetArray(F,&f)); 121 if (user->size == 1) { 122 /* Single processor */ 123 for (i=0;i<NOBSERVATIONS;i++) { 124 PetscCall(RunSimulation(x,i,&f[i],user)); 125 } 126 } else { 127 /* Multiprocessor main */ 128 PetscMPIInt tag; 129 PetscInt finishedtasks,next_task,checkedin; 130 PetscReal f_i=0.0; 131 MPI_Status status; 132 133 next_task=0; 134 finishedtasks=0; 135 checkedin=0; 136 137 while (finishedtasks < NOBSERVATIONS || checkedin < user->size-1) { 138 PetscCallMPI(MPI_Recv(&f_i,1,MPIU_REAL,MPI_ANY_SOURCE,MPI_ANY_TAG,PETSC_COMM_WORLD,&status)); 139 if (status.MPI_TAG == IDLE_TAG) { 140 checkedin++; 141 } else { 142 143 tag = status.MPI_TAG; 144 f[tag] = (PetscReal)f_i; 145 finishedtasks++; 146 } 147 148 if (next_task<NOBSERVATIONS) { 149 PetscCallMPI(MPI_Send(x,NPARAMETERS,MPIU_REAL,status.MPI_SOURCE,next_task,PETSC_COMM_WORLD)); 150 next_task++; 151 152 } else { 153 /* Send idle message */ 154 PetscCallMPI(MPI_Send(x,NPARAMETERS,MPIU_REAL,status.MPI_SOURCE,IDLE_TAG,PETSC_COMM_WORLD)); 155 } 156 } 157 } 158 PetscCall(VecRestoreArray(X,&x)); 159 PetscCall(VecRestoreArray(F,&f)); 160 PetscLogFlops(6*NOBSERVATIONS); 161 PetscFunctionReturn(0); 162 } 163 164 /* ------------------------------------------------------------ */ 165 PetscErrorCode FormStartingPoint(Vec X) 166 { 167 PetscReal *x; 168 169 PetscFunctionBegin; 170 PetscCall(VecGetArray(X,&x)); 171 x[0] = 0.15; 172 x[1] = 0.008; 173 x[2] = 0.010; 174 PetscCall(VecRestoreArray(X,&x)); 175 PetscFunctionReturn(0); 176 } 177 178 /* ---------------------------------------------------------------------- */ 179 PetscErrorCode InitializeData(AppCtx *user) 180 { 181 PetscReal *t=user->t,*y=user->y; 182 PetscInt i=0; 183 184 PetscFunctionBegin; 185 y[i] = 92.9000; t[i++] = 0.5000; 186 y[i] = 78.7000; t[i++] = 0.6250; 187 y[i] = 64.2000; t[i++] = 0.7500; 188 y[i] = 64.9000; t[i++] = 0.8750; 189 y[i] = 57.1000; t[i++] = 1.0000; 190 y[i] = 43.3000; t[i++] = 1.2500; 191 y[i] = 31.1000; t[i++] = 1.7500; 192 y[i] = 23.6000; t[i++] = 2.2500; 193 y[i] = 31.0500; t[i++] = 1.7500; 194 y[i] = 23.7750; t[i++] = 2.2500; 195 y[i] = 17.7375; t[i++] = 2.7500; 196 y[i] = 13.8000; t[i++] = 3.2500; 197 y[i] = 11.5875; t[i++] = 3.7500; 198 y[i] = 9.4125; t[i++] = 4.2500; 199 y[i] = 7.7250; t[i++] = 4.7500; 200 y[i] = 7.3500; t[i++] = 5.2500; 201 y[i] = 8.0250; t[i++] = 5.7500; 202 y[i] = 90.6000; t[i++] = 0.5000; 203 y[i] = 76.9000; t[i++] = 0.6250; 204 y[i] = 71.6000; t[i++] = 0.7500; 205 y[i] = 63.6000; t[i++] = 0.8750; 206 y[i] = 54.0000; t[i++] = 1.0000; 207 y[i] = 39.2000; t[i++] = 1.2500; 208 y[i] = 29.3000; t[i++] = 1.7500; 209 y[i] = 21.4000; t[i++] = 2.2500; 210 y[i] = 29.1750; t[i++] = 1.7500; 211 y[i] = 22.1250; t[i++] = 2.2500; 212 y[i] = 17.5125; t[i++] = 2.7500; 213 y[i] = 14.2500; t[i++] = 3.2500; 214 y[i] = 9.4500; t[i++] = 3.7500; 215 y[i] = 9.1500; t[i++] = 4.2500; 216 y[i] = 7.9125; t[i++] = 4.7500; 217 y[i] = 8.4750; t[i++] = 5.2500; 218 y[i] = 6.1125; t[i++] = 5.7500; 219 y[i] = 80.0000; t[i++] = 0.5000; 220 y[i] = 79.0000; t[i++] = 0.6250; 221 y[i] = 63.8000; t[i++] = 0.7500; 222 y[i] = 57.2000; t[i++] = 0.8750; 223 y[i] = 53.2000; t[i++] = 1.0000; 224 y[i] = 42.5000; t[i++] = 1.2500; 225 y[i] = 26.8000; t[i++] = 1.7500; 226 y[i] = 20.4000; t[i++] = 2.2500; 227 y[i] = 26.8500; t[i++] = 1.7500; 228 y[i] = 21.0000; t[i++] = 2.2500; 229 y[i] = 16.4625; t[i++] = 2.7500; 230 y[i] = 12.5250; t[i++] = 3.2500; 231 y[i] = 10.5375; t[i++] = 3.7500; 232 y[i] = 8.5875; t[i++] = 4.2500; 233 y[i] = 7.1250; t[i++] = 4.7500; 234 y[i] = 6.1125; t[i++] = 5.2500; 235 y[i] = 5.9625; t[i++] = 5.7500; 236 y[i] = 74.1000; t[i++] = 0.5000; 237 y[i] = 67.3000; t[i++] = 0.6250; 238 y[i] = 60.8000; t[i++] = 0.7500; 239 y[i] = 55.5000; t[i++] = 0.8750; 240 y[i] = 50.3000; t[i++] = 1.0000; 241 y[i] = 41.0000; t[i++] = 1.2500; 242 y[i] = 29.4000; t[i++] = 1.7500; 243 y[i] = 20.4000; t[i++] = 2.2500; 244 y[i] = 29.3625; t[i++] = 1.7500; 245 y[i] = 21.1500; t[i++] = 2.2500; 246 y[i] = 16.7625; t[i++] = 2.7500; 247 y[i] = 13.2000; t[i++] = 3.2500; 248 y[i] = 10.8750; t[i++] = 3.7500; 249 y[i] = 8.1750; t[i++] = 4.2500; 250 y[i] = 7.3500; t[i++] = 4.7500; 251 y[i] = 5.9625; t[i++] = 5.2500; 252 y[i] = 5.6250; t[i++] = 5.7500; 253 y[i] = 81.5000; t[i++] = .5000; 254 y[i] = 62.4000; t[i++] = .7500; 255 y[i] = 32.5000; t[i++] = 1.5000; 256 y[i] = 12.4100; t[i++] = 3.0000; 257 y[i] = 13.1200; t[i++] = 3.0000; 258 y[i] = 15.5600; t[i++] = 3.0000; 259 y[i] = 5.6300; t[i++] = 6.0000; 260 y[i] = 78.0000; t[i++] = .5000; 261 y[i] = 59.9000; t[i++] = .7500; 262 y[i] = 33.2000; t[i++] = 1.5000; 263 y[i] = 13.8400; t[i++] = 3.0000; 264 y[i] = 12.7500; t[i++] = 3.0000; 265 y[i] = 14.6200; t[i++] = 3.0000; 266 y[i] = 3.9400; t[i++] = 6.0000; 267 y[i] = 76.8000; t[i++] = .5000; 268 y[i] = 61.0000; t[i++] = .7500; 269 y[i] = 32.9000; t[i++] = 1.5000; 270 y[i] = 13.8700; t[i++] = 3.0000; 271 y[i] = 11.8100; t[i++] = 3.0000; 272 y[i] = 13.3100; t[i++] = 3.0000; 273 y[i] = 5.4400; t[i++] = 6.0000; 274 y[i] = 78.0000; t[i++] = .5000; 275 y[i] = 63.5000; t[i++] = .7500; 276 y[i] = 33.8000; t[i++] = 1.5000; 277 y[i] = 12.5600; t[i++] = 3.0000; 278 y[i] = 5.6300; t[i++] = 6.0000; 279 y[i] = 12.7500; t[i++] = 3.0000; 280 y[i] = 13.1200; t[i++] = 3.0000; 281 y[i] = 5.4400; t[i++] = 6.0000; 282 y[i] = 76.8000; t[i++] = .5000; 283 y[i] = 60.0000; t[i++] = .7500; 284 y[i] = 47.8000; t[i++] = 1.0000; 285 y[i] = 32.0000; t[i++] = 1.5000; 286 y[i] = 22.2000; t[i++] = 2.0000; 287 y[i] = 22.5700; t[i++] = 2.0000; 288 y[i] = 18.8200; t[i++] = 2.5000; 289 y[i] = 13.9500; t[i++] = 3.0000; 290 y[i] = 11.2500; t[i++] = 4.0000; 291 y[i] = 9.0000; t[i++] = 5.0000; 292 y[i] = 6.6700; t[i++] = 6.0000; 293 y[i] = 75.8000; t[i++] = .5000; 294 y[i] = 62.0000; t[i++] = .7500; 295 y[i] = 48.8000; t[i++] = 1.0000; 296 y[i] = 35.2000; t[i++] = 1.5000; 297 y[i] = 20.0000; t[i++] = 2.0000; 298 y[i] = 20.3200; t[i++] = 2.0000; 299 y[i] = 19.3100; t[i++] = 2.5000; 300 y[i] = 12.7500; t[i++] = 3.0000; 301 y[i] = 10.4200; t[i++] = 4.0000; 302 y[i] = 7.3100; t[i++] = 5.0000; 303 y[i] = 7.4200; t[i++] = 6.0000; 304 y[i] = 70.5000; t[i++] = .5000; 305 y[i] = 59.5000; t[i++] = .7500; 306 y[i] = 48.5000; t[i++] = 1.0000; 307 y[i] = 35.8000; t[i++] = 1.5000; 308 y[i] = 21.0000; t[i++] = 2.0000; 309 y[i] = 21.6700; t[i++] = 2.0000; 310 y[i] = 21.0000; t[i++] = 2.5000; 311 y[i] = 15.6400; t[i++] = 3.0000; 312 y[i] = 8.1700; t[i++] = 4.0000; 313 y[i] = 8.5500; t[i++] = 5.0000; 314 y[i] = 10.1200; t[i++] = 6.0000; 315 y[i] = 78.0000; t[i++] = .5000; 316 y[i] = 66.0000; t[i++] = .6250; 317 y[i] = 62.0000; t[i++] = .7500; 318 y[i] = 58.0000; t[i++] = .8750; 319 y[i] = 47.7000; t[i++] = 1.0000; 320 y[i] = 37.8000; t[i++] = 1.2500; 321 y[i] = 20.2000; t[i++] = 2.2500; 322 y[i] = 21.0700; t[i++] = 2.2500; 323 y[i] = 13.8700; t[i++] = 2.7500; 324 y[i] = 9.6700; t[i++] = 3.2500; 325 y[i] = 7.7600; t[i++] = 3.7500; 326 y[i] = 5.4400; t[i++] = 4.2500; 327 y[i] = 4.8700; t[i++] = 4.7500; 328 y[i] = 4.0100; t[i++] = 5.2500; 329 y[i] = 3.7500; t[i++] = 5.7500; 330 y[i] = 24.1900; t[i++] = 3.0000; 331 y[i] = 25.7600; t[i++] = 3.0000; 332 y[i] = 18.0700; t[i++] = 3.0000; 333 y[i] = 11.8100; t[i++] = 3.0000; 334 y[i] = 12.0700; t[i++] = 3.0000; 335 y[i] = 16.1200; t[i++] = 3.0000; 336 y[i] = 70.8000; t[i++] = .5000; 337 y[i] = 54.7000; t[i++] = .7500; 338 y[i] = 48.0000; t[i++] = 1.0000; 339 y[i] = 39.8000; t[i++] = 1.5000; 340 y[i] = 29.8000; t[i++] = 2.0000; 341 y[i] = 23.7000; t[i++] = 2.5000; 342 y[i] = 29.6200; t[i++] = 2.0000; 343 y[i] = 23.8100; t[i++] = 2.5000; 344 y[i] = 17.7000; t[i++] = 3.0000; 345 y[i] = 11.5500; t[i++] = 4.0000; 346 y[i] = 12.0700; t[i++] = 5.0000; 347 y[i] = 8.7400; t[i++] = 6.0000; 348 y[i] = 80.7000; t[i++] = .5000; 349 y[i] = 61.3000; t[i++] = .7500; 350 y[i] = 47.5000; t[i++] = 1.0000; 351 y[i] = 29.0000; t[i++] = 1.5000; 352 y[i] = 24.0000; t[i++] = 2.0000; 353 y[i] = 17.7000; t[i++] = 2.5000; 354 y[i] = 24.5600; t[i++] = 2.0000; 355 y[i] = 18.6700; t[i++] = 2.5000; 356 y[i] = 16.2400; t[i++] = 3.0000; 357 y[i] = 8.7400; t[i++] = 4.0000; 358 y[i] = 7.8700; t[i++] = 5.0000; 359 y[i] = 8.5100; t[i++] = 6.0000; 360 y[i] = 66.7000; t[i++] = .5000; 361 y[i] = 59.2000; t[i++] = .7500; 362 y[i] = 40.8000; t[i++] = 1.0000; 363 y[i] = 30.7000; t[i++] = 1.5000; 364 y[i] = 25.7000; t[i++] = 2.0000; 365 y[i] = 16.3000; t[i++] = 2.5000; 366 y[i] = 25.9900; t[i++] = 2.0000; 367 y[i] = 16.9500; t[i++] = 2.5000; 368 y[i] = 13.3500; t[i++] = 3.0000; 369 y[i] = 8.6200; t[i++] = 4.0000; 370 y[i] = 7.2000; t[i++] = 5.0000; 371 y[i] = 6.6400; t[i++] = 6.0000; 372 y[i] = 13.6900; t[i++] = 3.0000; 373 y[i] = 81.0000; t[i++] = .5000; 374 y[i] = 64.5000; t[i++] = .7500; 375 y[i] = 35.5000; t[i++] = 1.5000; 376 y[i] = 13.3100; t[i++] = 3.0000; 377 y[i] = 4.8700; t[i++] = 6.0000; 378 y[i] = 12.9400; t[i++] = 3.0000; 379 y[i] = 5.0600; t[i++] = 6.0000; 380 y[i] = 15.1900; t[i++] = 3.0000; 381 y[i] = 14.6200; t[i++] = 3.0000; 382 y[i] = 15.6400; t[i++] = 3.0000; 383 y[i] = 25.5000; t[i++] = 1.7500; 384 y[i] = 25.9500; t[i++] = 1.7500; 385 y[i] = 81.7000; t[i++] = .5000; 386 y[i] = 61.6000; t[i++] = .7500; 387 y[i] = 29.8000; t[i++] = 1.7500; 388 y[i] = 29.8100; t[i++] = 1.7500; 389 y[i] = 17.1700; t[i++] = 2.7500; 390 y[i] = 10.3900; t[i++] = 3.7500; 391 y[i] = 28.4000; t[i++] = 1.7500; 392 y[i] = 28.6900; t[i++] = 1.7500; 393 y[i] = 81.3000; t[i++] = .5000; 394 y[i] = 60.9000; t[i++] = .7500; 395 y[i] = 16.6500; t[i++] = 2.7500; 396 y[i] = 10.0500; t[i++] = 3.7500; 397 y[i] = 28.9000; t[i++] = 1.7500; 398 y[i] = 28.9500; t[i++] = 1.7500; 399 PetscFunctionReturn(0); 400 } 401 402 PetscErrorCode TaskWorker(AppCtx *user) 403 { 404 PetscReal x[NPARAMETERS],f = 0.0; 405 PetscMPIInt tag=IDLE_TAG; 406 PetscInt index; 407 MPI_Status status; 408 409 PetscFunctionBegin; 410 /* Send check-in message to rank-0 */ 411 412 PetscCallMPI(MPI_Send(&f,1,MPIU_REAL,0,IDLE_TAG,PETSC_COMM_WORLD)); 413 while (tag != DIE_TAG) { 414 PetscCallMPI(MPI_Recv(x,NPARAMETERS,MPIU_REAL,0,MPI_ANY_TAG,PETSC_COMM_WORLD,&status)); 415 tag = status.MPI_TAG; 416 if (tag == IDLE_TAG) { 417 PetscCallMPI(MPI_Send(&f,1,MPIU_REAL,0,IDLE_TAG,PETSC_COMM_WORLD)); 418 } else if (tag != DIE_TAG) { 419 index = (PetscInt)tag; 420 PetscCall(RunSimulation(x,index,&f,user)); 421 PetscCallMPI(MPI_Send(&f,1,MPIU_REAL,0,tag,PETSC_COMM_WORLD)); 422 } 423 } 424 PetscFunctionReturn(0); 425 } 426 427 PetscErrorCode RunSimulation(PetscReal *x, PetscInt i, PetscReal*f, AppCtx *user) 428 { 429 PetscReal *t = user->t; 430 PetscReal *y = user->y; 431 #if defined(PETSC_USE_REAL_SINGLE) 432 *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 */ 433 #else 434 *f = y[i] - PetscExpScalar(-x[0]*t[i])/(x[1] + x[2]*t[i]); 435 #endif 436 return(0); 437 } 438 439 PetscErrorCode StopWorkers(AppCtx *user) 440 { 441 PetscInt checkedin; 442 MPI_Status status; 443 PetscReal f,x[NPARAMETERS]; 444 445 PetscFunctionBegin; 446 checkedin=0; 447 while (checkedin < user->size-1) { 448 PetscCallMPI(MPI_Recv(&f,1,MPIU_REAL,MPI_ANY_SOURCE,MPI_ANY_TAG,PETSC_COMM_WORLD,&status)); 449 checkedin++; 450 PetscCall(PetscArrayzero(x,NPARAMETERS)); 451 PetscCallMPI(MPI_Send(x,NPARAMETERS,MPIU_REAL,status.MPI_SOURCE,DIE_TAG,PETSC_COMM_WORLD)); 452 } 453 PetscFunctionReturn(0); 454 } 455 456 /*TEST 457 458 build: 459 requires: !complex 460 461 test: 462 nsize: 3 463 requires: !single 464 args: -tao_smonitor -tao_max_it 100 -tao_type pounders -tao_gatol 1.e-5 465 466 TEST*/ 467