1c4762a1bSJed Brown /* Program usage: mpiexec -n <proc> eptorsion2 [-help] [all TAO options] */ 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* ---------------------------------------------------------------------- 4c4762a1bSJed Brown 5c4762a1bSJed Brown Elastic-plastic torsion problem. 6c4762a1bSJed Brown 7c4762a1bSJed Brown The elastic plastic torsion problem arises from the determination 8c4762a1bSJed Brown of the stress field on an infinitely long cylindrical bar, which is 9c4762a1bSJed Brown equivalent to the solution of the following problem: 10c4762a1bSJed Brown 11c4762a1bSJed Brown min{ .5 * integral(||gradient(v(x))||^2 dx) - C * integral(v(x) dx)} 12c4762a1bSJed Brown 13c4762a1bSJed Brown where C is the torsion angle per unit length. 14c4762a1bSJed Brown 15c4762a1bSJed Brown The uniprocessor version of this code is eptorsion1.c; the Fortran 16c4762a1bSJed Brown version of this code is eptorsion2f.F. 17c4762a1bSJed Brown 18c4762a1bSJed Brown This application solves the problem without calculating hessians 19c4762a1bSJed Brown ---------------------------------------------------------------------- */ 20c4762a1bSJed Brown 21c4762a1bSJed Brown /* 22c4762a1bSJed Brown Include "petsctao.h" so that we can use TAO solvers. Note that this 23c4762a1bSJed Brown file automatically includes files for lower-level support, such as those 24c4762a1bSJed Brown provided by the PETSc library: 25c4762a1bSJed Brown petsc.h - base PETSc routines petscvec.h - vectors 26a5b23f4aSJose E. Roman petscsys.h - system routines petscmat.h - matrices 27c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 28c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 29c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMs) for managing 30c4762a1bSJed Brown the parallel mesh. 31c4762a1bSJed Brown */ 32c4762a1bSJed Brown 33c4762a1bSJed Brown #include <petsctao.h> 34c4762a1bSJed Brown #include <petscdmda.h> 35c4762a1bSJed Brown 36c4762a1bSJed Brown static char help[] = 37c4762a1bSJed Brown "Demonstrates use of the TAO package to solve \n\ 38c4762a1bSJed Brown unconstrained minimization problems in parallel. This example is based on \n\ 39c4762a1bSJed Brown the Elastic-Plastic Torsion (dept) problem from the MINPACK-2 test suite.\n\ 40c4762a1bSJed Brown The command line options are:\n\ 41c4762a1bSJed Brown -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\ 42c4762a1bSJed Brown -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\ 43c4762a1bSJed Brown -par <param>, where <param> = angle of twist per unit length\n\n"; 44c4762a1bSJed Brown 45c4762a1bSJed Brown /*T 46c4762a1bSJed Brown Concepts: TAO^Solving an unconstrained minimization problem 47c4762a1bSJed Brown Routines: TaoCreate(); TaoSetType(); 48a82e8c82SStefano Zampini Routines: TaoSetSolution(); 49a82e8c82SStefano Zampini Routines: TaoSetObjectiveAndGradient(); 50a82e8c82SStefano Zampini Routines: TaoSetHessian(); TaoSetFromOptions(); 51c4762a1bSJed Brown Routines: TaoSolve(); 52c4762a1bSJed Brown Routines: TaoDestroy(); 53c4762a1bSJed Brown Processors: n 54c4762a1bSJed Brown T*/ 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* 57c4762a1bSJed Brown User-defined application context - contains data needed by the 58c4762a1bSJed Brown application-provided call-back routines, FormFunction() and 59c4762a1bSJed Brown FormGradient(). 60c4762a1bSJed Brown */ 61c4762a1bSJed Brown typedef struct { 62c4762a1bSJed Brown /* parameters */ 63c4762a1bSJed Brown PetscInt mx, my; /* global discretization in x- and y-directions */ 64c4762a1bSJed Brown PetscReal param; /* nonlinearity parameter */ 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* work space */ 67c4762a1bSJed Brown Vec localX; /* local vectors */ 68c4762a1bSJed Brown DM dm; /* distributed array data structure */ 69c4762a1bSJed Brown } AppCtx; 70c4762a1bSJed Brown 71c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx*, Vec); 72c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*); 73c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*); 74c4762a1bSJed Brown 75c4762a1bSJed Brown int main(int argc, char **argv) 76c4762a1bSJed Brown { 77c4762a1bSJed Brown Vec x; 78c4762a1bSJed Brown Mat H; 79c4762a1bSJed Brown PetscInt Nx, Ny; 80c4762a1bSJed Brown Tao tao; 81c4762a1bSJed Brown PetscBool flg; 82c4762a1bSJed Brown KSP ksp; 83c4762a1bSJed Brown PC pc; 84c4762a1bSJed Brown AppCtx user; 85c4762a1bSJed Brown 86c4762a1bSJed Brown PetscInitialize(&argc, &argv, (char *)0, help); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* Specify default dimension of the problem */ 89c4762a1bSJed Brown user.param = 5.0; user.mx = 10; user.my = 10; 90c4762a1bSJed Brown Nx = Ny = PETSC_DECIDE; 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* Check for any command line arguments that override defaults */ 93*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-par",&user.param,&flg)); 94*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg)); 95*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg)); 96c4762a1bSJed Brown 97*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n---- Elastic-Plastic Torsion Problem -----\n")); 98*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mx: %D my: %D \n\n",user.mx,user.my)); 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* Set up distributed array */ 101*9566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.mx,user.my,Nx,Ny,1,1,NULL,NULL,&user.dm)); 102*9566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.dm)); 103*9566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.dm)); 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* Create vectors */ 106*9566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.dm,&x)); 107c4762a1bSJed Brown 108*9566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(user.dm,&user.localX)); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* Create Hessian */ 111*9566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(user.dm,&H)); 112*9566063dSJacob Faibussowitsch PetscCall(MatSetOption(H,MAT_SYMMETRIC,PETSC_TRUE)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* The TAO code begins here */ 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 117*9566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao)); 118*9566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao,TAOCG)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* Set initial solution guess */ 121*9566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(&user,x)); 122*9566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao,x)); 123c4762a1bSJed Brown 124c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 125*9566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user)); 126c4762a1bSJed Brown 127*9566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao,H,H,FormHessian,(void*)&user)); 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* Check for any TAO command line options */ 130*9566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 131c4762a1bSJed Brown 132*9566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao,&ksp)); 133c4762a1bSJed Brown if (ksp) { 134*9566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 135*9566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCNONE)); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 139*9566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* Free TAO data structures */ 142*9566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* Free PETSc data structures */ 145*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 146*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&H)); 147c4762a1bSJed Brown 148*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.localX)); 149*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 150c4762a1bSJed Brown 151c4762a1bSJed Brown PetscFinalize(); 152c4762a1bSJed Brown return 0; 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 155c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 156c4762a1bSJed Brown /* 157c4762a1bSJed Brown FormInitialGuess - Computes an initial approximation to the solution. 158c4762a1bSJed Brown 159c4762a1bSJed Brown Input Parameters: 160c4762a1bSJed Brown . user - user-defined application context 161c4762a1bSJed Brown . X - vector 162c4762a1bSJed Brown 163c4762a1bSJed Brown Output Parameters: 164c4762a1bSJed Brown X - vector 165c4762a1bSJed Brown */ 166c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx *user,Vec X) 167c4762a1bSJed Brown { 168c4762a1bSJed Brown PetscInt i, j, k, mx = user->mx, my = user->my; 169c4762a1bSJed Brown PetscInt xs, ys, xm, ym, gxm, gym, gxs, gys, xe, ye; 170c4762a1bSJed Brown PetscReal hx = 1.0/(mx+1), hy = 1.0/(my+1), temp, val; 171c4762a1bSJed Brown 172c4762a1bSJed Brown PetscFunctionBegin; 173c4762a1bSJed Brown /* Get local mesh boundaries */ 174*9566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 175*9566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL)); 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* Compute initial guess over locally owned part of mesh */ 178c4762a1bSJed Brown xe = xs+xm; 179c4762a1bSJed Brown ye = ys+ym; 180c4762a1bSJed Brown for (j=ys; j<ye; j++) { /* for (j=0; j<my; j++) */ 181c4762a1bSJed Brown temp = PetscMin(j+1,my-j)*hy; 182c4762a1bSJed Brown for (i=xs; i<xe; i++) { /* for (i=0; i<mx; i++) */ 183c4762a1bSJed Brown k = (j-gys)*gxm + i-gxs; 184c4762a1bSJed Brown val = PetscMin((PetscMin(i+1,mx-i))*hx,temp); 185*9566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(X,1,&k,&val,ADD_VALUES)); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown } 188*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X)); 189*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X)); 190c4762a1bSJed Brown PetscFunctionReturn(0); 191c4762a1bSJed Brown } 192c4762a1bSJed Brown 193c4762a1bSJed Brown /* ------------------------------------------------------------------ */ 194c4762a1bSJed Brown /* 195c4762a1bSJed Brown FormFunctionGradient - Evaluates the function and corresponding gradient. 196c4762a1bSJed Brown 197c4762a1bSJed Brown Input Parameters: 198c4762a1bSJed Brown tao - the Tao context 199c4762a1bSJed Brown X - the input vector 200a82e8c82SStefano Zampini ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient() 201c4762a1bSJed Brown 202c4762a1bSJed Brown Output Parameters: 203c4762a1bSJed Brown f - the newly evaluated function 204c4762a1bSJed Brown G - the newly evaluated gradient 205c4762a1bSJed Brown */ 20670a7d78aSStefano Zampini PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f,Vec G,void *ptr) 20770a7d78aSStefano Zampini { 208c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 209c4762a1bSJed Brown PetscInt i,j,k,ind; 210c4762a1bSJed Brown PetscInt xe,ye,xsm,ysm,xep,yep; 211c4762a1bSJed Brown PetscInt xs, ys, xm, ym, gxm, gym, gxs, gys; 212c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my; 213c4762a1bSJed Brown PetscReal three = 3.0, zero = 0.0, *x, floc, cdiv3 = user->param/three; 214c4762a1bSJed Brown PetscReal p5 = 0.5, area, val, flin, fquad; 215c4762a1bSJed Brown PetscReal v,vb,vl,vr,vt,dvdx,dvdy; 216c4762a1bSJed Brown PetscReal hx = 1.0/(user->mx + 1); 217c4762a1bSJed Brown PetscReal hy = 1.0/(user->my + 1); 218c4762a1bSJed Brown Vec localX = user->localX; 219c4762a1bSJed Brown 220c4762a1bSJed Brown PetscFunctionBegin; 221c4762a1bSJed Brown /* Initialize */ 222c4762a1bSJed Brown flin = fquad = zero; 223c4762a1bSJed Brown 224*9566063dSJacob Faibussowitsch PetscCall(VecSet(G, zero)); 225c4762a1bSJed Brown /* 226c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 227c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 228c4762a1bSJed Brown By placing code between these two statements, computations can be 229c4762a1bSJed Brown done while messages are in transition. 230c4762a1bSJed Brown */ 231*9566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX)); 232*9566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX)); 233c4762a1bSJed Brown 234c4762a1bSJed Brown /* Get pointer to vector data */ 235*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX,&x)); 236c4762a1bSJed Brown 237c4762a1bSJed Brown /* Get local mesh boundaries */ 238*9566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 239*9566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL)); 240c4762a1bSJed Brown 241c4762a1bSJed Brown /* Set local loop dimensions */ 242c4762a1bSJed Brown xe = xs+xm; 243c4762a1bSJed Brown ye = ys+ym; 244c4762a1bSJed Brown if (xs == 0) xsm = xs-1; 245c4762a1bSJed Brown else xsm = xs; 246c4762a1bSJed Brown if (ys == 0) ysm = ys-1; 247c4762a1bSJed Brown else ysm = ys; 248c4762a1bSJed Brown if (xe == mx) xep = xe+1; 249c4762a1bSJed Brown else xep = xe; 250c4762a1bSJed Brown if (ye == my) yep = ye+1; 251c4762a1bSJed Brown else yep = ye; 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* Compute local gradient contributions over the lower triangular elements */ 254c4762a1bSJed Brown for (j=ysm; j<ye; j++) { /* for (j=-1; j<my; j++) */ 255c4762a1bSJed Brown for (i=xsm; i<xe; i++) { /* for (i=-1; i<mx; i++) */ 256c4762a1bSJed Brown k = (j-gys)*gxm + i-gxs; 257c4762a1bSJed Brown v = zero; 258c4762a1bSJed Brown vr = zero; 259c4762a1bSJed Brown vt = zero; 260c4762a1bSJed Brown if (i >= 0 && j >= 0) v = x[k]; 261c4762a1bSJed Brown if (i < mx-1 && j > -1) vr = x[k+1]; 262c4762a1bSJed Brown if (i > -1 && j < my-1) vt = x[k+gxm]; 263c4762a1bSJed Brown dvdx = (vr-v)/hx; 264c4762a1bSJed Brown dvdy = (vt-v)/hy; 265c4762a1bSJed Brown if (i != -1 && j != -1) { 266c4762a1bSJed Brown ind = k; val = - dvdx/hx - dvdy/hy - cdiv3; 267*9566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(G,1,&k,&val,ADD_VALUES)); 268c4762a1bSJed Brown } 269c4762a1bSJed Brown if (i != mx-1 && j != -1) { 270c4762a1bSJed Brown ind = k+1; val = dvdx/hx - cdiv3; 271*9566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES)); 272c4762a1bSJed Brown } 273c4762a1bSJed Brown if (i != -1 && j != my-1) { 274c4762a1bSJed Brown ind = k+gxm; val = dvdy/hy - cdiv3; 275*9566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES)); 276c4762a1bSJed Brown } 277c4762a1bSJed Brown fquad += dvdx*dvdx + dvdy*dvdy; 278c4762a1bSJed Brown flin -= cdiv3 * (v + vr + vt); 279c4762a1bSJed Brown } 280c4762a1bSJed Brown } 281c4762a1bSJed Brown 282c4762a1bSJed Brown /* Compute local gradient contributions over the upper triangular elements */ 283c4762a1bSJed Brown for (j=ys; j<yep; j++) { /* for (j=0; j<=my; j++) */ 284c4762a1bSJed Brown for (i=xs; i<xep; i++) { /* for (i=0; i<=mx; i++) */ 285c4762a1bSJed Brown k = (j-gys)*gxm + i-gxs; 286c4762a1bSJed Brown vb = zero; 287c4762a1bSJed Brown vl = zero; 288c4762a1bSJed Brown v = zero; 289c4762a1bSJed Brown if (i < mx && j > 0) vb = x[k-gxm]; 290c4762a1bSJed Brown if (i > 0 && j < my) vl = x[k-1]; 291c4762a1bSJed Brown if (i < mx && j < my) v = x[k]; 292c4762a1bSJed Brown dvdx = (v-vl)/hx; 293c4762a1bSJed Brown dvdy = (v-vb)/hy; 294c4762a1bSJed Brown if (i != mx && j != 0) { 295c4762a1bSJed Brown ind = k-gxm; val = - dvdy/hy - cdiv3; 296*9566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES)); 297c4762a1bSJed Brown } 298c4762a1bSJed Brown if (i != 0 && j != my) { 299c4762a1bSJed Brown ind = k-1; val = - dvdx/hx - cdiv3; 300*9566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES)); 301c4762a1bSJed Brown } 302c4762a1bSJed Brown if (i != mx && j != my) { 303c4762a1bSJed Brown ind = k; val = dvdx/hx + dvdy/hy - cdiv3; 304*9566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(G,1,&ind,&val,ADD_VALUES)); 305c4762a1bSJed Brown } 306c4762a1bSJed Brown fquad += dvdx*dvdx + dvdy*dvdy; 307c4762a1bSJed Brown flin -= cdiv3 * (vb + vl + v); 308c4762a1bSJed Brown } 309c4762a1bSJed Brown } 310c4762a1bSJed Brown 311c4762a1bSJed Brown /* Restore vector */ 312*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX,&x)); 313c4762a1bSJed Brown 314c4762a1bSJed Brown /* Assemble gradient vector */ 315*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(G)); 316*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(G)); 317c4762a1bSJed Brown 318c4762a1bSJed Brown /* Scale the gradient */ 319c4762a1bSJed Brown area = p5*hx*hy; 320c4762a1bSJed Brown floc = area * (p5 * fquad + flin); 321*9566063dSJacob Faibussowitsch PetscCall(VecScale(G, area)); 322c4762a1bSJed Brown 323*9566063dSJacob Faibussowitsch /* Sum function contributions from all processes */ /* TODO: Change to PetscCallMPI() */ 324*9566063dSJacob Faibussowitsch PetscCall((PetscErrorCode)MPI_Allreduce((void*)&floc,(void*)f,1,MPIU_REAL,MPIU_SUM,MPI_COMM_WORLD)); 325c4762a1bSJed Brown 326*9566063dSJacob Faibussowitsch PetscCall(PetscLogFlops((ye-ysm)*(xe-xsm)*20+(xep-xs)*(yep-ys)*16)); 327c4762a1bSJed Brown PetscFunctionReturn(0); 328c4762a1bSJed Brown } 329c4762a1bSJed Brown 330c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao, Vec X, Mat A, Mat Hpre, void*ctx) 331c4762a1bSJed Brown { 332c4762a1bSJed Brown AppCtx *user= (AppCtx*) ctx; 333c4762a1bSJed Brown PetscInt i,j,k; 334c4762a1bSJed Brown PetscInt col[5],row; 335c4762a1bSJed Brown PetscInt xs,xm,gxs,gxm,ys,ym,gys,gym; 336c4762a1bSJed Brown PetscReal v[5]; 337c4762a1bSJed Brown PetscReal hx=1.0/(user->mx+1), hy=1.0/(user->my+1), hxhx=1.0/(hx*hx), hyhy=1.0/(hy*hy), area=0.5*hx*hy; 338c4762a1bSJed Brown 339c4762a1bSJed Brown /* Compute the quadratic term in the objective function */ 340c4762a1bSJed Brown 341c4762a1bSJed Brown /* 342c4762a1bSJed Brown Get local grid boundaries 343c4762a1bSJed Brown */ 344c4762a1bSJed Brown 345c4762a1bSJed Brown PetscFunctionBegin; 346*9566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 347*9566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL)); 348c4762a1bSJed Brown 349c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 350c4762a1bSJed Brown 351c4762a1bSJed Brown for (i=xs; i< xs+xm; i++) { 352c4762a1bSJed Brown 353c4762a1bSJed Brown row=(j-gys)*gxm + (i-gxs); 354c4762a1bSJed Brown 355c4762a1bSJed Brown k=0; 356c4762a1bSJed Brown if (j>gys) { 357c4762a1bSJed Brown v[k]=-2*hyhy; col[k]=row - gxm; k++; 358c4762a1bSJed Brown } 359c4762a1bSJed Brown 360c4762a1bSJed Brown if (i>gxs) { 361c4762a1bSJed Brown v[k]= -2*hxhx; col[k]=row - 1; k++; 362c4762a1bSJed Brown } 363c4762a1bSJed Brown 364c4762a1bSJed Brown v[k]= 4.0*(hxhx+hyhy); col[k]=row; k++; 365c4762a1bSJed Brown 366c4762a1bSJed Brown if (i+1 < gxs+gxm) { 367c4762a1bSJed Brown v[k]= -2.0*hxhx; col[k]=row+1; k++; 368c4762a1bSJed Brown } 369c4762a1bSJed Brown 370c4762a1bSJed Brown if (j+1 <gys+gym) { 371c4762a1bSJed Brown v[k]= -2*hyhy; col[k] = row+gxm; k++; 372c4762a1bSJed Brown } 373c4762a1bSJed Brown 374*9566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A,1,&row,k,col,v,INSERT_VALUES)); 375c4762a1bSJed Brown 376c4762a1bSJed Brown } 377c4762a1bSJed Brown } 378c4762a1bSJed Brown /* 379c4762a1bSJed Brown Assemble matrix, using the 2-step process: 380c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd(). 381c4762a1bSJed Brown By placing code between these two statements, computations can be 382c4762a1bSJed Brown done while messages are in transition. 383c4762a1bSJed Brown */ 384*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 385*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 386c4762a1bSJed Brown /* 387c4762a1bSJed Brown Tell the matrix we will never add a new nonzero location to the 388c4762a1bSJed Brown matrix. If we do it will generate an error. 389c4762a1bSJed Brown */ 390*9566063dSJacob Faibussowitsch PetscCall(MatScale(A,area)); 391*9566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 392*9566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE)); 393*9566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(9*xm*ym+49*xm)); 394c4762a1bSJed Brown PetscFunctionReturn(0); 395c4762a1bSJed Brown } 396c4762a1bSJed Brown 397c4762a1bSJed Brown /*TEST 398c4762a1bSJed Brown 399c4762a1bSJed Brown build: 400c4762a1bSJed Brown requires: !complex 401c4762a1bSJed Brown 402c4762a1bSJed Brown test: 403c4762a1bSJed Brown suffix: 1 404c4762a1bSJed Brown nsize: 2 405c4762a1bSJed Brown args: -tao_smonitor -tao_type nls -mx 16 -my 16 -tao_gatol 1.e-4 406c4762a1bSJed Brown 407c4762a1bSJed Brown TEST*/ 408