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