1 #include <petscdmda.h> 2 #include <petsctao.h> 3 4 static char help[] = 5 "This example demonstrates use of the TAO package to \n\ 6 solve an unconstrained minimization problem. This example is based on a \n\ 7 problem from the MINPACK-2 test suite. Given a rectangular 2-D domain, \n\ 8 boundary values along the edges of the domain, and a plate represented by \n\ 9 lower boundary conditions, the objective is to find the\n\ 10 surface with the minimal area that satisfies the boundary conditions.\n\ 11 The command line options are:\n\ 12 -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\ 13 -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\ 14 -bmx <bxg>, where <bxg> = number of grid points under plate in 1st direction\n\ 15 -bmy <byg>, where <byg> = number of grid points under plate in 2nd direction\n\ 16 -bheight <ht>, where <ht> = height of the plate\n\ 17 -start <st>, where <st> =0 for zero vector, <st> >0 for random start, and <st> <0 \n\ 18 for an average of the boundary conditions\n\n"; 19 20 /*T 21 Concepts: TAO^Solving a bound constrained minimization problem 22 Routines: TaoCreate(); 23 Routines: TaoSetType(); TaoSetObjectiveAndGradient(); 24 Routines: TaoSetHessian(); 25 Routines: TaoSetSolution(); 26 Routines: TaoSetVariableBounds(); 27 Routines: TaoSetFromOptions(); 28 Routines: TaoSolve(); TaoView(); 29 Routines: TaoDestroy(); 30 Processors: n 31 T*/ 32 33 /* 34 User-defined application context - contains data needed by the 35 application-provided call-back routines, FormFunctionGradient(), 36 FormHessian(). 37 */ 38 typedef struct { 39 /* problem parameters */ 40 PetscReal bheight; /* Height of plate under the surface */ 41 PetscInt mx, my; /* discretization in x, y directions */ 42 PetscInt bmx,bmy; /* Size of plate under the surface */ 43 Vec Bottom, Top, Left, Right; /* boundary values */ 44 45 /* Working space */ 46 Vec localX, localV; /* ghosted local vector */ 47 DM dm; /* distributed array data structure */ 48 Mat H; 49 } AppCtx; 50 51 /* -------- User-defined Routines --------- */ 52 53 static PetscErrorCode MSA_BoundaryConditions(AppCtx*); 54 static PetscErrorCode MSA_InitialPoint(AppCtx*,Vec); 55 static PetscErrorCode MSA_Plate(Vec,Vec,void*); 56 PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*); 57 PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*); 58 59 /* For testing matrix free submatrices */ 60 PetscErrorCode MatrixFreeHessian(Tao,Vec,Mat, Mat,void*); 61 PetscErrorCode MyMatMult(Mat,Vec,Vec); 62 63 int main(int argc, char **argv) 64 { 65 PetscErrorCode ierr; /* used to check for functions returning nonzeros */ 66 PetscInt Nx, Ny; /* number of processors in x- and y- directions */ 67 PetscInt m, N; /* number of local and global elements in vectors */ 68 Vec x,xl,xu; /* solution vector and bounds*/ 69 PetscBool flg; /* A return variable when checking for user options */ 70 Tao tao; /* Tao solver context */ 71 ISLocalToGlobalMapping isltog; /* local-to-global mapping object */ 72 Mat H_shell; /* to test matrix-free submatrices */ 73 AppCtx user; /* user-defined work context */ 74 75 /* Initialize PETSc, TAO */ 76 ierr = PetscInitialize(&argc, &argv,(char *)0,help);if (ierr) return ierr; 77 78 /* Specify default dimension of the problem */ 79 user.mx = 10; user.my = 10; user.bheight=0.1; 80 81 /* Check for any command line arguments that override defaults */ 82 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg)); 83 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg)); 84 CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-bheight",&user.bheight,&flg)); 85 86 user.bmx = user.mx/2; user.bmy = user.my/2; 87 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bmx",&user.bmx,&flg)); 88 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bmy",&user.bmy,&flg)); 89 90 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n---- Minimum Surface Area With Plate Problem -----\n")); 91 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"mx:%D, my:%D, bmx:%D, bmy:%D, height:%g\n",user.mx,user.my,user.bmx,user.bmy,(double)user.bheight)); 92 93 /* Calculate any derived values from parameters */ 94 N = user.mx*user.my; 95 96 /* Let Petsc determine the dimensions of the local vectors */ 97 Nx = PETSC_DECIDE; Ny = PETSC_DECIDE; 98 99 /* 100 A two dimensional distributed array will help define this problem, 101 which derives from an elliptic PDE on two dimensional domain. From 102 the distributed array, Create the vectors. 103 */ 104 CHKERRQ(DMDACreate2d(MPI_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,user.mx,user.my,Nx,Ny,1,1,NULL,NULL,&user.dm)); 105 CHKERRQ(DMSetFromOptions(user.dm)); 106 CHKERRQ(DMSetUp(user.dm)); 107 /* 108 Extract global and local vectors from DM; The local vectors are 109 used solely as work space for the evaluation of the function, 110 gradient, and Hessian. Duplicate for remaining vectors that are 111 the same types. 112 */ 113 CHKERRQ(DMCreateGlobalVector(user.dm,&x)); /* Solution */ 114 CHKERRQ(DMCreateLocalVector(user.dm,&user.localX)); 115 CHKERRQ(VecDuplicate(user.localX,&user.localV)); 116 117 CHKERRQ(VecDuplicate(x,&xl)); 118 CHKERRQ(VecDuplicate(x,&xu)); 119 120 /* The TAO code begins here */ 121 122 /* 123 Create TAO solver and set desired solution method 124 The method must either be TAOTRON or TAOBLMVM 125 If TAOBLMVM is used, then hessian function is not called. 126 */ 127 CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao)); 128 CHKERRQ(TaoSetType(tao,TAOBLMVM)); 129 130 /* Set initial solution guess; */ 131 CHKERRQ(MSA_BoundaryConditions(&user)); 132 CHKERRQ(MSA_InitialPoint(&user,x)); 133 CHKERRQ(TaoSetSolution(tao,x)); 134 135 /* Set routines for function, gradient and hessian evaluation */ 136 CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void*) &user)); 137 138 CHKERRQ(VecGetLocalSize(x,&m)); 139 CHKERRQ(MatCreateAIJ(MPI_COMM_WORLD,m,m,N,N,7,NULL,3,NULL,&(user.H))); 140 CHKERRQ(MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE)); 141 142 CHKERRQ(DMGetLocalToGlobalMapping(user.dm,&isltog)); 143 CHKERRQ(MatSetLocalToGlobalMapping(user.H,isltog,isltog)); 144 CHKERRQ(PetscOptionsHasName(NULL,NULL,"-matrixfree",&flg)); 145 if (flg) { 146 CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,m,m,N,N,(void*)&user,&H_shell)); 147 CHKERRQ(MatShellSetOperation(H_shell,MATOP_MULT,(void(*)(void))MyMatMult)); 148 CHKERRQ(MatSetOption(H_shell,MAT_SYMMETRIC,PETSC_TRUE)); 149 CHKERRQ(TaoSetHessian(tao,H_shell,H_shell,MatrixFreeHessian,(void*)&user)); 150 } else { 151 CHKERRQ(TaoSetHessian(tao,user.H,user.H,FormHessian,(void*)&user)); 152 } 153 154 /* Set Variable bounds */ 155 CHKERRQ(MSA_Plate(xl,xu,(void*)&user)); 156 CHKERRQ(TaoSetVariableBounds(tao,xl,xu)); 157 158 /* Check for any tao command line options */ 159 CHKERRQ(TaoSetFromOptions(tao)); 160 161 /* SOLVE THE APPLICATION */ 162 CHKERRQ(TaoSolve(tao)); 163 164 CHKERRQ(TaoView(tao,PETSC_VIEWER_STDOUT_WORLD)); 165 166 /* Free TAO data structures */ 167 CHKERRQ(TaoDestroy(&tao)); 168 169 /* Free PETSc data structures */ 170 CHKERRQ(VecDestroy(&x)); 171 CHKERRQ(VecDestroy(&xl)); 172 CHKERRQ(VecDestroy(&xu)); 173 CHKERRQ(MatDestroy(&user.H)); 174 CHKERRQ(VecDestroy(&user.localX)); 175 CHKERRQ(VecDestroy(&user.localV)); 176 CHKERRQ(VecDestroy(&user.Bottom)); 177 CHKERRQ(VecDestroy(&user.Top)); 178 CHKERRQ(VecDestroy(&user.Left)); 179 CHKERRQ(VecDestroy(&user.Right)); 180 CHKERRQ(DMDestroy(&user.dm)); 181 if (flg) { 182 CHKERRQ(MatDestroy(&H_shell)); 183 } 184 ierr = PetscFinalize(); 185 return ierr; 186 } 187 188 /* FormFunctionGradient - Evaluates f(x) and gradient g(x). 189 190 Input Parameters: 191 . tao - the Tao context 192 . X - input vector 193 . userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradient() 194 195 Output Parameters: 196 . fcn - the function value 197 . G - vector containing the newly evaluated gradient 198 199 Notes: 200 In this case, we discretize the domain and Create triangles. The 201 surface of each triangle is planar, whose surface area can be easily 202 computed. The total surface area is found by sweeping through the grid 203 and computing the surface area of the two triangles that have their 204 right angle at the grid point. The diagonal line segments on the 205 grid that define the triangles run from top left to lower right. 206 The numbering of points starts at the lower left and runs left to 207 right, then bottom to top. 208 */ 209 PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G,void *userCtx) 210 { 211 AppCtx *user = (AppCtx *) userCtx; 212 PetscInt i,j,row; 213 PetscInt mx=user->mx, my=user->my; 214 PetscInt xs,xm,gxs,gxm,ys,ym,gys,gym; 215 PetscReal ft=0; 216 PetscReal zero=0.0; 217 PetscReal hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy, area=0.5*hx*hy; 218 PetscReal rhx=mx+1, rhy=my+1; 219 PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 220 PetscReal df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc; 221 PetscReal *g, *x,*left,*right,*bottom,*top; 222 Vec localX = user->localX, localG = user->localV; 223 224 /* Get local mesh boundaries */ 225 CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 226 CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL)); 227 228 /* Scatter ghost points to local vector */ 229 CHKERRQ(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX)); 230 CHKERRQ(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX)); 231 232 /* Initialize vector to zero */ 233 CHKERRQ(VecSet(localG, zero)); 234 235 /* Get pointers to vector data */ 236 CHKERRQ(VecGetArray(localX,&x)); 237 CHKERRQ(VecGetArray(localG,&g)); 238 CHKERRQ(VecGetArray(user->Top,&top)); 239 CHKERRQ(VecGetArray(user->Bottom,&bottom)); 240 CHKERRQ(VecGetArray(user->Left,&left)); 241 CHKERRQ(VecGetArray(user->Right,&right)); 242 243 /* Compute function over the locally owned part of the mesh */ 244 for (j=ys; j<ys+ym; j++) { 245 for (i=xs; i< xs+xm; i++) { 246 row=(j-gys)*gxm + (i-gxs); 247 248 xc = x[row]; 249 xlt=xrb=xl=xr=xb=xt=xc; 250 251 if (i==0) { /* left side */ 252 xl= left[j-ys+1]; 253 xlt = left[j-ys+2]; 254 } else { 255 xl = x[row-1]; 256 } 257 258 if (j==0) { /* bottom side */ 259 xb=bottom[i-xs+1]; 260 xrb = bottom[i-xs+2]; 261 } else { 262 xb = x[row-gxm]; 263 } 264 265 if (i+1 == gxs+gxm) { /* right side */ 266 xr=right[j-ys+1]; 267 xrb = right[j-ys]; 268 } else { 269 xr = x[row+1]; 270 } 271 272 if (j+1==gys+gym) { /* top side */ 273 xt=top[i-xs+1]; 274 xlt = top[i-xs]; 275 }else { 276 xt = x[row+gxm]; 277 } 278 279 if (i>gxs && j+1<gys+gym) { 280 xlt = x[row-1+gxm]; 281 } 282 if (j>gys && i+1<gxs+gxm) { 283 xrb = x[row+1-gxm]; 284 } 285 286 d1 = (xc-xl); 287 d2 = (xc-xr); 288 d3 = (xc-xt); 289 d4 = (xc-xb); 290 d5 = (xr-xrb); 291 d6 = (xrb-xb); 292 d7 = (xlt-xl); 293 d8 = (xt-xlt); 294 295 df1dxc = d1*hydhx; 296 df2dxc = (d1*hydhx + d4*hxdhy); 297 df3dxc = d3*hxdhy; 298 df4dxc = (d2*hydhx + d3*hxdhy); 299 df5dxc = d2*hydhx; 300 df6dxc = d4*hxdhy; 301 302 d1 *= rhx; 303 d2 *= rhx; 304 d3 *= rhy; 305 d4 *= rhy; 306 d5 *= rhy; 307 d6 *= rhx; 308 d7 *= rhy; 309 d8 *= rhx; 310 311 f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7); 312 f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 313 f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8); 314 f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 315 f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5); 316 f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6); 317 318 ft = ft + (f2 + f4); 319 320 df1dxc /= f1; 321 df2dxc /= f2; 322 df3dxc /= f3; 323 df4dxc /= f4; 324 df5dxc /= f5; 325 df6dxc /= f6; 326 327 g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc) * 0.5; 328 329 } 330 } 331 332 /* Compute triangular areas along the border of the domain. */ 333 if (xs==0) { /* left side */ 334 for (j=ys; j<ys+ym; j++) { 335 d3=(left[j-ys+1] - left[j-ys+2])*rhy; 336 d2=(left[j-ys+1] - x[(j-gys)*gxm])*rhx; 337 ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 338 } 339 } 340 if (ys==0) { /* bottom side */ 341 for (i=xs; i<xs+xm; i++) { 342 d2=(bottom[i+1-xs]-bottom[i-xs+2])*rhx; 343 d3=(bottom[i-xs+1]-x[i-gxs])*rhy; 344 ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 345 } 346 } 347 348 if (xs+xm==mx) { /* right side */ 349 for (j=ys; j< ys+ym; j++) { 350 d1=(x[(j+1-gys)*gxm-1]-right[j-ys+1])*rhx; 351 d4=(right[j-ys]-right[j-ys+1])*rhy; 352 ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 353 } 354 } 355 if (ys+ym==my) { /* top side */ 356 for (i=xs; i<xs+xm; i++) { 357 d1=(x[(gym-1)*gxm + i-gxs] - top[i-xs+1])*rhy; 358 d4=(top[i-xs+1] - top[i-xs])*rhx; 359 ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 360 } 361 } 362 363 if (ys==0 && xs==0) { 364 d1=(left[0]-left[1])*rhy; 365 d2=(bottom[0]-bottom[1])*rhx; 366 ft +=PetscSqrtScalar(1.0 + d1*d1 + d2*d2); 367 } 368 if (ys+ym == my && xs+xm == mx) { 369 d1=(right[ym+1] - right[ym])*rhy; 370 d2=(top[xm+1] - top[xm])*rhx; 371 ft +=PetscSqrtScalar(1.0 + d1*d1 + d2*d2); 372 } 373 374 ft=ft*area; 375 CHKERRMPI(MPI_Allreduce(&ft,fcn,1,MPIU_REAL,MPIU_SUM,MPI_COMM_WORLD)); 376 377 /* Restore vectors */ 378 CHKERRQ(VecRestoreArray(localX,&x)); 379 CHKERRQ(VecRestoreArray(localG,&g)); 380 CHKERRQ(VecRestoreArray(user->Left,&left)); 381 CHKERRQ(VecRestoreArray(user->Top,&top)); 382 CHKERRQ(VecRestoreArray(user->Bottom,&bottom)); 383 CHKERRQ(VecRestoreArray(user->Right,&right)); 384 385 /* Scatter values to global vector */ 386 CHKERRQ(DMLocalToGlobalBegin(user->dm,localG,INSERT_VALUES,G)); 387 CHKERRQ(DMLocalToGlobalEnd(user->dm,localG,INSERT_VALUES,G)); 388 389 CHKERRQ(PetscLogFlops(70.0*xm*ym)); 390 391 return 0; 392 } 393 394 /* ------------------------------------------------------------------- */ 395 /* 396 FormHessian - Evaluates Hessian matrix. 397 398 Input Parameters: 399 . tao - the Tao context 400 . x - input vector 401 . ptr - optional user-defined context, as set by TaoSetHessian() 402 403 Output Parameters: 404 . A - Hessian matrix 405 . B - optionally different preconditioning matrix 406 407 Notes: 408 Due to mesh point reordering with DMs, we must always work 409 with the local mesh points, and then transform them to the new 410 global numbering with the local-to-global mapping. We cannot work 411 directly with the global numbers for the original uniprocessor mesh! 412 413 Two methods are available for imposing this transformation 414 when setting matrix entries: 415 (A) MatSetValuesLocal(), using the local ordering (including 416 ghost points!) 417 - Do the following two steps once, before calling TaoSolve() 418 - Use DMGetISLocalToGlobalMapping() to extract the 419 local-to-global map from the DM 420 - Associate this map with the matrix by calling 421 MatSetLocalToGlobalMapping() 422 - Then set matrix entries using the local ordering 423 by calling MatSetValuesLocal() 424 (B) MatSetValues(), using the global ordering 425 - Use DMGetGlobalIndices() to extract the local-to-global map 426 - Then apply this map explicitly yourself 427 - Set matrix entries using the global ordering by calling 428 MatSetValues() 429 Option (A) seems cleaner/easier in many cases, and is the procedure 430 used in this example. 431 */ 432 PetscErrorCode FormHessian(Tao tao,Vec X,Mat Hptr, Mat Hessian, void *ptr) 433 { 434 AppCtx *user = (AppCtx *) ptr; 435 PetscInt i,j,k,row; 436 PetscInt mx=user->mx, my=user->my; 437 PetscInt xs,xm,gxs,gxm,ys,ym,gys,gym,col[7]; 438 PetscReal hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy; 439 PetscReal rhx=mx+1, rhy=my+1; 440 PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 441 PetscReal hl,hr,ht,hb,hc,htl,hbr; 442 PetscReal *x,*left,*right,*bottom,*top; 443 PetscReal v[7]; 444 Vec localX = user->localX; 445 PetscBool assembled; 446 447 /* Set various matrix options */ 448 CHKERRQ(MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE)); 449 450 /* Initialize matrix entries to zero */ 451 CHKERRQ(MatAssembled(Hessian,&assembled)); 452 if (assembled) CHKERRQ(MatZeroEntries(Hessian)); 453 454 /* Get local mesh boundaries */ 455 CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 456 CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL)); 457 458 /* Scatter ghost points to local vector */ 459 CHKERRQ(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX)); 460 CHKERRQ(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX)); 461 462 /* Get pointers to vector data */ 463 CHKERRQ(VecGetArray(localX,&x)); 464 CHKERRQ(VecGetArray(user->Top,&top)); 465 CHKERRQ(VecGetArray(user->Bottom,&bottom)); 466 CHKERRQ(VecGetArray(user->Left,&left)); 467 CHKERRQ(VecGetArray(user->Right,&right)); 468 469 /* Compute Hessian over the locally owned part of the mesh */ 470 471 for (i=xs; i< xs+xm; i++) { 472 473 for (j=ys; j<ys+ym; j++) { 474 475 row=(j-gys)*gxm + (i-gxs); 476 477 xc = x[row]; 478 xlt=xrb=xl=xr=xb=xt=xc; 479 480 /* Left side */ 481 if (i==gxs) { 482 xl= left[j-ys+1]; 483 xlt = left[j-ys+2]; 484 } else { 485 xl = x[row-1]; 486 } 487 488 if (j==gys) { 489 xb=bottom[i-xs+1]; 490 xrb = bottom[i-xs+2]; 491 } else { 492 xb = x[row-gxm]; 493 } 494 495 if (i+1 == gxs+gxm) { 496 xr=right[j-ys+1]; 497 xrb = right[j-ys]; 498 } else { 499 xr = x[row+1]; 500 } 501 502 if (j+1==gys+gym) { 503 xt=top[i-xs+1]; 504 xlt = top[i-xs]; 505 }else { 506 xt = x[row+gxm]; 507 } 508 509 if (i>gxs && j+1<gys+gym) { 510 xlt = x[row-1+gxm]; 511 } 512 if (j>gys && i+1<gxs+gxm) { 513 xrb = x[row+1-gxm]; 514 } 515 516 d1 = (xc-xl)*rhx; 517 d2 = (xc-xr)*rhx; 518 d3 = (xc-xt)*rhy; 519 d4 = (xc-xb)*rhy; 520 d5 = (xrb-xr)*rhy; 521 d6 = (xrb-xb)*rhx; 522 d7 = (xlt-xl)*rhy; 523 d8 = (xlt-xt)*rhx; 524 525 f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7); 526 f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 527 f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8); 528 f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 529 f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5); 530 f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6); 531 532 hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+ 533 (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2); 534 hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+ 535 (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4); 536 ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+ 537 (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4); 538 hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+ 539 (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2); 540 541 hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6); 542 htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3); 543 544 hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) + 545 hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) + 546 (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) + 547 (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4); 548 549 hl*=0.5; hr*=0.5; ht*=0.5; hb*=0.5; hbr*=0.5; htl*=0.5; hc*=0.5; 550 551 k=0; 552 if (j>0) { 553 v[k]=hb; col[k]=row - gxm; k++; 554 } 555 556 if (j>0 && i < mx -1) { 557 v[k]=hbr; col[k]=row - gxm+1; k++; 558 } 559 560 if (i>0) { 561 v[k]= hl; col[k]=row - 1; k++; 562 } 563 564 v[k]= hc; col[k]=row; k++; 565 566 if (i < mx-1) { 567 v[k]= hr; col[k]=row+1; k++; 568 } 569 570 if (i>0 && j < my-1) { 571 v[k]= htl; col[k] = row+gxm-1; k++; 572 } 573 574 if (j < my-1) { 575 v[k]= ht; col[k] = row+gxm; k++; 576 } 577 578 /* 579 Set matrix values using local numbering, which was defined 580 earlier, in the main routine. 581 */ 582 CHKERRQ(MatSetValuesLocal(Hessian,1,&row,k,col,v,INSERT_VALUES)); 583 584 } 585 } 586 587 /* Restore vectors */ 588 CHKERRQ(VecRestoreArray(localX,&x)); 589 CHKERRQ(VecRestoreArray(user->Left,&left)); 590 CHKERRQ(VecRestoreArray(user->Top,&top)); 591 CHKERRQ(VecRestoreArray(user->Bottom,&bottom)); 592 CHKERRQ(VecRestoreArray(user->Right,&right)); 593 594 /* Assemble the matrix */ 595 CHKERRQ(MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY)); 596 CHKERRQ(MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY)); 597 598 CHKERRQ(PetscLogFlops(199.0*xm*ym)); 599 return 0; 600 } 601 602 /* ------------------------------------------------------------------- */ 603 /* 604 MSA_BoundaryConditions - Calculates the boundary conditions for 605 the region. 606 607 Input Parameter: 608 . user - user-defined application context 609 610 Output Parameter: 611 . user - user-defined application context 612 */ 613 static PetscErrorCode MSA_BoundaryConditions(AppCtx * user) 614 { 615 PetscInt i,j,k,maxits=5,limit=0; 616 PetscInt xs,ys,xm,ym,gxs,gys,gxm,gym; 617 PetscInt mx=user->mx,my=user->my; 618 PetscInt bsize=0, lsize=0, tsize=0, rsize=0; 619 PetscReal one=1.0, two=2.0, three=3.0, scl=1.0, tol=1e-10; 620 PetscReal fnorm,det,hx,hy,xt=0,yt=0; 621 PetscReal u1,u2,nf1,nf2,njac11,njac12,njac21,njac22; 622 PetscReal b=-0.5, t=0.5, l=-0.5, r=0.5; 623 PetscReal *boundary; 624 PetscBool flg; 625 Vec Bottom,Top,Right,Left; 626 627 /* Get local mesh boundaries */ 628 CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 629 CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL)); 630 631 bsize=xm+2; 632 lsize=ym+2; 633 rsize=ym+2; 634 tsize=xm+2; 635 636 CHKERRQ(VecCreateMPI(MPI_COMM_WORLD,bsize,PETSC_DECIDE,&Bottom)); 637 CHKERRQ(VecCreateMPI(MPI_COMM_WORLD,tsize,PETSC_DECIDE,&Top)); 638 CHKERRQ(VecCreateMPI(MPI_COMM_WORLD,lsize,PETSC_DECIDE,&Left)); 639 CHKERRQ(VecCreateMPI(MPI_COMM_WORLD,rsize,PETSC_DECIDE,&Right)); 640 641 user->Top=Top; 642 user->Left=Left; 643 user->Bottom=Bottom; 644 user->Right=Right; 645 646 hx= (r-l)/(mx+1); hy=(t-b)/(my+1); 647 648 for (j=0; j<4; j++) { 649 if (j==0) { 650 yt=b; 651 xt=l+hx*xs; 652 limit=bsize; 653 VecGetArray(Bottom,&boundary); 654 } else if (j==1) { 655 yt=t; 656 xt=l+hx*xs; 657 limit=tsize; 658 VecGetArray(Top,&boundary); 659 } else if (j==2) { 660 yt=b+hy*ys; 661 xt=l; 662 limit=lsize; 663 VecGetArray(Left,&boundary); 664 } else if (j==3) { 665 yt=b+hy*ys; 666 xt=r; 667 limit=rsize; 668 VecGetArray(Right,&boundary); 669 } 670 671 for (i=0; i<limit; i++) { 672 u1=xt; 673 u2=-yt; 674 for (k=0; k<maxits; k++) { 675 nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt; 676 nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt; 677 fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2); 678 if (fnorm <= tol) break; 679 njac11=one+u2*u2-u1*u1; 680 njac12=two*u1*u2; 681 njac21=-two*u1*u2; 682 njac22=-one - u1*u1 + u2*u2; 683 det = njac11*njac22-njac21*njac12; 684 u1 = u1-(njac22*nf1-njac12*nf2)/det; 685 u2 = u2-(njac11*nf2-njac21*nf1)/det; 686 } 687 688 boundary[i]=u1*u1-u2*u2; 689 if (j==0 || j==1) { 690 xt=xt+hx; 691 } else if (j==2 || j==3) { 692 yt=yt+hy; 693 } 694 } 695 if (j==0) { 696 CHKERRQ(VecRestoreArray(Bottom,&boundary)); 697 } else if (j==1) { 698 CHKERRQ(VecRestoreArray(Top,&boundary)); 699 } else if (j==2) { 700 CHKERRQ(VecRestoreArray(Left,&boundary)); 701 } else if (j==3) { 702 CHKERRQ(VecRestoreArray(Right,&boundary)); 703 } 704 } 705 706 /* Scale the boundary if desired */ 707 708 CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-bottom",&scl,&flg)); 709 if (flg) { 710 CHKERRQ(VecScale(Bottom, scl)); 711 } 712 CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-top",&scl,&flg)); 713 if (flg) { 714 CHKERRQ(VecScale(Top, scl)); 715 } 716 CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-right",&scl,&flg)); 717 if (flg) { 718 CHKERRQ(VecScale(Right, scl)); 719 } 720 721 CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-left",&scl,&flg)); 722 if (flg) { 723 CHKERRQ(VecScale(Left, scl)); 724 } 725 return 0; 726 } 727 728 /* ------------------------------------------------------------------- */ 729 /* 730 MSA_Plate - Calculates an obstacle for surface to stretch over. 731 732 Input Parameter: 733 . user - user-defined application context 734 735 Output Parameter: 736 . user - user-defined application context 737 */ 738 static PetscErrorCode MSA_Plate(Vec XL,Vec XU,void *ctx) 739 { 740 AppCtx *user=(AppCtx *)ctx; 741 PetscInt i,j,row; 742 PetscInt xs,ys,xm,ym; 743 PetscInt mx=user->mx, my=user->my, bmy, bmx; 744 PetscReal t1,t2,t3; 745 PetscReal *xl, lb=PETSC_NINFINITY, ub=PETSC_INFINITY; 746 PetscBool cylinder; 747 748 user->bmy = PetscMax(0,user->bmy);user->bmy = PetscMin(my,user->bmy); 749 user->bmx = PetscMax(0,user->bmx);user->bmx = PetscMin(mx,user->bmx); 750 bmy=user->bmy; bmx=user->bmx; 751 752 CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 753 754 CHKERRQ(VecSet(XL, lb)); 755 CHKERRQ(VecSet(XU, ub)); 756 757 CHKERRQ(VecGetArray(XL,&xl)); 758 759 CHKERRQ(PetscOptionsHasName(NULL,NULL,"-cylinder",&cylinder)); 760 /* Compute the optional lower box */ 761 if (cylinder) { 762 for (i=xs; i< xs+xm; i++) { 763 for (j=ys; j<ys+ym; j++) { 764 row=(j-ys)*xm + (i-xs); 765 t1=(2.0*i-mx)*bmy; 766 t2=(2.0*j-my)*bmx; 767 t3=bmx*bmx*bmy*bmy; 768 if (t1*t1 + t2*t2 <= t3) { 769 xl[row] = user->bheight; 770 } 771 } 772 } 773 } else { 774 /* Compute the optional lower box */ 775 for (i=xs; i< xs+xm; i++) { 776 for (j=ys; j<ys+ym; j++) { 777 row=(j-ys)*xm + (i-xs); 778 if (i>=(mx-bmx)/2 && i<mx-(mx-bmx)/2 && 779 j>=(my-bmy)/2 && j<my-(my-bmy)/2) { 780 xl[row] = user->bheight; 781 } 782 } 783 } 784 } 785 CHKERRQ(VecRestoreArray(XL,&xl)); 786 787 return 0; 788 } 789 790 /* ------------------------------------------------------------------- */ 791 /* 792 MSA_InitialPoint - Calculates the initial guess in one of three ways. 793 794 Input Parameters: 795 . user - user-defined application context 796 . X - vector for initial guess 797 798 Output Parameters: 799 . X - newly computed initial guess 800 */ 801 static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X) 802 { 803 PetscInt start=-1,i,j; 804 PetscReal zero=0.0; 805 PetscBool flg; 806 807 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg)); 808 if (flg && start==0) { /* The zero vector is reasonable */ 809 CHKERRQ(VecSet(X, zero)); 810 } else if (flg && start>0) { /* Try a random start between -0.5 and 0.5 */ 811 PetscRandom rctx; PetscReal np5=-0.5; 812 813 CHKERRQ(PetscRandomCreate(MPI_COMM_WORLD,&rctx)); 814 for (i=0; i<start; i++) { 815 CHKERRQ(VecSetRandom(X, rctx)); 816 } 817 CHKERRQ(PetscRandomDestroy(&rctx)); 818 CHKERRQ(VecShift(X, np5)); 819 820 } else { /* Take an average of the boundary conditions */ 821 822 PetscInt row,xs,xm,gxs,gxm,ys,ym,gys,gym; 823 PetscInt mx=user->mx,my=user->my; 824 PetscReal *x,*left,*right,*bottom,*top; 825 Vec localX = user->localX; 826 827 /* Get local mesh boundaries */ 828 CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 829 CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL)); 830 831 /* Get pointers to vector data */ 832 CHKERRQ(VecGetArray(user->Top,&top)); 833 CHKERRQ(VecGetArray(user->Bottom,&bottom)); 834 CHKERRQ(VecGetArray(user->Left,&left)); 835 CHKERRQ(VecGetArray(user->Right,&right)); 836 837 CHKERRQ(VecGetArray(localX,&x)); 838 /* Perform local computations */ 839 for (j=ys; j<ys+ym; j++) { 840 for (i=xs; i< xs+xm; i++) { 841 row=(j-gys)*gxm + (i-gxs); 842 x[row] = ((j+1)*bottom[i-xs+1]/my + (my-j+1)*top[i-xs+1]/(my+2)+(i+1)*left[j-ys+1]/mx + (mx-i+1)*right[j-ys+1]/(mx+2))/2.0; 843 } 844 } 845 846 /* Restore vectors */ 847 CHKERRQ(VecRestoreArray(localX,&x)); 848 849 CHKERRQ(VecRestoreArray(user->Left,&left)); 850 CHKERRQ(VecRestoreArray(user->Top,&top)); 851 CHKERRQ(VecRestoreArray(user->Bottom,&bottom)); 852 CHKERRQ(VecRestoreArray(user->Right,&right)); 853 854 /* Scatter values into global vector */ 855 CHKERRQ(DMLocalToGlobalBegin(user->dm,localX,INSERT_VALUES,X)); 856 CHKERRQ(DMLocalToGlobalEnd(user->dm,localX,INSERT_VALUES,X)); 857 858 } 859 return 0; 860 } 861 862 /* For testing matrix free submatrices */ 863 PetscErrorCode MatrixFreeHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) 864 { 865 AppCtx *user = (AppCtx*)ptr; 866 PetscFunctionBegin; 867 CHKERRQ(FormHessian(tao,x,user->H,user->H,ptr)); 868 PetscFunctionReturn(0); 869 } 870 PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y) 871 { 872 void *ptr; 873 AppCtx *user; 874 PetscFunctionBegin; 875 CHKERRQ(MatShellGetContext(H_shell,&ptr)); 876 user = (AppCtx*)ptr; 877 CHKERRQ(MatMult(user->H,X,Y)); 878 PetscFunctionReturn(0); 879 } 880 881 /*TEST 882 883 build: 884 requires: !complex 885 886 test: 887 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type tron -tao_gatol 1.e-5 888 requires: !single 889 890 test: 891 suffix: 2 892 nsize: 2 893 args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_type blmvm -tao_gatol 1.e-4 894 requires: !single 895 896 test: 897 suffix: 3 898 nsize: 3 899 args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_type tron -tao_gatol 1.e-5 900 requires: !single 901 902 test: 903 suffix: 4 904 nsize: 3 905 args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type mask -tao_type tron -tao_gatol 1.e-5 906 requires: !single 907 908 test: 909 suffix: 5 910 nsize: 3 911 args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -matrixfree -pc_type none -tao_type tron -tao_gatol 1.e-5 912 requires: !single 913 914 test: 915 suffix: 6 916 nsize: 3 917 args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -matrixfree -tao_type blmvm -tao_gatol 1.e-4 918 requires: !single 919 920 test: 921 suffix: 7 922 nsize: 3 923 args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -pc_type none -tao_type gpcg -tao_gatol 1.e-5 924 requires: !single 925 926 test: 927 suffix: 8 928 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4 929 requires: !single 930 931 test: 932 suffix: 9 933 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_gatol 1e-4 934 requires: !single 935 936 test: 937 suffix: 10 938 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5 939 requires: !single 940 941 test: 942 suffix: 11 943 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5 944 requires: !single 945 946 test: 947 suffix: 12 948 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5 949 requires: !single 950 951 test: 952 suffix: 13 953 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 954 requires: !single 955 956 test: 957 suffix: 14 958 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 959 requires: !single 960 961 test: 962 suffix: 15 963 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 964 requires: !single 965 966 test: 967 suffix: 16 968 args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnls 969 requires: !single 970 971 test: 972 suffix: 17 973 args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs 974 requires: !single 975 976 test: 977 suffix: 18 978 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5 -tao_mf_hessian 979 requires: !single 980 981 test: 982 suffix: 19 983 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5 -tao_mf_hessian 984 requires: !single 985 986 test: 987 suffix: 20 988 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5 -tao_mf_hessian 989 990 TEST*/ 991