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 /* 23c4762a1bSJed Brown User-defined application context - contains data needed by the 24c4762a1bSJed Brown application-provided call-back routines, FormFunctionGradient() 25c4762a1bSJed Brown and FormHessian(). 26c4762a1bSJed Brown */ 27c4762a1bSJed Brown typedef struct { 28c4762a1bSJed Brown PetscInt mx, my; /* discretization in x, y directions */ 29c4762a1bSJed Brown PetscReal *bottom, *top, *left, *right; /* boundary values */ 30c4762a1bSJed Brown DM dm; /* distributed array data structure */ 31c4762a1bSJed Brown Mat H; /* Hessian */ 32c4762a1bSJed Brown } AppCtx; 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* -------- User-defined Routines --------- */ 35c4762a1bSJed Brown 36c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx*); 37c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx*,Vec); 38c4762a1bSJed Brown PetscErrorCode QuadraticH(AppCtx*,Vec,Mat); 39c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal *,Vec,void*); 40c4762a1bSJed Brown PetscErrorCode FormGradient(Tao,Vec,Vec,void*); 41c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*); 42c4762a1bSJed Brown PetscErrorCode My_Monitor(Tao, void *); 43c4762a1bSJed Brown 44c4762a1bSJed Brown int main(int argc, char **argv) 45c4762a1bSJed Brown { 46c4762a1bSJed Brown PetscInt Nx, Ny; /* number of processors in x- and y- directions */ 47c4762a1bSJed Brown Vec x; /* solution, gradient vectors */ 48c4762a1bSJed Brown PetscBool flg, viewmat; /* flags */ 49c4762a1bSJed Brown PetscBool fddefault, fdcoloring; /* flags */ 50c4762a1bSJed Brown Tao tao; /* TAO solver context */ 51c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 52c4762a1bSJed Brown ISColoring iscoloring; 53c4762a1bSJed Brown MatFDColoring matfdcoloring; 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Initialize TAO */ 569566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv,(char *)0,help)); 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* Specify dimension of the problem */ 59c4762a1bSJed Brown user.mx = 10; user.my = 10; 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* Check for any command line arguments that override defaults */ 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg)); 639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg)); 64c4762a1bSJed Brown 659566063dSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD,"\n---- Minimum Surface Area Problem -----\n")); 66*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD,"mx: %" PetscInt_FMT " my: %" PetscInt_FMT " \n\n",user.mx,user.my)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* Let PETSc determine the vector distribution */ 69c4762a1bSJed Brown Nx = PETSC_DECIDE; Ny = PETSC_DECIDE; 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* Create distributed array (DM) to manage parallel grid and vectors */ 729566063dSJacob Faibussowitsch PetscCall(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)); 739566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.dm)); 749566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.dm)); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* Create TAO solver and set desired solution method.*/ 779566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao)); 789566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao,TAOCG)); 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* 81c4762a1bSJed Brown Extract global vector from DA for the vector of variables -- PETSC routine 82c4762a1bSJed Brown Compute the initial solution -- application specific, see below 83c4762a1bSJed Brown Set this vector for use by TAO -- TAO routine 84c4762a1bSJed Brown */ 859566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.dm,&x)); 869566063dSJacob Faibussowitsch PetscCall(MSA_BoundaryConditions(&user)); 879566063dSJacob Faibussowitsch PetscCall(MSA_InitialPoint(&user,x)); 889566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao,x)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* 91c4762a1bSJed Brown Initialize the Application context for use in function evaluations -- application specific, see below. 92c4762a1bSJed Brown Set routines for function and gradient evaluation 93c4762a1bSJed Brown */ 949566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user)); 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* 97c4762a1bSJed Brown Given the command line arguments, calculate the hessian with either the user- 98c4762a1bSJed Brown provided function FormHessian, or the default finite-difference driven Hessian 99c4762a1bSJed Brown functions 100c4762a1bSJed Brown */ 1019566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-tao_fddefault",&fddefault)); 1029566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-tao_fdcoloring",&fdcoloring)); 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* 105c4762a1bSJed Brown Create a matrix data structure to store the Hessian and set 106c4762a1bSJed Brown the Hessian evalution routine. 107c4762a1bSJed Brown Set the matrix structure to be used for Hessian evalutions 108c4762a1bSJed Brown */ 1099566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(user.dm,&user.H)); 1109566063dSJacob Faibussowitsch PetscCall(MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE)); 111c4762a1bSJed Brown 112c4762a1bSJed Brown if (fdcoloring) { 1139566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(user.dm,IS_COLORING_GLOBAL,&iscoloring)); 1149566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(user.H,iscoloring,&matfdcoloring)); 1159566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring)); 1169566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))FormGradient,(void*)&user)); 1179566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFromOptions(matfdcoloring)); 1189566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao,user.H,user.H,TaoDefaultComputeHessianColor,(void *)matfdcoloring)); 119c4762a1bSJed Brown } else if (fddefault) { 1209566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao,user.H,user.H,TaoDefaultComputeHessian,(void *)NULL)); 121c4762a1bSJed Brown } else { 1229566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao,user.H,user.H,FormHessian,(void *)&user)); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* 126c4762a1bSJed Brown If my_monitor option is in command line, then use the user-provided 127c4762a1bSJed Brown monitoring function 128c4762a1bSJed Brown */ 1299566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-my_monitor",&viewmat)); 130c4762a1bSJed Brown if (viewmat) { 1319566063dSJacob Faibussowitsch PetscCall(TaoSetMonitor(tao,My_Monitor,NULL,NULL)); 132c4762a1bSJed Brown } 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* Check for any tao command line options */ 1359566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 1389566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 139c4762a1bSJed Brown 1409566063dSJacob Faibussowitsch PetscCall(TaoView(tao,PETSC_VIEWER_STDOUT_WORLD)); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* Free TAO data structures */ 1439566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Free PETSc data structures */ 1469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.H)); 148c4762a1bSJed Brown if (fdcoloring) { 1499566063dSJacob Faibussowitsch PetscCall(MatFDColoringDestroy(&matfdcoloring)); 150c4762a1bSJed Brown } 1519566063dSJacob Faibussowitsch PetscCall(PetscFree(user.bottom)); 1529566063dSJacob Faibussowitsch PetscCall(PetscFree(user.top)); 1539566063dSJacob Faibussowitsch PetscCall(PetscFree(user.left)); 1549566063dSJacob Faibussowitsch PetscCall(PetscFree(user.right)); 1559566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 1569566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 157b122ec5aSJacob Faibussowitsch return 0; 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown PetscErrorCode FormGradient(Tao tao, Vec X, Vec G,void *userCtx) 161c4762a1bSJed Brown { 162c4762a1bSJed Brown PetscReal fcn; 163c4762a1bSJed Brown 164c4762a1bSJed Brown PetscFunctionBegin; 1659566063dSJacob Faibussowitsch PetscCall(FormFunctionGradient(tao,X,&fcn,G,userCtx)); 166c4762a1bSJed Brown PetscFunctionReturn(0); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown 169c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 170c4762a1bSJed Brown /* FormFunctionGradient - Evaluates the function and corresponding gradient. 171c4762a1bSJed Brown 172c4762a1bSJed Brown Input Parameters: 173c4762a1bSJed Brown . tao - the Tao context 174c4762a1bSJed Brown . XX - input vector 175a82e8c82SStefano Zampini . userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradient() 176c4762a1bSJed Brown 177c4762a1bSJed Brown Output Parameters: 178c4762a1bSJed Brown . fcn - the newly evaluated function 179c4762a1bSJed Brown . GG - vector containing the newly evaluated gradient 180c4762a1bSJed Brown */ 18170a7d78aSStefano Zampini PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn,Vec G,void *userCtx) 18270a7d78aSStefano Zampini { 183c4762a1bSJed Brown AppCtx *user = (AppCtx *) userCtx; 184c4762a1bSJed Brown PetscInt i,j; 185c4762a1bSJed Brown PetscInt mx=user->mx, my=user->my; 186c4762a1bSJed Brown PetscInt xs,xm,gxs,gxm,ys,ym,gys,gym; 187c4762a1bSJed Brown PetscReal ft=0.0; 188c4762a1bSJed Brown PetscReal hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy, area=0.5*hx*hy; 189c4762a1bSJed Brown PetscReal rhx=mx+1, rhy=my+1; 190c4762a1bSJed Brown PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 191c4762a1bSJed Brown PetscReal df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc; 192c4762a1bSJed Brown PetscReal **g, **x; 193c4762a1bSJed Brown Vec localX; 194c4762a1bSJed Brown 195c4762a1bSJed Brown PetscFunctionBegin; 196c4762a1bSJed Brown /* Get local mesh boundaries */ 1979566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->dm,&localX)); 1989566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 1999566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL)); 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* Scatter ghost points to local vector */ 2029566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX)); 2039566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX)); 204c4762a1bSJed Brown 205c4762a1bSJed Brown /* Get pointers to vector data */ 2069566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->dm,localX,(void**)&x)); 2079566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->dm,G,(void**)&g)); 208c4762a1bSJed Brown 209c4762a1bSJed Brown /* Compute function and gradient over the locally owned part of the mesh */ 210c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 211c4762a1bSJed Brown for (i=xs; i< xs+xm; i++) { 212c4762a1bSJed Brown 213c4762a1bSJed Brown xc = x[j][i]; 214c4762a1bSJed Brown xlt=xrb=xl=xr=xb=xt=xc; 215c4762a1bSJed Brown 216c4762a1bSJed Brown if (i==0) { /* left side */ 217c4762a1bSJed Brown xl= user->left[j-ys+1]; 218c4762a1bSJed Brown xlt = user->left[j-ys+2]; 219c4762a1bSJed Brown } else { 220c4762a1bSJed Brown xl = x[j][i-1]; 221c4762a1bSJed Brown } 222c4762a1bSJed Brown 223c4762a1bSJed Brown if (j==0) { /* bottom side */ 224c4762a1bSJed Brown xb=user->bottom[i-xs+1]; 225c4762a1bSJed Brown xrb =user->bottom[i-xs+2]; 226c4762a1bSJed Brown } else { 227c4762a1bSJed Brown xb = x[j-1][i]; 228c4762a1bSJed Brown } 229c4762a1bSJed Brown 230c4762a1bSJed Brown if (i+1 == gxs+gxm) { /* right side */ 231c4762a1bSJed Brown xr=user->right[j-ys+1]; 232c4762a1bSJed Brown xrb = user->right[j-ys]; 233c4762a1bSJed Brown } else { 234c4762a1bSJed Brown xr = x[j][i+1]; 235c4762a1bSJed Brown } 236c4762a1bSJed Brown 237c4762a1bSJed Brown if (j+1==gys+gym) { /* top side */ 238c4762a1bSJed Brown xt=user->top[i-xs+1]; 239c4762a1bSJed Brown xlt = user->top[i-xs]; 240c4762a1bSJed Brown }else { 241c4762a1bSJed Brown xt = x[j+1][i]; 242c4762a1bSJed Brown } 243c4762a1bSJed Brown 244c4762a1bSJed Brown if (i>gxs && j+1<gys+gym) { 245c4762a1bSJed Brown xlt = x[j+1][i-1]; 246c4762a1bSJed Brown } 247c4762a1bSJed Brown if (j>gys && i+1<gxs+gxm) { 248c4762a1bSJed Brown xrb = x[j-1][i+1]; 249c4762a1bSJed Brown } 250c4762a1bSJed Brown 251c4762a1bSJed Brown d1 = (xc-xl); 252c4762a1bSJed Brown d2 = (xc-xr); 253c4762a1bSJed Brown d3 = (xc-xt); 254c4762a1bSJed Brown d4 = (xc-xb); 255c4762a1bSJed Brown d5 = (xr-xrb); 256c4762a1bSJed Brown d6 = (xrb-xb); 257c4762a1bSJed Brown d7 = (xlt-xl); 258c4762a1bSJed Brown d8 = (xt-xlt); 259c4762a1bSJed Brown 260c4762a1bSJed Brown df1dxc = d1*hydhx; 261c4762a1bSJed Brown df2dxc = (d1*hydhx + d4*hxdhy); 262c4762a1bSJed Brown df3dxc = d3*hxdhy; 263c4762a1bSJed Brown df4dxc = (d2*hydhx + d3*hxdhy); 264c4762a1bSJed Brown df5dxc = d2*hydhx; 265c4762a1bSJed Brown df6dxc = d4*hxdhy; 266c4762a1bSJed Brown 267c4762a1bSJed Brown d1 *= rhx; 268c4762a1bSJed Brown d2 *= rhx; 269c4762a1bSJed Brown d3 *= rhy; 270c4762a1bSJed Brown d4 *= rhy; 271c4762a1bSJed Brown d5 *= rhy; 272c4762a1bSJed Brown d6 *= rhx; 273c4762a1bSJed Brown d7 *= rhy; 274c4762a1bSJed Brown d8 *= rhx; 275c4762a1bSJed Brown 276c4762a1bSJed Brown f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7); 277c4762a1bSJed Brown f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4); 278c4762a1bSJed Brown f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8); 279c4762a1bSJed Brown f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2); 280c4762a1bSJed Brown f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5); 281c4762a1bSJed Brown f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6); 282c4762a1bSJed Brown 283c4762a1bSJed Brown ft = ft + (f2 + f4); 284c4762a1bSJed Brown 285c4762a1bSJed Brown df1dxc /= f1; 286c4762a1bSJed Brown df2dxc /= f2; 287c4762a1bSJed Brown df3dxc /= f3; 288c4762a1bSJed Brown df4dxc /= f4; 289c4762a1bSJed Brown df5dxc /= f5; 290c4762a1bSJed Brown df6dxc /= f6; 291c4762a1bSJed Brown 292c4762a1bSJed Brown g[j][i] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc) * 0.5; 293c4762a1bSJed Brown 294c4762a1bSJed Brown } 295c4762a1bSJed Brown } 296c4762a1bSJed Brown 297c4762a1bSJed Brown /* Compute triangular areas along the border of the domain. */ 298c4762a1bSJed Brown if (xs==0) { /* left side */ 299c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 300c4762a1bSJed Brown d3=(user->left[j-ys+1] - user->left[j-ys+2])*rhy; 301c4762a1bSJed Brown d2=(user->left[j-ys+1] - x[j][0]) *rhx; 302c4762a1bSJed Brown ft = ft+PetscSqrtReal(1.0 + d3*d3 + d2*d2); 303c4762a1bSJed Brown } 304c4762a1bSJed Brown } 305c4762a1bSJed Brown if (ys==0) { /* bottom side */ 306c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 307c4762a1bSJed Brown d2=(user->bottom[i+1-xs]-user->bottom[i-xs+2])*rhx; 308c4762a1bSJed Brown d3=(user->bottom[i-xs+1]-x[0][i])*rhy; 309c4762a1bSJed Brown ft = ft+PetscSqrtReal(1.0 + d3*d3 + d2*d2); 310c4762a1bSJed Brown } 311c4762a1bSJed Brown } 312c4762a1bSJed Brown 313c4762a1bSJed Brown if (xs+xm==mx) { /* right side */ 314c4762a1bSJed Brown for (j=ys; j< ys+ym; j++) { 315c4762a1bSJed Brown d1=(x[j][mx-1] - user->right[j-ys+1])*rhx; 316c4762a1bSJed Brown d4=(user->right[j-ys]-user->right[j-ys+1])*rhy; 317c4762a1bSJed Brown ft = ft+PetscSqrtReal(1.0 + d1*d1 + d4*d4); 318c4762a1bSJed Brown } 319c4762a1bSJed Brown } 320c4762a1bSJed Brown if (ys+ym==my) { /* top side */ 321c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 322c4762a1bSJed Brown d1=(x[my-1][i] - user->top[i-xs+1])*rhy; 323c4762a1bSJed Brown d4=(user->top[i-xs+1] - user->top[i-xs])*rhx; 324c4762a1bSJed Brown ft = ft+PetscSqrtReal(1.0 + d1*d1 + d4*d4); 325c4762a1bSJed Brown } 326c4762a1bSJed Brown } 327c4762a1bSJed Brown 328c4762a1bSJed Brown if (ys==0 && xs==0) { 329c4762a1bSJed Brown d1=(user->left[0]-user->left[1])/hy; 330c4762a1bSJed Brown d2=(user->bottom[0]-user->bottom[1])*rhx; 331c4762a1bSJed Brown ft +=PetscSqrtReal(1.0 + d1*d1 + d2*d2); 332c4762a1bSJed Brown } 333c4762a1bSJed Brown if (ys+ym == my && xs+xm == mx) { 334c4762a1bSJed Brown d1=(user->right[ym+1] - user->right[ym])*rhy; 335c4762a1bSJed Brown d2=(user->top[xm+1] - user->top[xm])*rhx; 336c4762a1bSJed Brown ft +=PetscSqrtReal(1.0 + d1*d1 + d2*d2); 337c4762a1bSJed Brown } 338c4762a1bSJed Brown 339c4762a1bSJed Brown ft=ft*area; 3409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&ft,fcn,1,MPIU_REAL,MPIU_SUM,MPI_COMM_WORLD)); 341c4762a1bSJed Brown 342c4762a1bSJed Brown /* Restore vectors */ 3439566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->dm,localX,(void**)&x)); 3449566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->dm,G,(void**)&g)); 345c4762a1bSJed Brown 346c4762a1bSJed Brown /* Scatter values to global vector */ 3479566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->dm,&localX)); 3489566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(67.0*xm*ym)); 349c4762a1bSJed Brown PetscFunctionReturn(0); 350c4762a1bSJed Brown } 351c4762a1bSJed Brown 352c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 353c4762a1bSJed Brown /* 354c4762a1bSJed Brown FormHessian - Evaluates Hessian matrix. 355c4762a1bSJed Brown 356c4762a1bSJed Brown Input Parameters: 357c4762a1bSJed Brown . tao - the Tao context 358c4762a1bSJed Brown . x - input vector 359a82e8c82SStefano Zampini . ptr - optional user-defined context, as set by TaoSetHessian() 360c4762a1bSJed Brown 361c4762a1bSJed Brown Output Parameters: 362c4762a1bSJed Brown . H - Hessian matrix 363c4762a1bSJed Brown . Hpre - optionally different preconditioning matrix 364c4762a1bSJed Brown . flg - flag indicating matrix structure 365c4762a1bSJed Brown 366c4762a1bSJed Brown */ 367c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr) 368c4762a1bSJed Brown { 369c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 370c4762a1bSJed Brown 371c4762a1bSJed Brown PetscFunctionBegin; 372c4762a1bSJed Brown /* Evaluate the Hessian entries*/ 3739566063dSJacob Faibussowitsch PetscCall(QuadraticH(user,X,H)); 374c4762a1bSJed Brown PetscFunctionReturn(0); 375c4762a1bSJed Brown } 376c4762a1bSJed Brown 377c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 378c4762a1bSJed Brown /* 379c4762a1bSJed Brown QuadraticH - Evaluates Hessian matrix. 380c4762a1bSJed Brown 381c4762a1bSJed Brown Input Parameters: 382c4762a1bSJed Brown . user - user-defined context, as set by TaoSetHessian() 383c4762a1bSJed Brown . X - input vector 384c4762a1bSJed Brown 385c4762a1bSJed Brown Output Parameter: 386c4762a1bSJed Brown . H - Hessian matrix 387c4762a1bSJed Brown */ 388c4762a1bSJed Brown PetscErrorCode QuadraticH(AppCtx *user, Vec X, Mat Hessian) 389c4762a1bSJed Brown { 390c4762a1bSJed Brown PetscInt i,j,k; 391c4762a1bSJed Brown PetscInt mx=user->mx, my=user->my; 392c4762a1bSJed Brown PetscInt xs,xm,gxs,gxm,ys,ym,gys,gym; 393c4762a1bSJed Brown PetscReal hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy; 394c4762a1bSJed Brown PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 395c4762a1bSJed Brown PetscReal hl,hr,ht,hb,hc,htl,hbr; 396c4762a1bSJed Brown PetscReal **x, v[7]; 397c4762a1bSJed Brown MatStencil col[7],row; 398c4762a1bSJed Brown Vec localX; 399c4762a1bSJed Brown PetscBool assembled; 400c4762a1bSJed Brown 401c4762a1bSJed Brown PetscFunctionBegin; 402c4762a1bSJed Brown /* Get local mesh boundaries */ 4039566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->dm,&localX)); 404c4762a1bSJed Brown 4059566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 4069566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL)); 407c4762a1bSJed Brown 408c4762a1bSJed Brown /* Scatter ghost points to local vector */ 4099566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX)); 4109566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX)); 411c4762a1bSJed Brown 412c4762a1bSJed Brown /* Get pointers to vector data */ 4139566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->dm,localX,(void**)&x)); 414c4762a1bSJed Brown 415c4762a1bSJed Brown /* Initialize matrix entries to zero */ 4169566063dSJacob Faibussowitsch PetscCall(MatAssembled(Hessian,&assembled)); 4179566063dSJacob Faibussowitsch if (assembled) PetscCall(MatZeroEntries(Hessian)); 418c4762a1bSJed Brown 419c4762a1bSJed Brown /* Set various matrix options */ 4209566063dSJacob Faibussowitsch PetscCall(MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE)); 421c4762a1bSJed Brown 422c4762a1bSJed Brown /* Compute Hessian over the locally owned part of the mesh */ 423c4762a1bSJed Brown 424c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 425c4762a1bSJed Brown 426c4762a1bSJed Brown for (i=xs; i< xs+xm; i++) { 427c4762a1bSJed Brown 428c4762a1bSJed Brown xc = x[j][i]; 429c4762a1bSJed Brown xlt=xrb=xl=xr=xb=xt=xc; 430c4762a1bSJed Brown 431c4762a1bSJed Brown /* Left side */ 432c4762a1bSJed Brown if (i==0) { 433c4762a1bSJed Brown xl = user->left[j-ys+1]; 434c4762a1bSJed Brown xlt = user->left[j-ys+2]; 435c4762a1bSJed Brown } else { 436c4762a1bSJed Brown xl = x[j][i-1]; 437c4762a1bSJed Brown } 438c4762a1bSJed Brown 439c4762a1bSJed Brown if (j==0) { 440c4762a1bSJed Brown xb = user->bottom[i-xs+1]; 441c4762a1bSJed Brown xrb = user->bottom[i-xs+2]; 442c4762a1bSJed Brown } else { 443c4762a1bSJed Brown xb = x[j-1][i]; 444c4762a1bSJed Brown } 445c4762a1bSJed Brown 446c4762a1bSJed Brown if (i+1 == mx) { 447c4762a1bSJed Brown xr = user->right[j-ys+1]; 448c4762a1bSJed Brown xrb = user->right[j-ys]; 449c4762a1bSJed Brown } else { 450c4762a1bSJed Brown xr = x[j][i+1]; 451c4762a1bSJed Brown } 452c4762a1bSJed Brown 453c4762a1bSJed Brown if (j+1==my) { 454c4762a1bSJed Brown xt = user->top[i-xs+1]; 455c4762a1bSJed Brown xlt = user->top[i-xs]; 456c4762a1bSJed Brown }else { 457c4762a1bSJed Brown xt = x[j+1][i]; 458c4762a1bSJed Brown } 459c4762a1bSJed Brown 460c4762a1bSJed Brown if (i>0 && j+1<my) { 461c4762a1bSJed Brown xlt = x[j+1][i-1]; 462c4762a1bSJed Brown } 463c4762a1bSJed Brown if (j>0 && i+1<mx) { 464c4762a1bSJed Brown xrb = x[j-1][i+1]; 465c4762a1bSJed Brown } 466c4762a1bSJed Brown 467c4762a1bSJed Brown d1 = (xc-xl)/hx; 468c4762a1bSJed Brown d2 = (xc-xr)/hx; 469c4762a1bSJed Brown d3 = (xc-xt)/hy; 470c4762a1bSJed Brown d4 = (xc-xb)/hy; 471c4762a1bSJed Brown d5 = (xrb-xr)/hy; 472c4762a1bSJed Brown d6 = (xrb-xb)/hx; 473c4762a1bSJed Brown d7 = (xlt-xl)/hy; 474c4762a1bSJed Brown d8 = (xlt-xt)/hx; 475c4762a1bSJed Brown 476c4762a1bSJed Brown f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7); 477c4762a1bSJed Brown f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4); 478c4762a1bSJed Brown f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8); 479c4762a1bSJed Brown f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2); 480c4762a1bSJed Brown f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5); 481c4762a1bSJed Brown f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6); 482c4762a1bSJed Brown 483c4762a1bSJed Brown hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+ 484c4762a1bSJed Brown (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2); 485c4762a1bSJed Brown hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+ 486c4762a1bSJed Brown (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4); 487c4762a1bSJed Brown ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+ 488c4762a1bSJed Brown (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4); 489c4762a1bSJed Brown hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+ 490c4762a1bSJed Brown (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2); 491c4762a1bSJed Brown 492c4762a1bSJed Brown hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6); 493c4762a1bSJed Brown htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3); 494c4762a1bSJed Brown 495c4762a1bSJed Brown hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) + 496c4762a1bSJed Brown hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) + 497c4762a1bSJed Brown (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) + 498c4762a1bSJed Brown (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4); 499c4762a1bSJed Brown 500c4762a1bSJed Brown hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0; hc/=2.0; 501c4762a1bSJed Brown 502c4762a1bSJed Brown row.j = j; row.i = i; 503c4762a1bSJed Brown k=0; 504c4762a1bSJed Brown if (j>0) { 505c4762a1bSJed Brown v[k]=hb; 506c4762a1bSJed Brown col[k].j = j - 1; col[k].i = i; 507c4762a1bSJed Brown k++; 508c4762a1bSJed Brown } 509c4762a1bSJed Brown 510c4762a1bSJed Brown if (j>0 && i < mx -1) { 511c4762a1bSJed Brown v[k]=hbr; 512c4762a1bSJed Brown col[k].j = j - 1; col[k].i = i+1; 513c4762a1bSJed Brown k++; 514c4762a1bSJed Brown } 515c4762a1bSJed Brown 516c4762a1bSJed Brown if (i>0) { 517c4762a1bSJed Brown v[k]= hl; 518c4762a1bSJed Brown col[k].j = j; col[k].i = i-1; 519c4762a1bSJed Brown k++; 520c4762a1bSJed Brown } 521c4762a1bSJed Brown 522c4762a1bSJed Brown v[k]= hc; 523c4762a1bSJed Brown col[k].j = j; col[k].i = i; 524c4762a1bSJed Brown k++; 525c4762a1bSJed Brown 526c4762a1bSJed Brown if (i < mx-1) { 527c4762a1bSJed Brown v[k]= hr; 528c4762a1bSJed Brown col[k].j = j; col[k].i = i+1; 529c4762a1bSJed Brown k++; 530c4762a1bSJed Brown } 531c4762a1bSJed Brown 532c4762a1bSJed Brown if (i>0 && j < my-1) { 533c4762a1bSJed Brown v[k]= htl; 534c4762a1bSJed Brown col[k].j = j+1; col[k].i = i-1; 535c4762a1bSJed Brown k++; 536c4762a1bSJed Brown } 537c4762a1bSJed Brown 538c4762a1bSJed Brown if (j < my-1) { 539c4762a1bSJed Brown v[k]= ht; 540c4762a1bSJed Brown col[k].j = j+1; col[k].i = i; 541c4762a1bSJed Brown k++; 542c4762a1bSJed Brown } 543c4762a1bSJed Brown 544c4762a1bSJed Brown /* 545c4762a1bSJed Brown Set matrix values using local numbering, which was defined 546c4762a1bSJed Brown earlier, in the main routine. 547c4762a1bSJed Brown */ 5489566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(Hessian,1,&row,k,col,v,INSERT_VALUES)); 549c4762a1bSJed Brown } 550c4762a1bSJed Brown } 551c4762a1bSJed Brown 5529566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->dm,localX,(void**)&x)); 5539566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->dm,&localX)); 554c4762a1bSJed Brown 5559566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY)); 5569566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY)); 557c4762a1bSJed Brown 5589566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(199.0*xm*ym)); 559c4762a1bSJed Brown PetscFunctionReturn(0); 560c4762a1bSJed Brown } 561c4762a1bSJed Brown 562c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 563c4762a1bSJed Brown /* 564c4762a1bSJed Brown MSA_BoundaryConditions - Calculates the boundary conditions for 565c4762a1bSJed Brown the region. 566c4762a1bSJed Brown 567c4762a1bSJed Brown Input Parameter: 568c4762a1bSJed Brown . user - user-defined application context 569c4762a1bSJed Brown 570c4762a1bSJed Brown Output Parameter: 571c4762a1bSJed Brown . user - user-defined application context 572c4762a1bSJed Brown */ 573c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx * user) 574c4762a1bSJed Brown { 575c4762a1bSJed Brown PetscInt i,j,k,limit=0,maxits=5; 576c4762a1bSJed Brown PetscInt xs,ys,xm,ym,gxs,gys,gxm,gym; 577c4762a1bSJed Brown PetscInt mx=user->mx,my=user->my; 578c4762a1bSJed Brown PetscInt bsize=0, lsize=0, tsize=0, rsize=0; 579c4762a1bSJed Brown PetscReal one=1.0, two=2.0, three=3.0, tol=1e-10; 580c4762a1bSJed Brown PetscReal fnorm,det,hx,hy,xt=0,yt=0; 581c4762a1bSJed Brown PetscReal u1,u2,nf1,nf2,njac11,njac12,njac21,njac22; 582c4762a1bSJed Brown PetscReal b=-0.5, t=0.5, l=-0.5, r=0.5; 583c4762a1bSJed Brown PetscReal *boundary; 584c4762a1bSJed Brown PetscBool flg; 585c4762a1bSJed Brown 586c4762a1bSJed Brown PetscFunctionBegin; 587c4762a1bSJed Brown /* Get local mesh boundaries */ 5889566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 5899566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL)); 590c4762a1bSJed Brown 591c4762a1bSJed Brown bsize=xm+2; 592c4762a1bSJed Brown lsize=ym+2; 593c4762a1bSJed Brown rsize=ym+2; 594c4762a1bSJed Brown tsize=xm+2; 595c4762a1bSJed Brown 5969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bsize,&user->bottom)); 5979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tsize,&user->top)); 5989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize,&user->left)); 5999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(rsize,&user->right)); 600c4762a1bSJed Brown 601c4762a1bSJed Brown hx= (r-l)/(mx+1); hy=(t-b)/(my+1); 602c4762a1bSJed Brown 603c4762a1bSJed Brown for (j=0; j<4; j++) { 604c4762a1bSJed Brown if (j==0) { 605c4762a1bSJed Brown yt=b; 606c4762a1bSJed Brown xt=l+hx*xs; 607c4762a1bSJed Brown limit=bsize; 608c4762a1bSJed Brown boundary=user->bottom; 609c4762a1bSJed Brown } else if (j==1) { 610c4762a1bSJed Brown yt=t; 611c4762a1bSJed Brown xt=l+hx*xs; 612c4762a1bSJed Brown limit=tsize; 613c4762a1bSJed Brown boundary=user->top; 614c4762a1bSJed Brown } else if (j==2) { 615c4762a1bSJed Brown yt=b+hy*ys; 616c4762a1bSJed Brown xt=l; 617c4762a1bSJed Brown limit=lsize; 618c4762a1bSJed Brown boundary=user->left; 619c4762a1bSJed Brown } else { /* if (j==3) */ 620c4762a1bSJed Brown yt=b+hy*ys; 621c4762a1bSJed Brown xt=r; 622c4762a1bSJed Brown limit=rsize; 623c4762a1bSJed Brown boundary=user->right; 624c4762a1bSJed Brown } 625c4762a1bSJed Brown 626c4762a1bSJed Brown for (i=0; i<limit; i++) { 627c4762a1bSJed Brown u1=xt; 628c4762a1bSJed Brown u2=-yt; 629c4762a1bSJed Brown for (k=0; k<maxits; k++) { 630c4762a1bSJed Brown nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt; 631c4762a1bSJed Brown nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt; 632c4762a1bSJed Brown fnorm=PetscSqrtReal(nf1*nf1+nf2*nf2); 633c4762a1bSJed Brown if (fnorm <= tol) break; 634c4762a1bSJed Brown njac11=one+u2*u2-u1*u1; 635c4762a1bSJed Brown njac12=two*u1*u2; 636c4762a1bSJed Brown njac21=-two*u1*u2; 637c4762a1bSJed Brown njac22=-one - u1*u1 + u2*u2; 638c4762a1bSJed Brown det = njac11*njac22-njac21*njac12; 639c4762a1bSJed Brown u1 = u1-(njac22*nf1-njac12*nf2)/det; 640c4762a1bSJed Brown u2 = u2-(njac11*nf2-njac21*nf1)/det; 641c4762a1bSJed Brown } 642c4762a1bSJed Brown 643c4762a1bSJed Brown boundary[i]=u1*u1-u2*u2; 644c4762a1bSJed Brown if (j==0 || j==1) { 645c4762a1bSJed Brown xt=xt+hx; 646c4762a1bSJed Brown } else { /* if (j==2 || j==3) */ 647c4762a1bSJed Brown yt=yt+hy; 648c4762a1bSJed Brown } 649c4762a1bSJed Brown 650c4762a1bSJed Brown } 651c4762a1bSJed Brown 652c4762a1bSJed Brown } 653c4762a1bSJed Brown 654c4762a1bSJed Brown /* Scale the boundary if desired */ 655c4762a1bSJed Brown if (1==1) { 656c4762a1bSJed Brown PetscReal scl = 1.0; 657c4762a1bSJed Brown 6589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-bottom",&scl,&flg)); 659c4762a1bSJed Brown if (flg) { 660c4762a1bSJed Brown for (i=0;i<bsize;i++) user->bottom[i]*=scl; 661c4762a1bSJed Brown } 662c4762a1bSJed Brown 6639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-top",&scl,&flg)); 664c4762a1bSJed Brown if (flg) { 665c4762a1bSJed Brown for (i=0;i<tsize;i++) user->top[i]*=scl; 666c4762a1bSJed Brown } 667c4762a1bSJed Brown 6689566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-right",&scl,&flg)); 669c4762a1bSJed Brown if (flg) { 670c4762a1bSJed Brown for (i=0;i<rsize;i++) user->right[i]*=scl; 671c4762a1bSJed Brown } 672c4762a1bSJed Brown 6739566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-left",&scl,&flg)); 674c4762a1bSJed Brown if (flg) { 675c4762a1bSJed Brown for (i=0;i<lsize;i++) user->left[i]*=scl; 676c4762a1bSJed Brown } 677c4762a1bSJed Brown } 678c4762a1bSJed Brown PetscFunctionReturn(0); 679c4762a1bSJed Brown } 680c4762a1bSJed Brown 681c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 682c4762a1bSJed Brown /* 683c4762a1bSJed Brown MSA_InitialPoint - Calculates the initial guess in one of three ways. 684c4762a1bSJed Brown 685c4762a1bSJed Brown Input Parameters: 686c4762a1bSJed Brown . user - user-defined application context 687c4762a1bSJed Brown . X - vector for initial guess 688c4762a1bSJed Brown 689c4762a1bSJed Brown Output Parameters: 690c4762a1bSJed Brown . X - newly computed initial guess 691c4762a1bSJed Brown */ 692c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X) 693c4762a1bSJed Brown { 694c4762a1bSJed Brown PetscInt start2=-1,i,j; 695c4762a1bSJed Brown PetscReal start1=0; 696c4762a1bSJed Brown PetscBool flg1,flg2; 697c4762a1bSJed Brown 698c4762a1bSJed Brown PetscFunctionBegin; 6999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-start",&start1,&flg1)); 7009566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-random",&start2,&flg2)); 701c4762a1bSJed Brown 702c4762a1bSJed Brown if (flg1) { /* The zero vector is reasonable */ 703c4762a1bSJed Brown 7049566063dSJacob Faibussowitsch PetscCall(VecSet(X, start1)); 705c4762a1bSJed Brown 706c4762a1bSJed Brown } else if (flg2 && start2>0) { /* Try a random start between -0.5 and 0.5 */ 707c4762a1bSJed Brown PetscRandom rctx; PetscReal np5=-0.5; 708c4762a1bSJed Brown 7099566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rctx)); 7109566063dSJacob Faibussowitsch PetscCall(VecSetRandom(X, rctx)); 7119566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 7129566063dSJacob Faibussowitsch PetscCall(VecShift(X, np5)); 713c4762a1bSJed Brown 714c4762a1bSJed Brown } else { /* Take an average of the boundary conditions */ 715c4762a1bSJed Brown PetscInt xs,xm,ys,ym; 716c4762a1bSJed Brown PetscInt mx=user->mx,my=user->my; 717c4762a1bSJed Brown PetscReal **x; 718c4762a1bSJed Brown 719c4762a1bSJed Brown /* Get local mesh boundaries */ 7209566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 721c4762a1bSJed Brown 722c4762a1bSJed Brown /* Get pointers to vector data */ 7239566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->dm,X,(void**)&x)); 724c4762a1bSJed Brown 725c4762a1bSJed Brown /* Perform local computations */ 726c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 727c4762a1bSJed Brown for (i=xs; i< xs+xm; i++) { 728c4762a1bSJed 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; 729c4762a1bSJed Brown } 730c4762a1bSJed Brown } 7319566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->dm,X,(void**)&x)); 7329566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(9.0*xm*ym)); 733c4762a1bSJed Brown } 734c4762a1bSJed Brown PetscFunctionReturn(0); 735c4762a1bSJed Brown } 736c4762a1bSJed Brown 737c4762a1bSJed Brown /*-----------------------------------------------------------------------*/ 738c4762a1bSJed Brown PetscErrorCode My_Monitor(Tao tao, void *ctx) 739c4762a1bSJed Brown { 740c4762a1bSJed Brown Vec X; 741c4762a1bSJed Brown 742c4762a1bSJed Brown PetscFunctionBegin; 7439566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao,&X)); 7449566063dSJacob Faibussowitsch PetscCall(VecView(X,PETSC_VIEWER_STDOUT_WORLD)); 745c4762a1bSJed Brown PetscFunctionReturn(0); 746c4762a1bSJed Brown } 747c4762a1bSJed Brown 748c4762a1bSJed Brown /*TEST 749c4762a1bSJed Brown 750c4762a1bSJed Brown build: 751c4762a1bSJed Brown requires: !complex 752c4762a1bSJed Brown 753c4762a1bSJed Brown test: 754c4762a1bSJed Brown args: -tao_smonitor -tao_type lmvm -mx 10 -my 8 -tao_gatol 1.e-3 755c4762a1bSJed Brown requires: !single 756c4762a1bSJed Brown 757c4762a1bSJed Brown test: 758c4762a1bSJed Brown suffix: 2 759c4762a1bSJed Brown nsize: 2 760c4762a1bSJed Brown args: -tao_smonitor -tao_type nls -tao_nls_ksp_max_it 15 -tao_gatol 1.e-4 761c4762a1bSJed Brown filter: grep -v "nls ksp" 762c4762a1bSJed Brown requires: !single 763c4762a1bSJed Brown 764c4762a1bSJed Brown test: 765c4762a1bSJed Brown suffix: 3 766c4762a1bSJed Brown nsize: 3 767c4762a1bSJed Brown args: -tao_smonitor -tao_type cg -tao_cg_type fr -mx 10 -my 10 -tao_gatol 1.e-3 768c4762a1bSJed Brown requires: !single 769c4762a1bSJed Brown 770c4762a1bSJed Brown test: 771c4762a1bSJed Brown suffix: 5 772c4762a1bSJed Brown nsize: 2 773c4762a1bSJed Brown args: -tao_smonitor -tao_type bmrm -mx 10 -my 8 -tao_gatol 1.e-3 774c4762a1bSJed Brown requires: !single 775c4762a1bSJed Brown 776c4762a1bSJed Brown TEST*/ 777