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