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