1 #include <petsctao.h> 2 3 static char help[] = 4 "This example demonstrates use of the TAO package to\n\ 5 solve an unconstrained system of equations. This example is based on a\n\ 6 problem from the MINPACK-2 test suite. Given a rectangular 2-D domain and\n\ 7 boundary values along the edges of the domain, the objective is to find the\n\ 8 surface with the minimal area that satisfies the boundary conditions.\n\ 9 This application solves this problem using complimentarity -- We are actually\n\ 10 solving the system (grad f)_i >= 0, if x_i == l_i \n\ 11 (grad f)_i = 0, if l_i < x_i < u_i \n\ 12 (grad f)_i <= 0, if x_i == u_i \n\ 13 where f is the function to be minimized. \n\ 14 \n\ 15 The command line options are:\n\ 16 -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\ 17 -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\ 18 -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise \n\n"; 19 20 /*T 21 Concepts: TAO^Solving a complementarity problem 22 Routines: TaoCreate(); TaoDestroy(); 23 24 Processors: 1 25 T*/ 26 27 /* 28 User-defined application context - contains data needed by the 29 application-provided call-back routines, FormFunctionGradient(), 30 FormHessian(). 31 */ 32 typedef struct { 33 PetscInt mx, my; 34 PetscReal *bottom, *top, *left, *right; 35 } AppCtx; 36 37 /* -------- User-defined Routines --------- */ 38 39 static PetscErrorCode MSA_BoundaryConditions(AppCtx *); 40 static PetscErrorCode MSA_InitialPoint(AppCtx *, Vec); 41 PetscErrorCode FormConstraints(Tao, Vec, Vec, void *); 42 PetscErrorCode FormJacobian(Tao, Vec, Mat, Mat, void *); 43 44 int main(int argc, char **argv) 45 { 46 PetscErrorCode ierr; /* used to check for functions returning nonzeros */ 47 Vec x; /* solution vector */ 48 Vec c; /* Constraints function vector */ 49 Vec xl,xu; /* Bounds on the variables */ 50 PetscBool flg; /* A return variable when checking for user options */ 51 Tao tao; /* TAO solver context */ 52 Mat J; /* Jacobian matrix */ 53 PetscInt N; /* Number of elements in vector */ 54 PetscScalar lb = PETSC_NINFINITY; /* lower bound constant */ 55 PetscScalar ub = PETSC_INFINITY; /* upper bound constant */ 56 AppCtx user; /* user-defined work context */ 57 58 /* Initialize PETSc, TAO */ 59 ierr = PetscInitialize(&argc, &argv, (char *)0, help);if (ierr) return ierr; 60 61 /* Specify default dimension of the problem */ 62 user.mx = 4; user.my = 4; 63 64 /* Check for any command line arguments that override defaults */ 65 CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-mx", &user.mx, &flg)); 66 CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-my", &user.my, &flg)); 67 68 /* Calculate any derived values from parameters */ 69 N = user.mx*user.my; 70 71 CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n---- Minimum Surface Area Problem -----\n")); 72 CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"mx:%D, my:%D\n", user.mx,user.my)); 73 74 /* Create appropriate vectors and matrices */ 75 CHKERRQ(VecCreateSeq(MPI_COMM_SELF, N, &x)); 76 CHKERRQ(VecDuplicate(x, &c)); 77 CHKERRQ(MatCreateSeqAIJ(MPI_COMM_SELF, N, N, 7, NULL, &J)); 78 79 /* The TAO code begins here */ 80 81 /* Create TAO solver and set desired solution method */ 82 CHKERRQ(TaoCreate(PETSC_COMM_SELF,&tao)); 83 CHKERRQ(TaoSetType(tao,TAOSSILS)); 84 85 /* Set data structure */ 86 CHKERRQ(TaoSetSolution(tao, x)); 87 88 /* Set routines for constraints function and Jacobian evaluation */ 89 CHKERRQ(TaoSetConstraintsRoutine(tao, c, FormConstraints, (void *)&user)); 90 CHKERRQ(TaoSetJacobianRoutine(tao, J, J, FormJacobian, (void *)&user)); 91 92 /* Set the variable bounds */ 93 CHKERRQ(MSA_BoundaryConditions(&user)); 94 95 /* Set initial solution guess */ 96 CHKERRQ(MSA_InitialPoint(&user, x)); 97 98 /* Set Bounds on variables */ 99 CHKERRQ(VecDuplicate(x, &xl)); 100 CHKERRQ(VecDuplicate(x, &xu)); 101 CHKERRQ(VecSet(xl, lb)); 102 CHKERRQ(VecSet(xu, ub)); 103 CHKERRQ(TaoSetVariableBounds(tao,xl,xu)); 104 105 /* Check for any tao command line options */ 106 CHKERRQ(TaoSetFromOptions(tao)); 107 108 /* Solve the application */ 109 CHKERRQ(TaoSolve(tao)); 110 111 /* Free Tao data structures */ 112 CHKERRQ(TaoDestroy(&tao)); 113 114 /* Free PETSc data structures */ 115 CHKERRQ(VecDestroy(&x)); 116 CHKERRQ(VecDestroy(&xl)); 117 CHKERRQ(VecDestroy(&xu)); 118 CHKERRQ(VecDestroy(&c)); 119 CHKERRQ(MatDestroy(&J)); 120 121 /* Free user-created data structures */ 122 CHKERRQ(PetscFree(user.bottom)); 123 CHKERRQ(PetscFree(user.top)); 124 CHKERRQ(PetscFree(user.left)); 125 CHKERRQ(PetscFree(user.right)); 126 127 ierr = PetscFinalize(); 128 return ierr; 129 } 130 131 /* -------------------------------------------------------------------- */ 132 133 /* FormConstraints - Evaluates gradient of f. 134 135 Input Parameters: 136 . tao - the TAO_APPLICATION context 137 . X - input vector 138 . ptr - optional user-defined context, as set by TaoSetConstraintsRoutine() 139 140 Output Parameters: 141 . G - vector containing the newly evaluated gradient 142 */ 143 PetscErrorCode FormConstraints(Tao tao, Vec X, Vec G, void *ptr) 144 { 145 AppCtx *user = (AppCtx *) ptr; 146 PetscInt i,j,row; 147 PetscInt mx=user->mx, my=user->my; 148 PetscReal hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy; 149 PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 150 PetscReal df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc; 151 PetscScalar zero=0.0; 152 PetscScalar *g, *x; 153 154 PetscFunctionBegin; 155 /* Initialize vector to zero */ 156 CHKERRQ(VecSet(G, zero)); 157 158 /* Get pointers to vector data */ 159 CHKERRQ(VecGetArray(X, &x)); 160 CHKERRQ(VecGetArray(G, &g)); 161 162 /* Compute function over the locally owned part of the mesh */ 163 for (j=0; j<my; j++) { 164 for (i=0; i< mx; i++) { 165 row= j*mx + i; 166 167 xc = x[row]; 168 xlt=xrb=xl=xr=xb=xt=xc; 169 170 if (i==0) { /* left side */ 171 xl= user->left[j+1]; 172 xlt = user->left[j+2]; 173 } else { 174 xl = x[row-1]; 175 } 176 177 if (j==0) { /* bottom side */ 178 xb=user->bottom[i+1]; 179 xrb = user->bottom[i+2]; 180 } else { 181 xb = x[row-mx]; 182 } 183 184 if (i+1 == mx) { /* right side */ 185 xr=user->right[j+1]; 186 xrb = user->right[j]; 187 } else { 188 xr = x[row+1]; 189 } 190 191 if (j+1==0+my) { /* top side */ 192 xt=user->top[i+1]; 193 xlt = user->top[i]; 194 }else { 195 xt = x[row+mx]; 196 } 197 198 if (i>0 && j+1<my) { 199 xlt = x[row-1+mx]; 200 } 201 if (j>0 && i+1<mx) { 202 xrb = x[row+1-mx]; 203 } 204 205 d1 = (xc-xl); 206 d2 = (xc-xr); 207 d3 = (xc-xt); 208 d4 = (xc-xb); 209 d5 = (xr-xrb); 210 d6 = (xrb-xb); 211 d7 = (xlt-xl); 212 d8 = (xt-xlt); 213 214 df1dxc = d1*hydhx; 215 df2dxc = (d1*hydhx + d4*hxdhy); 216 df3dxc = d3*hxdhy; 217 df4dxc = (d2*hydhx + d3*hxdhy); 218 df5dxc = d2*hydhx; 219 df6dxc = d4*hxdhy; 220 221 d1 /= hx; 222 d2 /= hx; 223 d3 /= hy; 224 d4 /= hy; 225 d5 /= hy; 226 d6 /= hx; 227 d7 /= hy; 228 d8 /= hx; 229 230 f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7); 231 f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 232 f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8); 233 f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 234 f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5); 235 f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6); 236 237 df1dxc /= f1; 238 df2dxc /= f2; 239 df3dxc /= f3; 240 df4dxc /= f4; 241 df5dxc /= f5; 242 df6dxc /= f6; 243 244 g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0; 245 } 246 } 247 248 /* Restore vectors */ 249 CHKERRQ(VecRestoreArray(X, &x)); 250 CHKERRQ(VecRestoreArray(G, &g)); 251 CHKERRQ(PetscLogFlops(67*mx*my)); 252 PetscFunctionReturn(0); 253 } 254 255 /* ------------------------------------------------------------------- */ 256 /* 257 FormJacobian - Evaluates Jacobian matrix. 258 259 Input Parameters: 260 . tao - the TAO_APPLICATION context 261 . X - input vector 262 . ptr - optional user-defined context, as set by TaoSetJacobian() 263 264 Output Parameters: 265 . tH - Jacobian matrix 266 267 */ 268 PetscErrorCode FormJacobian(Tao tao, Vec X, Mat H, Mat tHPre, void *ptr) 269 { 270 AppCtx *user = (AppCtx *) ptr; 271 PetscInt i,j,k,row; 272 PetscInt mx=user->mx, my=user->my; 273 PetscInt col[7]; 274 PetscReal hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy; 275 PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 276 PetscReal hl,hr,ht,hb,hc,htl,hbr; 277 const PetscScalar *x; 278 PetscScalar v[7]; 279 PetscBool assembled; 280 281 /* Set various matrix options */ 282 PetscFunctionBegin; 283 CHKERRQ(MatSetOption(H,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE)); 284 CHKERRQ(MatAssembled(H,&assembled)); 285 if (assembled) CHKERRQ(MatZeroEntries(H)); 286 287 /* Get pointers to vector data */ 288 CHKERRQ(VecGetArrayRead(X, &x)); 289 290 /* Compute Jacobian over the locally owned part of the mesh */ 291 for (i=0; i< mx; i++) { 292 for (j=0; j<my; j++) { 293 row= j*mx + i; 294 295 xc = x[row]; 296 xlt=xrb=xl=xr=xb=xt=xc; 297 298 /* Left side */ 299 if (i==0) { 300 xl = user->left[j+1]; 301 xlt = user->left[j+2]; 302 } else { 303 xl = x[row-1]; 304 } 305 306 if (j==0) { 307 xb = user->bottom[i+1]; 308 xrb = user->bottom[i+2]; 309 } else { 310 xb = x[row-mx]; 311 } 312 313 if (i+1 == mx) { 314 xr = user->right[j+1]; 315 xrb = user->right[j]; 316 } else { 317 xr = x[row+1]; 318 } 319 320 if (j+1==my) { 321 xt = user->top[i+1]; 322 xlt = user->top[i]; 323 }else { 324 xt = x[row+mx]; 325 } 326 327 if (i>0 && j+1<my) { 328 xlt = x[row-1+mx]; 329 } 330 if (j>0 && i+1<mx) { 331 xrb = x[row+1-mx]; 332 } 333 334 d1 = (xc-xl)/hx; 335 d2 = (xc-xr)/hx; 336 d3 = (xc-xt)/hy; 337 d4 = (xc-xb)/hy; 338 d5 = (xrb-xr)/hy; 339 d6 = (xrb-xb)/hx; 340 d7 = (xlt-xl)/hy; 341 d8 = (xlt-xt)/hx; 342 343 f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7); 344 f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 345 f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8); 346 f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 347 f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5); 348 f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6); 349 350 hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+(-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2); 351 hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+(-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4); 352 ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+(-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4); 353 hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+(-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2); 354 355 hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6); 356 htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3); 357 358 hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) + hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) + 359 (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) + (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4); 360 361 hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0; hc/=2.0; 362 363 k=0; 364 if (j>0) { 365 v[k]=hb; col[k]=row - mx; k++; 366 } 367 368 if (j>0 && i < mx -1) { 369 v[k]=hbr; col[k]=row - mx+1; k++; 370 } 371 372 if (i>0) { 373 v[k]= hl; col[k]=row - 1; k++; 374 } 375 376 v[k]= hc; col[k]=row; k++; 377 378 if (i < mx-1) { 379 v[k]= hr; col[k]=row+1; k++; 380 } 381 382 if (i>0 && j < my-1) { 383 v[k]= htl; col[k] = row+mx-1; k++; 384 } 385 386 if (j < my-1) { 387 v[k]= ht; col[k] = row+mx; k++; 388 } 389 390 /* 391 Set matrix values using local numbering, which was defined 392 earlier, in the main routine. 393 */ 394 CHKERRQ(MatSetValues(H,1,&row,k,col,v,INSERT_VALUES)); 395 } 396 } 397 398 /* Restore vectors */ 399 CHKERRQ(VecRestoreArrayRead(X,&x)); 400 401 /* Assemble the matrix */ 402 CHKERRQ(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY)); 403 CHKERRQ(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY)); 404 CHKERRQ(PetscLogFlops(199*mx*my)); 405 PetscFunctionReturn(0); 406 } 407 408 /* ------------------------------------------------------------------- */ 409 /* 410 MSA_BoundaryConditions - Calculates the boundary conditions for 411 the region. 412 413 Input Parameter: 414 . user - user-defined application context 415 416 Output Parameter: 417 . user - user-defined application context 418 */ 419 static PetscErrorCode MSA_BoundaryConditions(AppCtx * user) 420 { 421 PetscInt i,j,k,limit=0,maxits=5; 422 PetscInt mx=user->mx,my=user->my; 423 PetscInt bsize=0, lsize=0, tsize=0, rsize=0; 424 PetscReal one=1.0, two=2.0, three=3.0, tol=1e-10; 425 PetscReal fnorm,det,hx,hy,xt=0,yt=0; 426 PetscReal u1,u2,nf1,nf2,njac11,njac12,njac21,njac22; 427 PetscReal b=-0.5, t=0.5, l=-0.5, r=0.5; 428 PetscReal *boundary; 429 430 PetscFunctionBegin; 431 bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2; 432 433 CHKERRQ(PetscMalloc1(bsize, &user->bottom)); 434 CHKERRQ(PetscMalloc1(tsize, &user->top)); 435 CHKERRQ(PetscMalloc1(lsize, &user->left)); 436 CHKERRQ(PetscMalloc1(rsize, &user->right)); 437 438 hx= (r-l)/(mx+1); hy=(t-b)/(my+1); 439 440 for (j=0; j<4; j++) { 441 if (j==0) { 442 yt=b; 443 xt=l; 444 limit=bsize; 445 boundary=user->bottom; 446 } else if (j==1) { 447 yt=t; 448 xt=l; 449 limit=tsize; 450 boundary=user->top; 451 } else if (j==2) { 452 yt=b; 453 xt=l; 454 limit=lsize; 455 boundary=user->left; 456 } else { /* if (j==3) */ 457 yt=b; 458 xt=r; 459 limit=rsize; 460 boundary=user->right; 461 } 462 463 for (i=0; i<limit; i++) { 464 u1=xt; 465 u2=-yt; 466 for (k=0; k<maxits; k++) { 467 nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt; 468 nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt; 469 fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2); 470 if (fnorm <= tol) break; 471 njac11=one+u2*u2-u1*u1; 472 njac12=two*u1*u2; 473 njac21=-two*u1*u2; 474 njac22=-one - u1*u1 + u2*u2; 475 det = njac11*njac22-njac21*njac12; 476 u1 = u1-(njac22*nf1-njac12*nf2)/det; 477 u2 = u2-(njac11*nf2-njac21*nf1)/det; 478 } 479 480 boundary[i]=u1*u1-u2*u2; 481 if (j==0 || j==1) { 482 xt=xt+hx; 483 } else { /* if (j==2 || j==3) */ 484 yt=yt+hy; 485 } 486 } 487 } 488 PetscFunctionReturn(0); 489 } 490 491 /* ------------------------------------------------------------------- */ 492 /* 493 MSA_InitialPoint - Calculates the initial guess in one of three ways. 494 495 Input Parameters: 496 . user - user-defined application context 497 . X - vector for initial guess 498 499 Output Parameters: 500 . X - newly computed initial guess 501 */ 502 static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X) 503 { 504 PetscInt start=-1,i,j; 505 PetscScalar zero=0.0; 506 PetscBool flg; 507 508 PetscFunctionBegin; 509 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg)); 510 511 if (flg && start==0) { /* The zero vector is reasonable */ 512 CHKERRQ(VecSet(X, zero)); 513 } else { /* Take an average of the boundary conditions */ 514 PetscInt row; 515 PetscInt mx=user->mx,my=user->my; 516 PetscScalar *x; 517 518 /* Get pointers to vector data */ 519 CHKERRQ(VecGetArray(X,&x)); 520 521 /* Perform local computations */ 522 for (j=0; j<my; j++) { 523 for (i=0; i< mx; i++) { 524 row=(j)*mx + (i); 525 x[row] = (((j+1)*user->bottom[i+1]+(my-j+1)*user->top[i+1])/(my+2)+ ((i+1)*user->left[j+1]+(mx-i+1)*user->right[j+1])/(mx+2))/2.0; 526 } 527 } 528 529 /* Restore vectors */ 530 CHKERRQ(VecRestoreArray(X,&x)); 531 } 532 PetscFunctionReturn(0); 533 } 534 535 /*TEST 536 537 build: 538 requires: !complex 539 540 test: 541 args: -tao_monitor -tao_view -tao_type ssils -tao_gttol 1.e-5 542 requires: !single 543 544 test: 545 suffix: 2 546 args: -tao_monitor -tao_view -tao_type ssfls -tao_gttol 1.e-5 547 548 TEST*/ 549