1c4762a1bSJed Brown /* Program usage: mpiexec -n <proc> minsurf2 [-help] [all TAO options] */ 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Include "petsctao.h" so we can use TAO solvers. 5c4762a1bSJed Brown petscdmda.h for distributed array 6c4762a1bSJed Brown */ 7c4762a1bSJed Brown #include <petsctao.h> 8c4762a1bSJed Brown #include <petscdmda.h> 9c4762a1bSJed Brown 10c4762a1bSJed Brown static char help[] = 11c4762a1bSJed Brown "This example demonstrates use of the TAO package to \n\ 12c4762a1bSJed Brown solve an unconstrained minimization problem. This example is based on a \n\ 13c4762a1bSJed Brown problem from the MINPACK-2 test suite. Given a rectangular 2-D domain and \n\ 14c4762a1bSJed Brown boundary values along the edges of the domain, the objective is to find the\n\ 15c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\ 16c4762a1bSJed Brown The command line options are:\n\ 17c4762a1bSJed Brown -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\ 18c4762a1bSJed Brown -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\ 19c4762a1bSJed Brown -start <st>, where <st> =0 for zero vector, <st> >0 for random start, and <st> <0 \n\ 20c4762a1bSJed Brown for an average of the boundary conditions\n\n"; 21c4762a1bSJed Brown 22c4762a1bSJed Brown /*T 23c4762a1bSJed Brown Concepts: TAO^Solving an unconstrained minimization problem 24c4762a1bSJed Brown Routines: TaoCreate(); TaoSetType(); 25a82e8c82SStefano Zampini Routines: TaoSetSolution(); 26a82e8c82SStefano Zampini Routines: TaoSetObjectiveAndGradient(); 27a82e8c82SStefano Zampini Routines: TaoSetHessian(); TaoSetFromOptions(); 28c4762a1bSJed Brown Routines: TaoSetMonitor(); 29c4762a1bSJed Brown Routines: TaoSolve(); TaoView(); 30c4762a1bSJed Brown Routines: TaoDestroy(); 31c4762a1bSJed Brown Processors: n 32c4762a1bSJed Brown T*/ 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* 35c4762a1bSJed Brown User-defined application context - contains data needed by the 36c4762a1bSJed Brown application-provided call-back routines, FormFunctionGradient() 37c4762a1bSJed Brown and FormHessian(). 38c4762a1bSJed Brown */ 39c4762a1bSJed Brown typedef struct { 40c4762a1bSJed Brown PetscInt mx, my; /* discretization in x, y directions */ 41c4762a1bSJed Brown PetscReal *bottom, *top, *left, *right; /* boundary values */ 42c4762a1bSJed Brown DM dm; /* distributed array data structure */ 43c4762a1bSJed Brown Mat H; /* Hessian */ 44c4762a1bSJed Brown } AppCtx; 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* -------- User-defined Routines --------- */ 47c4762a1bSJed Brown 48c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx*); 49c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx*,Vec); 50c4762a1bSJed Brown PetscErrorCode QuadraticH(AppCtx*,Vec,Mat); 51c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal *,Vec,void*); 52c4762a1bSJed Brown PetscErrorCode FormGradient(Tao,Vec,Vec,void*); 53c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*); 54c4762a1bSJed Brown PetscErrorCode My_Monitor(Tao, void *); 55c4762a1bSJed Brown 56c4762a1bSJed Brown int main(int argc, char **argv) 57c4762a1bSJed Brown { 58c4762a1bSJed Brown PetscErrorCode ierr; /* used to check for functions returning nonzeros */ 59c4762a1bSJed Brown PetscInt Nx, Ny; /* number of processors in x- and y- directions */ 60c4762a1bSJed Brown Vec x; /* solution, gradient vectors */ 61c4762a1bSJed Brown PetscBool flg, viewmat; /* flags */ 62c4762a1bSJed Brown PetscBool fddefault, fdcoloring; /* flags */ 63c4762a1bSJed Brown Tao tao; /* TAO solver context */ 64c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 65c4762a1bSJed Brown ISColoring iscoloring; 66c4762a1bSJed Brown MatFDColoring matfdcoloring; 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* Initialize TAO */ 69c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv,(char *)0,help);if (ierr) return ierr; 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* Specify dimension of the problem */ 72c4762a1bSJed Brown user.mx = 10; user.my = 10; 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* Check for any command line arguments that override defaults */ 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg)); 77c4762a1bSJed Brown 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(MPI_COMM_WORLD,"\n---- Minimum Surface Area Problem -----\n")); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(MPI_COMM_WORLD,"mx: %D my: %D \n\n",user.mx,user.my)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* Let PETSc determine the vector distribution */ 82c4762a1bSJed Brown Nx = PETSC_DECIDE; Ny = PETSC_DECIDE; 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* Create distributed array (DM) to manage parallel grid and vectors */ 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,user.mx, user.my,Nx,Ny,1,1,NULL,NULL,&user.dm)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(user.dm)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(user.dm)); 88c4762a1bSJed Brown 89c4762a1bSJed Brown /* Create TAO solver and set desired solution method.*/ 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetType(tao,TAOCG)); 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* 94c4762a1bSJed Brown Extract global vector from DA for the vector of variables -- PETSC routine 95c4762a1bSJed Brown Compute the initial solution -- application specific, see below 96c4762a1bSJed Brown Set this vector for use by TAO -- TAO routine 97c4762a1bSJed Brown */ 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(user.dm,&x)); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(MSA_BoundaryConditions(&user)); 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(MSA_InitialPoint(&user,x)); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetSolution(tao,x)); 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* 104c4762a1bSJed Brown Initialize the Application context for use in function evaluations -- application specific, see below. 105c4762a1bSJed Brown Set routines for function and gradient evaluation 106c4762a1bSJed Brown */ 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user)); 108c4762a1bSJed Brown 109c4762a1bSJed Brown /* 110c4762a1bSJed Brown Given the command line arguments, calculate the hessian with either the user- 111c4762a1bSJed Brown provided function FormHessian, or the default finite-difference driven Hessian 112c4762a1bSJed Brown functions 113c4762a1bSJed Brown */ 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-tao_fddefault",&fddefault)); 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-tao_fdcoloring",&fdcoloring)); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* 118c4762a1bSJed Brown Create a matrix data structure to store the Hessian and set 119c4762a1bSJed Brown the Hessian evalution routine. 120c4762a1bSJed Brown Set the matrix structure to be used for Hessian evalutions 121c4762a1bSJed Brown */ 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(user.dm,&user.H)); 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown if (fdcoloring) { 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateColoring(user.dm,IS_COLORING_GLOBAL,&iscoloring)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatFDColoringCreate(user.H,iscoloring,&matfdcoloring)); 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISColoringDestroy(&iscoloring)); 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))FormGradient,(void*)&user)); 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatFDColoringSetFromOptions(matfdcoloring)); 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetHessian(tao,user.H,user.H,TaoDefaultComputeHessianColor,(void *)matfdcoloring)); 132c4762a1bSJed Brown } else if (fddefault) { 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetHessian(tao,user.H,user.H,TaoDefaultComputeHessian,(void *)NULL)); 134c4762a1bSJed Brown } else { 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetHessian(tao,user.H,user.H,FormHessian,(void *)&user)); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* 139c4762a1bSJed Brown If my_monitor option is in command line, then use the user-provided 140c4762a1bSJed Brown monitoring function 141c4762a1bSJed Brown */ 142*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-my_monitor",&viewmat)); 143c4762a1bSJed Brown if (viewmat) { 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetMonitor(tao,My_Monitor,NULL,NULL)); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* Check for any tao command line options */ 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetFromOptions(tao)); 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSolve(tao)); 152c4762a1bSJed Brown 153*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoView(tao,PETSC_VIEWER_STDOUT_WORLD)); 154c4762a1bSJed Brown 155c4762a1bSJed Brown /* Free TAO data structures */ 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoDestroy(&tao)); 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* Free PETSc data structures */ 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 160*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user.H)); 161c4762a1bSJed Brown if (fdcoloring) { 162*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatFDColoringDestroy(&matfdcoloring)); 163c4762a1bSJed Brown } 164*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user.bottom)); 165*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user.top)); 166*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user.left)); 167*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user.right)); 168*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&user.dm)); 169c4762a1bSJed Brown ierr = PetscFinalize(); 170c4762a1bSJed Brown return ierr; 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173c4762a1bSJed Brown PetscErrorCode FormGradient(Tao tao, Vec X, Vec G,void *userCtx) 174c4762a1bSJed Brown { 175c4762a1bSJed Brown PetscReal fcn; 176c4762a1bSJed Brown 177c4762a1bSJed Brown PetscFunctionBegin; 178*5f80ce2aSJacob Faibussowitsch CHKERRQ(FormFunctionGradient(tao,X,&fcn,G,userCtx)); 179c4762a1bSJed Brown PetscFunctionReturn(0); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 183c4762a1bSJed Brown /* FormFunctionGradient - Evaluates the function and corresponding gradient. 184c4762a1bSJed Brown 185c4762a1bSJed Brown Input Parameters: 186c4762a1bSJed Brown . tao - the Tao context 187c4762a1bSJed Brown . XX - input vector 188a82e8c82SStefano Zampini . userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradient() 189c4762a1bSJed Brown 190c4762a1bSJed Brown Output Parameters: 191c4762a1bSJed Brown . fcn - the newly evaluated function 192c4762a1bSJed Brown . GG - vector containing the newly evaluated gradient 193c4762a1bSJed Brown */ 19470a7d78aSStefano Zampini PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn,Vec G,void *userCtx) 19570a7d78aSStefano Zampini { 196c4762a1bSJed Brown AppCtx *user = (AppCtx *) userCtx; 197c4762a1bSJed Brown PetscInt i,j; 198c4762a1bSJed Brown PetscInt mx=user->mx, my=user->my; 199c4762a1bSJed Brown PetscInt xs,xm,gxs,gxm,ys,ym,gys,gym; 200c4762a1bSJed Brown PetscReal ft=0.0; 201c4762a1bSJed Brown PetscReal hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy, area=0.5*hx*hy; 202c4762a1bSJed Brown PetscReal rhx=mx+1, rhy=my+1; 203c4762a1bSJed Brown PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 204c4762a1bSJed Brown PetscReal df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc; 205c4762a1bSJed Brown PetscReal **g, **x; 206c4762a1bSJed Brown Vec localX; 207c4762a1bSJed Brown 208c4762a1bSJed Brown PetscFunctionBegin; 209c4762a1bSJed Brown /* Get local mesh boundaries */ 210*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(user->dm,&localX)); 211*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 212*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL)); 213c4762a1bSJed Brown 214c4762a1bSJed Brown /* Scatter ghost points to local vector */ 215*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX)); 216*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX)); 217c4762a1bSJed Brown 218c4762a1bSJed Brown /* Get pointers to vector data */ 219*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(user->dm,localX,(void**)&x)); 220*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(user->dm,G,(void**)&g)); 221c4762a1bSJed Brown 222c4762a1bSJed Brown /* Compute function and gradient over the locally owned part of the mesh */ 223c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 224c4762a1bSJed Brown for (i=xs; i< xs+xm; i++) { 225c4762a1bSJed Brown 226c4762a1bSJed Brown xc = x[j][i]; 227c4762a1bSJed Brown xlt=xrb=xl=xr=xb=xt=xc; 228c4762a1bSJed Brown 229c4762a1bSJed Brown if (i==0) { /* left side */ 230c4762a1bSJed Brown xl= user->left[j-ys+1]; 231c4762a1bSJed Brown xlt = user->left[j-ys+2]; 232c4762a1bSJed Brown } else { 233c4762a1bSJed Brown xl = x[j][i-1]; 234c4762a1bSJed Brown } 235c4762a1bSJed Brown 236c4762a1bSJed Brown if (j==0) { /* bottom side */ 237c4762a1bSJed Brown xb=user->bottom[i-xs+1]; 238c4762a1bSJed Brown xrb =user->bottom[i-xs+2]; 239c4762a1bSJed Brown } else { 240c4762a1bSJed Brown xb = x[j-1][i]; 241c4762a1bSJed Brown } 242c4762a1bSJed Brown 243c4762a1bSJed Brown if (i+1 == gxs+gxm) { /* right side */ 244c4762a1bSJed Brown xr=user->right[j-ys+1]; 245c4762a1bSJed Brown xrb = user->right[j-ys]; 246c4762a1bSJed Brown } else { 247c4762a1bSJed Brown xr = x[j][i+1]; 248c4762a1bSJed Brown } 249c4762a1bSJed Brown 250c4762a1bSJed Brown if (j+1==gys+gym) { /* top side */ 251c4762a1bSJed Brown xt=user->top[i-xs+1]; 252c4762a1bSJed Brown xlt = user->top[i-xs]; 253c4762a1bSJed Brown }else { 254c4762a1bSJed Brown xt = x[j+1][i]; 255c4762a1bSJed Brown } 256c4762a1bSJed Brown 257c4762a1bSJed Brown if (i>gxs && j+1<gys+gym) { 258c4762a1bSJed Brown xlt = x[j+1][i-1]; 259c4762a1bSJed Brown } 260c4762a1bSJed Brown if (j>gys && i+1<gxs+gxm) { 261c4762a1bSJed Brown xrb = x[j-1][i+1]; 262c4762a1bSJed Brown } 263c4762a1bSJed Brown 264c4762a1bSJed Brown d1 = (xc-xl); 265c4762a1bSJed Brown d2 = (xc-xr); 266c4762a1bSJed Brown d3 = (xc-xt); 267c4762a1bSJed Brown d4 = (xc-xb); 268c4762a1bSJed Brown d5 = (xr-xrb); 269c4762a1bSJed Brown d6 = (xrb-xb); 270c4762a1bSJed Brown d7 = (xlt-xl); 271c4762a1bSJed Brown d8 = (xt-xlt); 272c4762a1bSJed Brown 273c4762a1bSJed Brown df1dxc = d1*hydhx; 274c4762a1bSJed Brown df2dxc = (d1*hydhx + d4*hxdhy); 275c4762a1bSJed Brown df3dxc = d3*hxdhy; 276c4762a1bSJed Brown df4dxc = (d2*hydhx + d3*hxdhy); 277c4762a1bSJed Brown df5dxc = d2*hydhx; 278c4762a1bSJed Brown df6dxc = d4*hxdhy; 279c4762a1bSJed Brown 280c4762a1bSJed Brown d1 *= rhx; 281c4762a1bSJed Brown d2 *= rhx; 282c4762a1bSJed Brown d3 *= rhy; 283c4762a1bSJed Brown d4 *= rhy; 284c4762a1bSJed Brown d5 *= rhy; 285c4762a1bSJed Brown d6 *= rhx; 286c4762a1bSJed Brown d7 *= rhy; 287c4762a1bSJed Brown d8 *= rhx; 288c4762a1bSJed Brown 289c4762a1bSJed Brown f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7); 290c4762a1bSJed Brown f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4); 291c4762a1bSJed Brown f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8); 292c4762a1bSJed Brown f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2); 293c4762a1bSJed Brown f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5); 294c4762a1bSJed Brown f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6); 295c4762a1bSJed Brown 296c4762a1bSJed Brown ft = ft + (f2 + f4); 297c4762a1bSJed Brown 298c4762a1bSJed Brown df1dxc /= f1; 299c4762a1bSJed Brown df2dxc /= f2; 300c4762a1bSJed Brown df3dxc /= f3; 301c4762a1bSJed Brown df4dxc /= f4; 302c4762a1bSJed Brown df5dxc /= f5; 303c4762a1bSJed Brown df6dxc /= f6; 304c4762a1bSJed Brown 305c4762a1bSJed Brown g[j][i] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc) * 0.5; 306c4762a1bSJed Brown 307c4762a1bSJed Brown } 308c4762a1bSJed Brown } 309c4762a1bSJed Brown 310c4762a1bSJed Brown /* Compute triangular areas along the border of the domain. */ 311c4762a1bSJed Brown if (xs==0) { /* left side */ 312c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 313c4762a1bSJed Brown d3=(user->left[j-ys+1] - user->left[j-ys+2])*rhy; 314c4762a1bSJed Brown d2=(user->left[j-ys+1] - x[j][0]) *rhx; 315c4762a1bSJed Brown ft = ft+PetscSqrtReal(1.0 + d3*d3 + d2*d2); 316c4762a1bSJed Brown } 317c4762a1bSJed Brown } 318c4762a1bSJed Brown if (ys==0) { /* bottom side */ 319c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 320c4762a1bSJed Brown d2=(user->bottom[i+1-xs]-user->bottom[i-xs+2])*rhx; 321c4762a1bSJed Brown d3=(user->bottom[i-xs+1]-x[0][i])*rhy; 322c4762a1bSJed Brown ft = ft+PetscSqrtReal(1.0 + d3*d3 + d2*d2); 323c4762a1bSJed Brown } 324c4762a1bSJed Brown } 325c4762a1bSJed Brown 326c4762a1bSJed Brown if (xs+xm==mx) { /* right side */ 327c4762a1bSJed Brown for (j=ys; j< ys+ym; j++) { 328c4762a1bSJed Brown d1=(x[j][mx-1] - user->right[j-ys+1])*rhx; 329c4762a1bSJed Brown d4=(user->right[j-ys]-user->right[j-ys+1])*rhy; 330c4762a1bSJed Brown ft = ft+PetscSqrtReal(1.0 + d1*d1 + d4*d4); 331c4762a1bSJed Brown } 332c4762a1bSJed Brown } 333c4762a1bSJed Brown if (ys+ym==my) { /* top side */ 334c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 335c4762a1bSJed Brown d1=(x[my-1][i] - user->top[i-xs+1])*rhy; 336c4762a1bSJed Brown d4=(user->top[i-xs+1] - user->top[i-xs])*rhx; 337c4762a1bSJed Brown ft = ft+PetscSqrtReal(1.0 + d1*d1 + d4*d4); 338c4762a1bSJed Brown } 339c4762a1bSJed Brown } 340c4762a1bSJed Brown 341c4762a1bSJed Brown if (ys==0 && xs==0) { 342c4762a1bSJed Brown d1=(user->left[0]-user->left[1])/hy; 343c4762a1bSJed Brown d2=(user->bottom[0]-user->bottom[1])*rhx; 344c4762a1bSJed Brown ft +=PetscSqrtReal(1.0 + d1*d1 + d2*d2); 345c4762a1bSJed Brown } 346c4762a1bSJed Brown if (ys+ym == my && xs+xm == mx) { 347c4762a1bSJed Brown d1=(user->right[ym+1] - user->right[ym])*rhy; 348c4762a1bSJed Brown d2=(user->top[xm+1] - user->top[xm])*rhx; 349c4762a1bSJed Brown ft +=PetscSqrtReal(1.0 + d1*d1 + d2*d2); 350c4762a1bSJed Brown } 351c4762a1bSJed Brown 352c4762a1bSJed Brown ft=ft*area; 353*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(&ft,fcn,1,MPIU_REAL,MPIU_SUM,MPI_COMM_WORLD)); 354c4762a1bSJed Brown 355c4762a1bSJed Brown /* Restore vectors */ 356*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(user->dm,localX,(void**)&x)); 357*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(user->dm,G,(void**)&g)); 358c4762a1bSJed Brown 359c4762a1bSJed Brown /* Scatter values to global vector */ 360*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(user->dm,&localX)); 361*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(67.0*xm*ym)); 362c4762a1bSJed Brown PetscFunctionReturn(0); 363c4762a1bSJed Brown } 364c4762a1bSJed Brown 365c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 366c4762a1bSJed Brown /* 367c4762a1bSJed Brown FormHessian - Evaluates Hessian matrix. 368c4762a1bSJed Brown 369c4762a1bSJed Brown Input Parameters: 370c4762a1bSJed Brown . tao - the Tao context 371c4762a1bSJed Brown . x - input vector 372a82e8c82SStefano Zampini . ptr - optional user-defined context, as set by TaoSetHessian() 373c4762a1bSJed Brown 374c4762a1bSJed Brown Output Parameters: 375c4762a1bSJed Brown . H - Hessian matrix 376c4762a1bSJed Brown . Hpre - optionally different preconditioning matrix 377c4762a1bSJed Brown . flg - flag indicating matrix structure 378c4762a1bSJed Brown 379c4762a1bSJed Brown */ 380c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr) 381c4762a1bSJed Brown { 382c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 383c4762a1bSJed Brown 384c4762a1bSJed Brown PetscFunctionBegin; 385c4762a1bSJed Brown /* Evaluate the Hessian entries*/ 386*5f80ce2aSJacob Faibussowitsch CHKERRQ(QuadraticH(user,X,H)); 387c4762a1bSJed Brown PetscFunctionReturn(0); 388c4762a1bSJed Brown } 389c4762a1bSJed Brown 390c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 391c4762a1bSJed Brown /* 392c4762a1bSJed Brown QuadraticH - Evaluates Hessian matrix. 393c4762a1bSJed Brown 394c4762a1bSJed Brown Input Parameters: 395c4762a1bSJed Brown . user - user-defined context, as set by TaoSetHessian() 396c4762a1bSJed Brown . X - input vector 397c4762a1bSJed Brown 398c4762a1bSJed Brown Output Parameter: 399c4762a1bSJed Brown . H - Hessian matrix 400c4762a1bSJed Brown */ 401c4762a1bSJed Brown PetscErrorCode QuadraticH(AppCtx *user, Vec X, Mat Hessian) 402c4762a1bSJed Brown { 403c4762a1bSJed Brown PetscInt i,j,k; 404c4762a1bSJed Brown PetscInt mx=user->mx, my=user->my; 405c4762a1bSJed Brown PetscInt xs,xm,gxs,gxm,ys,ym,gys,gym; 406c4762a1bSJed Brown PetscReal hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy; 407c4762a1bSJed Brown PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 408c4762a1bSJed Brown PetscReal hl,hr,ht,hb,hc,htl,hbr; 409c4762a1bSJed Brown PetscReal **x, v[7]; 410c4762a1bSJed Brown MatStencil col[7],row; 411c4762a1bSJed Brown Vec localX; 412c4762a1bSJed Brown PetscBool assembled; 413c4762a1bSJed Brown 414c4762a1bSJed Brown PetscFunctionBegin; 415c4762a1bSJed Brown /* Get local mesh boundaries */ 416*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(user->dm,&localX)); 417c4762a1bSJed Brown 418*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 419*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL)); 420c4762a1bSJed Brown 421c4762a1bSJed Brown /* Scatter ghost points to local vector */ 422*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX)); 423*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX)); 424c4762a1bSJed Brown 425c4762a1bSJed Brown /* Get pointers to vector data */ 426*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(user->dm,localX,(void**)&x)); 427c4762a1bSJed Brown 428c4762a1bSJed Brown /* Initialize matrix entries to zero */ 429*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssembled(Hessian,&assembled)); 430*5f80ce2aSJacob Faibussowitsch if (assembled) CHKERRQ(MatZeroEntries(Hessian)); 431c4762a1bSJed Brown 432c4762a1bSJed Brown /* Set various matrix options */ 433*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE)); 434c4762a1bSJed Brown 435c4762a1bSJed Brown /* Compute Hessian over the locally owned part of the mesh */ 436c4762a1bSJed Brown 437c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 438c4762a1bSJed Brown 439c4762a1bSJed Brown for (i=xs; i< xs+xm; i++) { 440c4762a1bSJed Brown 441c4762a1bSJed Brown xc = x[j][i]; 442c4762a1bSJed Brown xlt=xrb=xl=xr=xb=xt=xc; 443c4762a1bSJed Brown 444c4762a1bSJed Brown /* Left side */ 445c4762a1bSJed Brown if (i==0) { 446c4762a1bSJed Brown xl = user->left[j-ys+1]; 447c4762a1bSJed Brown xlt = user->left[j-ys+2]; 448c4762a1bSJed Brown } else { 449c4762a1bSJed Brown xl = x[j][i-1]; 450c4762a1bSJed Brown } 451c4762a1bSJed Brown 452c4762a1bSJed Brown if (j==0) { 453c4762a1bSJed Brown xb = user->bottom[i-xs+1]; 454c4762a1bSJed Brown xrb = user->bottom[i-xs+2]; 455c4762a1bSJed Brown } else { 456c4762a1bSJed Brown xb = x[j-1][i]; 457c4762a1bSJed Brown } 458c4762a1bSJed Brown 459c4762a1bSJed Brown if (i+1 == mx) { 460c4762a1bSJed Brown xr = user->right[j-ys+1]; 461c4762a1bSJed Brown xrb = user->right[j-ys]; 462c4762a1bSJed Brown } else { 463c4762a1bSJed Brown xr = x[j][i+1]; 464c4762a1bSJed Brown } 465c4762a1bSJed Brown 466c4762a1bSJed Brown if (j+1==my) { 467c4762a1bSJed Brown xt = user->top[i-xs+1]; 468c4762a1bSJed Brown xlt = user->top[i-xs]; 469c4762a1bSJed Brown }else { 470c4762a1bSJed Brown xt = x[j+1][i]; 471c4762a1bSJed Brown } 472c4762a1bSJed Brown 473c4762a1bSJed Brown if (i>0 && j+1<my) { 474c4762a1bSJed Brown xlt = x[j+1][i-1]; 475c4762a1bSJed Brown } 476c4762a1bSJed Brown if (j>0 && i+1<mx) { 477c4762a1bSJed Brown xrb = x[j-1][i+1]; 478c4762a1bSJed Brown } 479c4762a1bSJed Brown 480c4762a1bSJed Brown d1 = (xc-xl)/hx; 481c4762a1bSJed Brown d2 = (xc-xr)/hx; 482c4762a1bSJed Brown d3 = (xc-xt)/hy; 483c4762a1bSJed Brown d4 = (xc-xb)/hy; 484c4762a1bSJed Brown d5 = (xrb-xr)/hy; 485c4762a1bSJed Brown d6 = (xrb-xb)/hx; 486c4762a1bSJed Brown d7 = (xlt-xl)/hy; 487c4762a1bSJed Brown d8 = (xlt-xt)/hx; 488c4762a1bSJed Brown 489c4762a1bSJed Brown f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7); 490c4762a1bSJed Brown f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4); 491c4762a1bSJed Brown f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8); 492c4762a1bSJed Brown f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2); 493c4762a1bSJed Brown f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5); 494c4762a1bSJed Brown f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6); 495c4762a1bSJed Brown 496c4762a1bSJed Brown hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+ 497c4762a1bSJed Brown (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2); 498c4762a1bSJed Brown hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+ 499c4762a1bSJed Brown (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4); 500c4762a1bSJed Brown ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+ 501c4762a1bSJed Brown (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4); 502c4762a1bSJed Brown hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+ 503c4762a1bSJed Brown (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2); 504c4762a1bSJed Brown 505c4762a1bSJed Brown hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6); 506c4762a1bSJed Brown htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3); 507c4762a1bSJed Brown 508c4762a1bSJed Brown hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) + 509c4762a1bSJed Brown hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) + 510c4762a1bSJed Brown (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) + 511c4762a1bSJed Brown (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4); 512c4762a1bSJed Brown 513c4762a1bSJed Brown hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0; hc/=2.0; 514c4762a1bSJed Brown 515c4762a1bSJed Brown row.j = j; row.i = i; 516c4762a1bSJed Brown k=0; 517c4762a1bSJed Brown if (j>0) { 518c4762a1bSJed Brown v[k]=hb; 519c4762a1bSJed Brown col[k].j = j - 1; col[k].i = i; 520c4762a1bSJed Brown k++; 521c4762a1bSJed Brown } 522c4762a1bSJed Brown 523c4762a1bSJed Brown if (j>0 && i < mx -1) { 524c4762a1bSJed Brown v[k]=hbr; 525c4762a1bSJed Brown col[k].j = j - 1; col[k].i = i+1; 526c4762a1bSJed Brown k++; 527c4762a1bSJed Brown } 528c4762a1bSJed Brown 529c4762a1bSJed Brown if (i>0) { 530c4762a1bSJed Brown v[k]= hl; 531c4762a1bSJed Brown col[k].j = j; col[k].i = i-1; 532c4762a1bSJed Brown k++; 533c4762a1bSJed Brown } 534c4762a1bSJed Brown 535c4762a1bSJed Brown v[k]= hc; 536c4762a1bSJed Brown col[k].j = j; col[k].i = i; 537c4762a1bSJed Brown k++; 538c4762a1bSJed Brown 539c4762a1bSJed Brown if (i < mx-1) { 540c4762a1bSJed Brown v[k]= hr; 541c4762a1bSJed Brown col[k].j = j; col[k].i = i+1; 542c4762a1bSJed Brown k++; 543c4762a1bSJed Brown } 544c4762a1bSJed Brown 545c4762a1bSJed Brown if (i>0 && j < my-1) { 546c4762a1bSJed Brown v[k]= htl; 547c4762a1bSJed Brown col[k].j = j+1; col[k].i = i-1; 548c4762a1bSJed Brown k++; 549c4762a1bSJed Brown } 550c4762a1bSJed Brown 551c4762a1bSJed Brown if (j < my-1) { 552c4762a1bSJed Brown v[k]= ht; 553c4762a1bSJed Brown col[k].j = j+1; col[k].i = i; 554c4762a1bSJed Brown k++; 555c4762a1bSJed Brown } 556c4762a1bSJed Brown 557c4762a1bSJed Brown /* 558c4762a1bSJed Brown Set matrix values using local numbering, which was defined 559c4762a1bSJed Brown earlier, in the main routine. 560c4762a1bSJed Brown */ 561*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesStencil(Hessian,1,&row,k,col,v,INSERT_VALUES)); 562c4762a1bSJed Brown } 563c4762a1bSJed Brown } 564c4762a1bSJed Brown 565*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(user->dm,localX,(void**)&x)); 566*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(user->dm,&localX)); 567c4762a1bSJed Brown 568*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY)); 569*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY)); 570c4762a1bSJed Brown 571*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(199.0*xm*ym)); 572c4762a1bSJed Brown PetscFunctionReturn(0); 573c4762a1bSJed Brown } 574c4762a1bSJed Brown 575c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 576c4762a1bSJed Brown /* 577c4762a1bSJed Brown MSA_BoundaryConditions - Calculates the boundary conditions for 578c4762a1bSJed Brown the region. 579c4762a1bSJed Brown 580c4762a1bSJed Brown Input Parameter: 581c4762a1bSJed Brown . user - user-defined application context 582c4762a1bSJed Brown 583c4762a1bSJed Brown Output Parameter: 584c4762a1bSJed Brown . user - user-defined application context 585c4762a1bSJed Brown */ 586c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx * user) 587c4762a1bSJed Brown { 588c4762a1bSJed Brown PetscInt i,j,k,limit=0,maxits=5; 589c4762a1bSJed Brown PetscInt xs,ys,xm,ym,gxs,gys,gxm,gym; 590c4762a1bSJed Brown PetscInt mx=user->mx,my=user->my; 591c4762a1bSJed Brown PetscInt bsize=0, lsize=0, tsize=0, rsize=0; 592c4762a1bSJed Brown PetscReal one=1.0, two=2.0, three=3.0, tol=1e-10; 593c4762a1bSJed Brown PetscReal fnorm,det,hx,hy,xt=0,yt=0; 594c4762a1bSJed Brown PetscReal u1,u2,nf1,nf2,njac11,njac12,njac21,njac22; 595c4762a1bSJed Brown PetscReal b=-0.5, t=0.5, l=-0.5, r=0.5; 596c4762a1bSJed Brown PetscReal *boundary; 597c4762a1bSJed Brown PetscBool flg; 598c4762a1bSJed Brown 599c4762a1bSJed Brown PetscFunctionBegin; 600c4762a1bSJed Brown /* Get local mesh boundaries */ 601*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 602*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL)); 603c4762a1bSJed Brown 604c4762a1bSJed Brown bsize=xm+2; 605c4762a1bSJed Brown lsize=ym+2; 606c4762a1bSJed Brown rsize=ym+2; 607c4762a1bSJed Brown tsize=xm+2; 608c4762a1bSJed Brown 609*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(bsize,&user->bottom)); 610*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(tsize,&user->top)); 611*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(lsize,&user->left)); 612*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(rsize,&user->right)); 613c4762a1bSJed Brown 614c4762a1bSJed Brown hx= (r-l)/(mx+1); hy=(t-b)/(my+1); 615c4762a1bSJed Brown 616c4762a1bSJed Brown for (j=0; j<4; j++) { 617c4762a1bSJed Brown if (j==0) { 618c4762a1bSJed Brown yt=b; 619c4762a1bSJed Brown xt=l+hx*xs; 620c4762a1bSJed Brown limit=bsize; 621c4762a1bSJed Brown boundary=user->bottom; 622c4762a1bSJed Brown } else if (j==1) { 623c4762a1bSJed Brown yt=t; 624c4762a1bSJed Brown xt=l+hx*xs; 625c4762a1bSJed Brown limit=tsize; 626c4762a1bSJed Brown boundary=user->top; 627c4762a1bSJed Brown } else if (j==2) { 628c4762a1bSJed Brown yt=b+hy*ys; 629c4762a1bSJed Brown xt=l; 630c4762a1bSJed Brown limit=lsize; 631c4762a1bSJed Brown boundary=user->left; 632c4762a1bSJed Brown } else { /* if (j==3) */ 633c4762a1bSJed Brown yt=b+hy*ys; 634c4762a1bSJed Brown xt=r; 635c4762a1bSJed Brown limit=rsize; 636c4762a1bSJed Brown boundary=user->right; 637c4762a1bSJed Brown } 638c4762a1bSJed Brown 639c4762a1bSJed Brown for (i=0; i<limit; i++) { 640c4762a1bSJed Brown u1=xt; 641c4762a1bSJed Brown u2=-yt; 642c4762a1bSJed Brown for (k=0; k<maxits; k++) { 643c4762a1bSJed Brown nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt; 644c4762a1bSJed Brown nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt; 645c4762a1bSJed Brown fnorm=PetscSqrtReal(nf1*nf1+nf2*nf2); 646c4762a1bSJed Brown if (fnorm <= tol) break; 647c4762a1bSJed Brown njac11=one+u2*u2-u1*u1; 648c4762a1bSJed Brown njac12=two*u1*u2; 649c4762a1bSJed Brown njac21=-two*u1*u2; 650c4762a1bSJed Brown njac22=-one - u1*u1 + u2*u2; 651c4762a1bSJed Brown det = njac11*njac22-njac21*njac12; 652c4762a1bSJed Brown u1 = u1-(njac22*nf1-njac12*nf2)/det; 653c4762a1bSJed Brown u2 = u2-(njac11*nf2-njac21*nf1)/det; 654c4762a1bSJed Brown } 655c4762a1bSJed Brown 656c4762a1bSJed Brown boundary[i]=u1*u1-u2*u2; 657c4762a1bSJed Brown if (j==0 || j==1) { 658c4762a1bSJed Brown xt=xt+hx; 659c4762a1bSJed Brown } else { /* if (j==2 || j==3) */ 660c4762a1bSJed Brown yt=yt+hy; 661c4762a1bSJed Brown } 662c4762a1bSJed Brown 663c4762a1bSJed Brown } 664c4762a1bSJed Brown 665c4762a1bSJed Brown } 666c4762a1bSJed Brown 667c4762a1bSJed Brown /* Scale the boundary if desired */ 668c4762a1bSJed Brown if (1==1) { 669c4762a1bSJed Brown PetscReal scl = 1.0; 670c4762a1bSJed Brown 671*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-bottom",&scl,&flg)); 672c4762a1bSJed Brown if (flg) { 673c4762a1bSJed Brown for (i=0;i<bsize;i++) user->bottom[i]*=scl; 674c4762a1bSJed Brown } 675c4762a1bSJed Brown 676*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-top",&scl,&flg)); 677c4762a1bSJed Brown if (flg) { 678c4762a1bSJed Brown for (i=0;i<tsize;i++) user->top[i]*=scl; 679c4762a1bSJed Brown } 680c4762a1bSJed Brown 681*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-right",&scl,&flg)); 682c4762a1bSJed Brown if (flg) { 683c4762a1bSJed Brown for (i=0;i<rsize;i++) user->right[i]*=scl; 684c4762a1bSJed Brown } 685c4762a1bSJed Brown 686*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-left",&scl,&flg)); 687c4762a1bSJed Brown if (flg) { 688c4762a1bSJed Brown for (i=0;i<lsize;i++) user->left[i]*=scl; 689c4762a1bSJed Brown } 690c4762a1bSJed Brown } 691c4762a1bSJed Brown PetscFunctionReturn(0); 692c4762a1bSJed Brown } 693c4762a1bSJed Brown 694c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 695c4762a1bSJed Brown /* 696c4762a1bSJed Brown MSA_InitialPoint - Calculates the initial guess in one of three ways. 697c4762a1bSJed Brown 698c4762a1bSJed Brown Input Parameters: 699c4762a1bSJed Brown . user - user-defined application context 700c4762a1bSJed Brown . X - vector for initial guess 701c4762a1bSJed Brown 702c4762a1bSJed Brown Output Parameters: 703c4762a1bSJed Brown . X - newly computed initial guess 704c4762a1bSJed Brown */ 705c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X) 706c4762a1bSJed Brown { 707c4762a1bSJed Brown PetscInt start2=-1,i,j; 708c4762a1bSJed Brown PetscReal start1=0; 709c4762a1bSJed Brown PetscBool flg1,flg2; 710c4762a1bSJed Brown 711c4762a1bSJed Brown PetscFunctionBegin; 712*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-start",&start1,&flg1)); 713*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-random",&start2,&flg2)); 714c4762a1bSJed Brown 715c4762a1bSJed Brown if (flg1) { /* The zero vector is reasonable */ 716c4762a1bSJed Brown 717*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(X, start1)); 718c4762a1bSJed Brown 719c4762a1bSJed Brown } else if (flg2 && start2>0) { /* Try a random start between -0.5 and 0.5 */ 720c4762a1bSJed Brown PetscRandom rctx; PetscReal np5=-0.5; 721c4762a1bSJed Brown 722*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rctx)); 723*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(X, rctx)); 724*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rctx)); 725*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecShift(X, np5)); 726c4762a1bSJed Brown 727c4762a1bSJed Brown } else { /* Take an average of the boundary conditions */ 728c4762a1bSJed Brown PetscInt xs,xm,ys,ym; 729c4762a1bSJed Brown PetscInt mx=user->mx,my=user->my; 730c4762a1bSJed Brown PetscReal **x; 731c4762a1bSJed Brown 732c4762a1bSJed Brown /* Get local mesh boundaries */ 733*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 734c4762a1bSJed Brown 735c4762a1bSJed Brown /* Get pointers to vector data */ 736*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(user->dm,X,(void**)&x)); 737c4762a1bSJed Brown 738c4762a1bSJed Brown /* Perform local computations */ 739c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 740c4762a1bSJed Brown for (i=xs; i< xs+xm; i++) { 741c4762a1bSJed Brown x[j][i] = (((j+1)*user->bottom[i-xs+1]+(my-j+1)*user->top[i-xs+1])/(my+2)+((i+1)*user->left[j-ys+1]+(mx-i+1)*user->right[j-ys+1])/(mx+2))/2.0; 742c4762a1bSJed Brown } 743c4762a1bSJed Brown } 744*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(user->dm,X,(void**)&x)); 745*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(9.0*xm*ym)); 746c4762a1bSJed Brown } 747c4762a1bSJed Brown PetscFunctionReturn(0); 748c4762a1bSJed Brown } 749c4762a1bSJed Brown 750c4762a1bSJed Brown /*-----------------------------------------------------------------------*/ 751c4762a1bSJed Brown PetscErrorCode My_Monitor(Tao tao, void *ctx) 752c4762a1bSJed Brown { 753c4762a1bSJed Brown Vec X; 754c4762a1bSJed Brown 755c4762a1bSJed Brown PetscFunctionBegin; 756*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoGetSolution(tao,&X)); 757*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(X,PETSC_VIEWER_STDOUT_WORLD)); 758c4762a1bSJed Brown PetscFunctionReturn(0); 759c4762a1bSJed Brown } 760c4762a1bSJed Brown 761c4762a1bSJed Brown /*TEST 762c4762a1bSJed Brown 763c4762a1bSJed Brown build: 764c4762a1bSJed Brown requires: !complex 765c4762a1bSJed Brown 766c4762a1bSJed Brown test: 767c4762a1bSJed Brown args: -tao_smonitor -tao_type lmvm -mx 10 -my 8 -tao_gatol 1.e-3 768c4762a1bSJed Brown requires: !single 769c4762a1bSJed Brown 770c4762a1bSJed Brown test: 771c4762a1bSJed Brown suffix: 2 772c4762a1bSJed Brown nsize: 2 773c4762a1bSJed Brown args: -tao_smonitor -tao_type nls -tao_nls_ksp_max_it 15 -tao_gatol 1.e-4 774c4762a1bSJed Brown filter: grep -v "nls ksp" 775c4762a1bSJed Brown requires: !single 776c4762a1bSJed Brown 777c4762a1bSJed Brown test: 778c4762a1bSJed Brown suffix: 3 779c4762a1bSJed Brown nsize: 3 780c4762a1bSJed Brown args: -tao_smonitor -tao_type cg -tao_cg_type fr -mx 10 -my 10 -tao_gatol 1.e-3 781c4762a1bSJed Brown requires: !single 782c4762a1bSJed Brown 783c4762a1bSJed Brown test: 784c4762a1bSJed Brown suffix: 5 785c4762a1bSJed Brown nsize: 2 786c4762a1bSJed Brown args: -tao_smonitor -tao_type bmrm -mx 10 -my 8 -tao_gatol 1.e-3 787c4762a1bSJed Brown requires: !single 788c4762a1bSJed Brown 789c4762a1bSJed Brown TEST*/ 790