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