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